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 5a7748839SVaclav Hapla const char * const DMPlexInterpolatedFlags[] = {"none", "partial", "mixed", "full", "DMPlexInterpolatedFlag", "DMPLEX_INTERPOLATED_", 0}; 6a7748839SVaclav Hapla 7e8f14785SLisandro Dalcin /* HashIJKL */ 8e8f14785SLisandro Dalcin 9e8f14785SLisandro Dalcin #include <petsc/private/hashmap.h> 10e8f14785SLisandro Dalcin 11e8f14785SLisandro Dalcin typedef struct _PetscHashIJKLKey { PetscInt i, j, k, l; } PetscHashIJKLKey; 12e8f14785SLisandro Dalcin 13e8f14785SLisandro Dalcin #define PetscHashIJKLKeyHash(key) \ 14e8f14785SLisandro Dalcin PetscHashCombine(PetscHashCombine(PetscHashInt((key).i),PetscHashInt((key).j)), \ 15e8f14785SLisandro Dalcin PetscHashCombine(PetscHashInt((key).k),PetscHashInt((key).l))) 16e8f14785SLisandro Dalcin 17e8f14785SLisandro Dalcin #define PetscHashIJKLKeyEqual(k1,k2) \ 18e8f14785SLisandro Dalcin (((k1).i==(k2).i) ? ((k1).j==(k2).j) ? ((k1).k==(k2).k) ? ((k1).l==(k2).l) : 0 : 0 : 0) 19e8f14785SLisandro Dalcin 20e8f14785SLisandro Dalcin PETSC_HASH_MAP(HashIJKL, PetscHashIJKLKey, PetscInt, PetscHashIJKLKeyHash, PetscHashIJKLKeyEqual, -1) 21e8f14785SLisandro Dalcin 223be36e83SMatthew G. Knepley static PetscSFNode _PetscInvalidSFNode = {-1, -1}; 233be36e83SMatthew G. Knepley 243be36e83SMatthew G. Knepley typedef struct _PetscHashIJKLRemoteKey { PetscSFNode i, j, k, l; } PetscHashIJKLRemoteKey; 253be36e83SMatthew G. Knepley 263be36e83SMatthew G. Knepley #define PetscHashIJKLRemoteKeyHash(key) \ 273be36e83SMatthew G. Knepley PetscHashCombine(PetscHashCombine(PetscHashInt((key).i.rank + (key).i.index),PetscHashInt((key).j.rank + (key).j.index)), \ 283be36e83SMatthew G. Knepley PetscHashCombine(PetscHashInt((key).k.rank + (key).k.index),PetscHashInt((key).l.rank + (key).l.index))) 293be36e83SMatthew G. Knepley 303be36e83SMatthew G. Knepley #define PetscHashIJKLRemoteKeyEqual(k1,k2) \ 313be36e83SMatthew G. Knepley (((k1).i.rank==(k2).i.rank) ? ((k1).i.index==(k2).i.index) ? ((k1).j.rank==(k2).j.rank) ? ((k1).j.index==(k2).j.index) ? ((k1).k.rank==(k2).k.rank) ? ((k1).k.index==(k2).k.index) ? ((k1).l.rank==(k2).l.rank) ? ((k1).l.index==(k2).l.index) : 0 : 0 : 0 : 0 : 0 : 0 : 0) 323be36e83SMatthew G. Knepley 333be36e83SMatthew G. Knepley PETSC_HASH_MAP(HashIJKLRemote, PetscHashIJKLRemoteKey, PetscSFNode, PetscHashIJKLRemoteKeyHash, PetscHashIJKLRemoteKeyEqual, _PetscInvalidSFNode) 343be36e83SMatthew G. Knepley 353be36e83SMatthew G. Knepley static PetscErrorCode PetscSortSFNode(PetscInt n, PetscSFNode A[]) 363be36e83SMatthew G. Knepley { 373be36e83SMatthew G. Knepley PetscInt i; 383be36e83SMatthew G. Knepley 393be36e83SMatthew G. Knepley PetscFunctionBegin; 403be36e83SMatthew G. Knepley for (i = 1; i < n; ++i) { 413be36e83SMatthew G. Knepley PetscSFNode x = A[i]; 423be36e83SMatthew G. Knepley PetscInt j; 433be36e83SMatthew G. Knepley 443be36e83SMatthew G. Knepley for (j = i-1; j >= 0; --j) { 453be36e83SMatthew G. Knepley if ((A[j].rank > x.rank) || (A[j].rank == x.rank && A[j].index > x.index)) break; 463be36e83SMatthew G. Knepley A[j+1] = A[j]; 473be36e83SMatthew G. Knepley } 483be36e83SMatthew G. Knepley A[j+1] = x; 493be36e83SMatthew G. Knepley } 503be36e83SMatthew G. Knepley PetscFunctionReturn(0); 513be36e83SMatthew G. Knepley } 529f074e33SMatthew G Knepley 539f074e33SMatthew G Knepley /* 54439ece16SMatthew G. Knepley DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone 55439ece16SMatthew G. Knepley */ 56*412e9a14SMatthew G. Knepley PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, DMPolytopeType ct, const PetscInt cone[], PetscInt *numFaces, const DMPolytopeType *faceTypes[], const PetscInt *faceSizes[], const PetscInt *faces[]) 57439ece16SMatthew G. Knepley { 58*412e9a14SMatthew G. Knepley DMPolytopeType *typesTmp; 59*412e9a14SMatthew G. Knepley PetscInt *sizesTmp, *facesTmp; 60439ece16SMatthew G. Knepley PetscInt maxConeSize, maxSupportSize; 61439ece16SMatthew G. Knepley PetscErrorCode ierr; 62439ece16SMatthew G. Knepley 63439ece16SMatthew G. Knepley PetscFunctionBegin; 64439ece16SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 65ba2698f1SMatthew G. Knepley if (cone) PetscValidIntPointer(cone, 3); 66439ece16SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 67*412e9a14SMatthew G. Knepley if (faceTypes) {ierr = DMGetWorkArray(dm, PetscMax(maxConeSize, maxSupportSize), MPIU_INT, &typesTmp);CHKERRQ(ierr);} 68*412e9a14SMatthew G. Knepley if (faceSizes) {ierr = DMGetWorkArray(dm, PetscMax(maxConeSize, maxSupportSize), MPIU_INT, &sizesTmp);CHKERRQ(ierr);} 6969291d52SBarry Smith if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);CHKERRQ(ierr);} 70ba2698f1SMatthew G. Knepley switch (ct) { 71ba2698f1SMatthew G. Knepley case DM_POLYTOPE_POINT: 72ba2698f1SMatthew G. Knepley if (numFaces) *numFaces = 0; 73ba2698f1SMatthew G. Knepley break; 74ba2698f1SMatthew G. Knepley case DM_POLYTOPE_SEGMENT: 75*412e9a14SMatthew G. Knepley if (numFaces) *numFaces = 2; 76*412e9a14SMatthew G. Knepley if (faceTypes) { 77*412e9a14SMatthew G. Knepley typesTmp[0] = DM_POLYTOPE_POINT; typesTmp[1] = DM_POLYTOPE_POINT; 78*412e9a14SMatthew G. Knepley *faceTypes = typesTmp; 79*412e9a14SMatthew G. Knepley } 80*412e9a14SMatthew G. Knepley if (faceSizes) { 81*412e9a14SMatthew G. Knepley sizesTmp[0] = 1; sizesTmp[1] = 1; 82*412e9a14SMatthew G. Knepley *faceSizes = sizesTmp; 83*412e9a14SMatthew G. Knepley } 84c49d9212SMatthew G. Knepley if (faces) { 85c49d9212SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 86c49d9212SMatthew G. Knepley *faces = facesTmp; 87c49d9212SMatthew G. Knepley } 88*412e9a14SMatthew G. Knepley break; 89*412e9a14SMatthew G. Knepley case DM_POLYTOPE_POINT_PRISM_TENSOR: 90c49d9212SMatthew G. Knepley if (numFaces) *numFaces = 2; 91*412e9a14SMatthew G. Knepley if (faceTypes) { 92*412e9a14SMatthew G. Knepley typesTmp[0] = DM_POLYTOPE_POINT; typesTmp[1] = DM_POLYTOPE_POINT; 93*412e9a14SMatthew G. Knepley *faceTypes = typesTmp; 94*412e9a14SMatthew G. Knepley } 95*412e9a14SMatthew G. Knepley if (faceSizes) { 96*412e9a14SMatthew G. Knepley sizesTmp[0] = 1; sizesTmp[1] = 1; 97*412e9a14SMatthew G. Knepley *faceSizes = sizesTmp; 98*412e9a14SMatthew G. Knepley } 99*412e9a14SMatthew G. Knepley if (faces) { 100*412e9a14SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 101*412e9a14SMatthew G. Knepley *faces = facesTmp; 102*412e9a14SMatthew G. Knepley } 103c49d9212SMatthew G. Knepley break; 104ba2698f1SMatthew G. Knepley case DM_POLYTOPE_TRIANGLE: 105*412e9a14SMatthew G. Knepley if (numFaces) *numFaces = 3; 106*412e9a14SMatthew G. Knepley if (faceTypes) { 107*412e9a14SMatthew G. Knepley typesTmp[0] = DM_POLYTOPE_SEGMENT; typesTmp[1] = DM_POLYTOPE_SEGMENT; typesTmp[2] = DM_POLYTOPE_SEGMENT; 108*412e9a14SMatthew G. Knepley *faceTypes = typesTmp; 109*412e9a14SMatthew G. Knepley } 110*412e9a14SMatthew G. Knepley if (faceSizes) { 111*412e9a14SMatthew G. Knepley sizesTmp[0] = 2; sizesTmp[1] = 2; sizesTmp[2] = 2; 112*412e9a14SMatthew G. Knepley *faceSizes = sizesTmp; 113*412e9a14SMatthew G. Knepley } 1149a5b13a2SMatthew G Knepley if (faces) { 1159a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 1169a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 1179a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 1189a5b13a2SMatthew G Knepley *faces = facesTmp; 1199a5b13a2SMatthew G Knepley } 1209f074e33SMatthew G Knepley break; 121ba2698f1SMatthew G. Knepley case DM_POLYTOPE_QUADRILATERAL: 1229a5b13a2SMatthew G Knepley /* Vertices follow right hand rule */ 123*412e9a14SMatthew G. Knepley if (numFaces) *numFaces = 4; 124*412e9a14SMatthew G. Knepley if (faceTypes) { 125*412e9a14SMatthew G. Knepley typesTmp[0] = DM_POLYTOPE_SEGMENT; typesTmp[1] = DM_POLYTOPE_SEGMENT; typesTmp[2] = DM_POLYTOPE_SEGMENT; typesTmp[3] = DM_POLYTOPE_SEGMENT; 126*412e9a14SMatthew G. Knepley *faceTypes = typesTmp; 127*412e9a14SMatthew G. Knepley } 128*412e9a14SMatthew G. Knepley if (faceSizes) { 129*412e9a14SMatthew G. Knepley sizesTmp[0] = 2; sizesTmp[1] = 2; sizesTmp[2] = 2; sizesTmp[3] = 2; 130*412e9a14SMatthew G. Knepley *faceSizes = sizesTmp; 131*412e9a14SMatthew G. Knepley } 1329a5b13a2SMatthew G Knepley if (faces) { 1339a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 1349a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 1359a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[3]; 1369a5b13a2SMatthew G Knepley facesTmp[6] = cone[3]; facesTmp[7] = cone[0]; 1379a5b13a2SMatthew G Knepley *faces = facesTmp; 1389a5b13a2SMatthew G Knepley } 139*412e9a14SMatthew G. Knepley break; 140*412e9a14SMatthew G. Knepley case DM_POLYTOPE_SEG_PRISM_TENSOR: 1419a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 4; 142*412e9a14SMatthew G. Knepley if (faceTypes) { 143*412e9a14SMatthew G. Knepley typesTmp[0] = DM_POLYTOPE_SEGMENT; typesTmp[1] = DM_POLYTOPE_SEGMENT; typesTmp[2] = DM_POLYTOPE_POINT_PRISM_TENSOR; typesTmp[3] = DM_POLYTOPE_POINT_PRISM_TENSOR; 144*412e9a14SMatthew G. Knepley *faceTypes = typesTmp; 145*412e9a14SMatthew G. Knepley } 146*412e9a14SMatthew G. Knepley if (faceSizes) { 147*412e9a14SMatthew G. Knepley sizesTmp[0] = 2; sizesTmp[1] = 2; sizesTmp[2] = 2; sizesTmp[3] = 2; 148*412e9a14SMatthew G. Knepley *faceSizes = sizesTmp; 149*412e9a14SMatthew G. Knepley } 150*412e9a14SMatthew G. Knepley if (faces) { 151*412e9a14SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 152*412e9a14SMatthew G. Knepley facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; 153*412e9a14SMatthew G. Knepley facesTmp[4] = cone[0]; facesTmp[5] = cone[2]; 154*412e9a14SMatthew G. Knepley facesTmp[6] = cone[1]; facesTmp[7] = cone[3]; 155*412e9a14SMatthew G. Knepley *faces = facesTmp; 156*412e9a14SMatthew G. Knepley } 1579f074e33SMatthew G Knepley break; 158ba2698f1SMatthew G. Knepley case DM_POLYTOPE_TETRAHEDRON: 1592e1b13c2SMatthew G. Knepley /* Vertices of first face follow right hand rule and normal points away from last vertex */ 160*412e9a14SMatthew G. Knepley if (numFaces) *numFaces = 4; 161*412e9a14SMatthew G. Knepley if (faceTypes) { 162*412e9a14SMatthew G. Knepley typesTmp[0] = DM_POLYTOPE_TRIANGLE; typesTmp[1] = DM_POLYTOPE_TRIANGLE; typesTmp[2] = DM_POLYTOPE_TRIANGLE; typesTmp[3] = DM_POLYTOPE_TRIANGLE; 163*412e9a14SMatthew G. Knepley *faceTypes = typesTmp; 164*412e9a14SMatthew G. Knepley } 165*412e9a14SMatthew G. Knepley if (faceSizes) { 166*412e9a14SMatthew G. Knepley sizesTmp[0] = 3; sizesTmp[1] = 3; sizesTmp[2] = 3; sizesTmp[3] = 3; 167*412e9a14SMatthew G. Knepley *faceSizes = sizesTmp; 168*412e9a14SMatthew G. Knepley } 1699a5b13a2SMatthew G Knepley if (faces) { 1702e1b13c2SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; 1712e1b13c2SMatthew G. Knepley facesTmp[3] = cone[0]; facesTmp[4] = cone[3]; facesTmp[5] = cone[1]; 1722e1b13c2SMatthew G. Knepley facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[3]; 1732e1b13c2SMatthew G. Knepley facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3]; 1749a5b13a2SMatthew G Knepley *faces = facesTmp; 1759a5b13a2SMatthew G Knepley } 1769a5b13a2SMatthew G Knepley break; 177ba2698f1SMatthew G. Knepley case DM_POLYTOPE_HEXAHEDRON: 178e368ef20SMatthew G. Knepley /* 7--------6 179e368ef20SMatthew G. Knepley /| /| 180e368ef20SMatthew G. Knepley / | / | 181e368ef20SMatthew G. Knepley 4--------5 | 182e368ef20SMatthew G. Knepley | | | | 183e368ef20SMatthew G. Knepley | | | | 184e368ef20SMatthew G. Knepley | 1--------2 185e368ef20SMatthew G. Knepley | / | / 186e368ef20SMatthew G. Knepley |/ |/ 187e368ef20SMatthew G. Knepley 0--------3 188e368ef20SMatthew G. Knepley */ 189*412e9a14SMatthew G. Knepley if (numFaces) *numFaces = 6; 190*412e9a14SMatthew G. Knepley if (faceTypes) { 191*412e9a14SMatthew G. Knepley typesTmp[0] = DM_POLYTOPE_QUADRILATERAL; typesTmp[1] = DM_POLYTOPE_QUADRILATERAL; typesTmp[2] = DM_POLYTOPE_QUADRILATERAL; 192*412e9a14SMatthew G. Knepley typesTmp[3] = DM_POLYTOPE_QUADRILATERAL; typesTmp[4] = DM_POLYTOPE_QUADRILATERAL; typesTmp[5] = DM_POLYTOPE_QUADRILATERAL; 193*412e9a14SMatthew G. Knepley *faceTypes = typesTmp; 194*412e9a14SMatthew G. Knepley } 195*412e9a14SMatthew G. Knepley if (faceSizes) { 196*412e9a14SMatthew G. Knepley sizesTmp[0] = 4; sizesTmp[1] = 4; sizesTmp[2] = 4; sizesTmp[3] = 4; sizesTmp[4] = 4; sizesTmp[5] = 4; 197*412e9a14SMatthew G. Knepley *faceSizes = sizesTmp; 198*412e9a14SMatthew G. Knepley } 1999a5b13a2SMatthew G Knepley if (faces) { 20051cfd6a4SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; /* Bottom */ 20151cfd6a4SMatthew G. Knepley facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; /* Top */ 20251cfd6a4SMatthew G. Knepley facesTmp[8] = cone[0]; facesTmp[9] = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; /* Front */ 20351cfd6a4SMatthew G. Knepley facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; /* Back */ 20451cfd6a4SMatthew G. Knepley facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; /* Right */ 20551cfd6a4SMatthew G. Knepley facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; /* Left */ 2069a5b13a2SMatthew G Knepley *faces = facesTmp; 2079a5b13a2SMatthew G Knepley } 20899836e0eSStefano Zampini break; 209ba2698f1SMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM: 210*412e9a14SMatthew G. Knepley if (numFaces) *numFaces = 5; 211*412e9a14SMatthew G. Knepley if (faceTypes) { 212*412e9a14SMatthew G. Knepley typesTmp[0] = DM_POLYTOPE_TRIANGLE; typesTmp[1] = DM_POLYTOPE_TRIANGLE; 213*412e9a14SMatthew G. Knepley typesTmp[2] = DM_POLYTOPE_QUADRILATERAL; typesTmp[3] = DM_POLYTOPE_QUADRILATERAL; typesTmp[4] = DM_POLYTOPE_QUADRILATERAL; 214*412e9a14SMatthew G. Knepley *faceTypes = typesTmp; 215*412e9a14SMatthew G. Knepley } 216*412e9a14SMatthew G. Knepley if (faceSizes) { 217*412e9a14SMatthew G. Knepley sizesTmp[0] = 3; sizesTmp[1] = 3; 218*412e9a14SMatthew G. Knepley sizesTmp[2] = 4; sizesTmp[3] = 4; sizesTmp[4] = 4; 219*412e9a14SMatthew G. Knepley *faceSizes = sizesTmp; 220*412e9a14SMatthew G. Knepley } 221ba2698f1SMatthew G. Knepley if (faces) { 222*412e9a14SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; /* Bottom */ 223*412e9a14SMatthew G. Knepley facesTmp[3] = cone[3]; facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; /* Top */ 224*412e9a14SMatthew G. Knepley facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[4]; facesTmp[9] = cone[3]; /* Back left */ 225*412e9a14SMatthew G. Knepley facesTmp[10] = cone[2]; facesTmp[11] = cone[1]; facesTmp[12] = cone[5]; facesTmp[13] = cone[4]; /* Front */ 226*412e9a14SMatthew G. Knepley facesTmp[14] = cone[1]; facesTmp[15] = cone[0]; facesTmp[16] = cone[3]; facesTmp[17] = cone[5]; /* Back right */ 227ba2698f1SMatthew G. Knepley *faces = facesTmp; 22899836e0eSStefano Zampini } 22999836e0eSStefano Zampini break; 230ba2698f1SMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM_TENSOR: 231*412e9a14SMatthew G. Knepley if (numFaces) *numFaces = 5; 232*412e9a14SMatthew G. Knepley if (faceTypes) { 233*412e9a14SMatthew G. Knepley typesTmp[0] = DM_POLYTOPE_TRIANGLE; typesTmp[1] = DM_POLYTOPE_TRIANGLE; 234*412e9a14SMatthew G. Knepley typesTmp[2] = DM_POLYTOPE_SEG_PRISM_TENSOR; typesTmp[3] = DM_POLYTOPE_SEG_PRISM_TENSOR; typesTmp[4] = DM_POLYTOPE_SEG_PRISM_TENSOR; 235*412e9a14SMatthew G. Knepley *faceTypes = typesTmp; 236*412e9a14SMatthew G. Knepley } 237*412e9a14SMatthew G. Knepley if (faceSizes) { 238*412e9a14SMatthew G. Knepley sizesTmp[0] = 3; sizesTmp[1] = 3; 239*412e9a14SMatthew G. Knepley sizesTmp[2] = 4; sizesTmp[3] = 4; sizesTmp[4] = 4; 240*412e9a14SMatthew G. Knepley *faceSizes = sizesTmp; 241*412e9a14SMatthew G. Knepley } 24299836e0eSStefano Zampini if (faces) { 243*412e9a14SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; /* Bottom */ 244*412e9a14SMatthew G. Knepley facesTmp[3] = cone[3]; facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; /* Top */ 245*412e9a14SMatthew G. Knepley facesTmp[6] = cone[0]; facesTmp[7] = cone[1]; facesTmp[8] = cone[3]; facesTmp[9] = cone[4]; /* Back left */ 246*412e9a14SMatthew G. Knepley facesTmp[10] = cone[1]; facesTmp[11] = cone[2]; facesTmp[12] = cone[4]; facesTmp[13] = cone[5]; /* Back right */ 247*412e9a14SMatthew G. Knepley facesTmp[14] = cone[2]; facesTmp[15] = cone[0]; facesTmp[16] = cone[5]; facesTmp[17] = cone[3]; /* Front */ 24899836e0eSStefano Zampini *faces = facesTmp; 24999836e0eSStefano Zampini } 250*412e9a14SMatthew G. Knepley break; 251*412e9a14SMatthew G. Knepley case DM_POLYTOPE_QUAD_PRISM_TENSOR: 252*412e9a14SMatthew G. Knepley /* 7--------6 253*412e9a14SMatthew G. Knepley /| /| 254*412e9a14SMatthew G. Knepley / | / | 255*412e9a14SMatthew G. Knepley 4--------5 | 256*412e9a14SMatthew G. Knepley | | | | 257*412e9a14SMatthew G. Knepley | | | | 258*412e9a14SMatthew G. Knepley | 3--------2 259*412e9a14SMatthew G. Knepley | / | / 260*412e9a14SMatthew G. Knepley |/ |/ 261*412e9a14SMatthew G. Knepley 0--------1 262*412e9a14SMatthew G. Knepley */ 263*412e9a14SMatthew G. Knepley if (numFaces) *numFaces = 6; 264*412e9a14SMatthew G. Knepley if (faceTypes) { 265*412e9a14SMatthew G. Knepley typesTmp[0] = DM_POLYTOPE_QUADRILATERAL; typesTmp[1] = DM_POLYTOPE_QUADRILATERAL; typesTmp[2] = DM_POLYTOPE_SEG_PRISM_TENSOR; 266*412e9a14SMatthew G. Knepley typesTmp[3] = DM_POLYTOPE_SEG_PRISM_TENSOR; typesTmp[4] = DM_POLYTOPE_SEG_PRISM_TENSOR; typesTmp[5] = DM_POLYTOPE_SEG_PRISM_TENSOR; 267*412e9a14SMatthew G. Knepley *faceTypes = typesTmp; 268*412e9a14SMatthew G. Knepley } 269*412e9a14SMatthew G. Knepley if (faceSizes) { 270*412e9a14SMatthew G. Knepley sizesTmp[0] = 4; sizesTmp[1] = 4; sizesTmp[2] = 4; sizesTmp[3] = 4; sizesTmp[4] = 4; sizesTmp[5] = 4; 271*412e9a14SMatthew G. Knepley *faceSizes = sizesTmp; 272*412e9a14SMatthew G. Knepley } 273*412e9a14SMatthew G. Knepley if (faces) { 274*412e9a14SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; /* Bottom */ 275*412e9a14SMatthew G. Knepley facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; /* Top */ 276*412e9a14SMatthew G. Knepley facesTmp[8] = cone[0]; facesTmp[9] = cone[1]; facesTmp[10] = cone[4]; facesTmp[11] = cone[5]; /* Front */ 277*412e9a14SMatthew G. Knepley facesTmp[12] = cone[1]; facesTmp[13] = cone[2]; facesTmp[14] = cone[5]; facesTmp[15] = cone[6]; /* Right */ 278*412e9a14SMatthew G. Knepley facesTmp[16] = cone[2]; facesTmp[17] = cone[3]; facesTmp[18] = cone[6]; facesTmp[19] = cone[7]; /* Back */ 279*412e9a14SMatthew G. Knepley facesTmp[20] = cone[3]; facesTmp[21] = cone[0]; facesTmp[22] = cone[7]; facesTmp[23] = cone[4]; /* Left */ 280*412e9a14SMatthew G. Knepley *faces = facesTmp; 281*412e9a14SMatthew G. Knepley } 28299836e0eSStefano Zampini break; 283ba2698f1SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No face description for cell type %s", DMPolytopeTypes[ct]); 28499836e0eSStefano Zampini } 28599836e0eSStefano Zampini PetscFunctionReturn(0); 28699836e0eSStefano Zampini } 28799836e0eSStefano Zampini 288*412e9a14SMatthew G. Knepley PetscErrorCode DMPlexRestoreRawFaces_Internal(DM dm, DMPolytopeType ct, const PetscInt cone[], PetscInt *numFaces, const DMPolytopeType *faceTypes[], const PetscInt *faceSizes[], const PetscInt *faces[]) 28999836e0eSStefano Zampini { 29099836e0eSStefano Zampini PetscErrorCode ierr; 29199836e0eSStefano Zampini 29299836e0eSStefano Zampini PetscFunctionBegin; 293*412e9a14SMatthew G. Knepley if (faceTypes) {ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faceTypes);CHKERRQ(ierr);} 294*412e9a14SMatthew G. Knepley if (faceSizes) {ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faceSizes);CHKERRQ(ierr);} 29599836e0eSStefano Zampini if (faces) {ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);CHKERRQ(ierr);} 29699836e0eSStefano Zampini PetscFunctionReturn(0); 29799836e0eSStefano Zampini } 29899836e0eSStefano Zampini 2999a5b13a2SMatthew G Knepley /* This interpolates faces for cells at some stratum */ 3009a5b13a2SMatthew G Knepley static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm) 3019f074e33SMatthew G Knepley { 302*412e9a14SMatthew G. Knepley DMLabel ctLabel; 3039a5b13a2SMatthew G Knepley PetscHashIJKL faceTable; 304*412e9a14SMatthew G. Knepley PetscInt faceTypeNum[DM_NUM_POLYTOPES]; 305*412e9a14SMatthew G. Knepley PetscInt depth, d, Np, cStart, cEnd, c, fStart, fEnd; 3069f074e33SMatthew G Knepley PetscErrorCode ierr; 3079f074e33SMatthew G Knepley 3089f074e33SMatthew G Knepley PetscFunctionBegin; 3099a5b13a2SMatthew G Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 31099836e0eSStefano Zampini ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 311*412e9a14SMatthew G. Knepley ierr = PetscArrayzero(faceTypeNum, DM_NUM_POLYTOPES);CHKERRQ(ierr); 312*412e9a14SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, cellDepth, &cStart, &cEnd);CHKERRQ(ierr); 313*412e9a14SMatthew G. Knepley /* Number new faces and save face vertices in hash table */ 314*412e9a14SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, depth > cellDepth ? cellDepth : 0, NULL, &fStart);CHKERRQ(ierr); 315*412e9a14SMatthew G. Knepley fEnd = fStart; 316*412e9a14SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 317*412e9a14SMatthew G. Knepley const PetscInt *cone, *faceSizes, *faces; 318*412e9a14SMatthew G. Knepley const DMPolytopeType *faceTypes; 319ba2698f1SMatthew G. Knepley DMPolytopeType ct; 320*412e9a14SMatthew G. Knepley PetscInt numFaces, cf, foff = 0; 32199836e0eSStefano Zampini 322ba2698f1SMatthew G. Knepley ierr = DMPlexGetCellType(dm, c, &ct);CHKERRQ(ierr); 32399836e0eSStefano Zampini ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 324*412e9a14SMatthew G. Knepley ierr = DMPlexGetRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);CHKERRQ(ierr); 325*412e9a14SMatthew G. Knepley for (cf = 0; cf < numFaces; foff += faceSizes[cf], ++cf) { 326*412e9a14SMatthew G. Knepley const PetscInt faceSize = faceSizes[cf]; 327*412e9a14SMatthew G. Knepley const DMPolytopeType faceType = faceTypes[cf]; 328*412e9a14SMatthew G. Knepley const PetscInt *face = &faces[foff]; 3299a5b13a2SMatthew G Knepley PetscHashIJKLKey key; 330e8f14785SLisandro Dalcin PetscHashIter iter; 331e8f14785SLisandro Dalcin PetscBool missing; 3329f074e33SMatthew G Knepley 333*412e9a14SMatthew G. Knepley if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not support faces of size %D > 4", faceSize); 334*412e9a14SMatthew G. Knepley key.i = face[0]; 335*412e9a14SMatthew G. Knepley key.j = faceSize > 1 ? face[1] : PETSC_MAX_INT; 336*412e9a14SMatthew G. Knepley key.k = faceSize > 2 ? face[2] : PETSC_MAX_INT; 337*412e9a14SMatthew G. Knepley key.l = faceSize > 3 ? face[3] : PETSC_MAX_INT; 338302440fdSBarry Smith ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 339e8f14785SLisandro Dalcin ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 340e9fa77a1SMatthew G. Knepley if (missing) { 341*412e9a14SMatthew G. Knepley ierr = PetscHashIJKLIterSet(faceTable, iter, fEnd++);CHKERRQ(ierr); 342*412e9a14SMatthew G. Knepley ++faceTypeNum[faceType]; 343e9fa77a1SMatthew G. Knepley } 3449a5b13a2SMatthew G Knepley } 345*412e9a14SMatthew G. Knepley ierr = DMPlexRestoreRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);CHKERRQ(ierr); 34699836e0eSStefano Zampini } 347*412e9a14SMatthew G. Knepley /* We need to number faces contiguously among types */ 348*412e9a14SMatthew G. Knepley { 349*412e9a14SMatthew G. Knepley PetscInt faceTypeStart[DM_NUM_POLYTOPES], ct, numFT = 0; 35099836e0eSStefano Zampini 351*412e9a14SMatthew G. Knepley for (ct = 0; ct < DM_NUM_POLYTOPES; ++ct) {if (faceTypeNum[ct]) ++numFT; faceTypeStart[ct] = 0;} 352*412e9a14SMatthew G. Knepley if (numFT > 1) { 353*412e9a14SMatthew G. Knepley ierr = PetscHashIJKLClear(faceTable);CHKERRQ(ierr); 354*412e9a14SMatthew G. Knepley faceTypeStart[0] = fStart; 355*412e9a14SMatthew G. Knepley for (ct = 1; ct < DM_NUM_POLYTOPES; ++ct) faceTypeStart[ct] = faceTypeStart[ct-1] + faceTypeNum[ct-1]; 356*412e9a14SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 357*412e9a14SMatthew G. Knepley const PetscInt *cone, *faceSizes, *faces; 358*412e9a14SMatthew G. Knepley const DMPolytopeType *faceTypes; 359ba2698f1SMatthew G. Knepley DMPolytopeType ct; 360*412e9a14SMatthew G. Knepley PetscInt numFaces, cf, foff = 0; 36199836e0eSStefano Zampini 362ba2698f1SMatthew G. Knepley ierr = DMPlexGetCellType(dm, c, &ct);CHKERRQ(ierr); 36399836e0eSStefano Zampini ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 364*412e9a14SMatthew G. Knepley ierr = DMPlexGetRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);CHKERRQ(ierr); 365*412e9a14SMatthew G. Knepley for (cf = 0; cf < numFaces; foff += faceSizes[cf], ++cf) { 366*412e9a14SMatthew G. Knepley const PetscInt faceSize = faceSizes[cf]; 367*412e9a14SMatthew G. Knepley const DMPolytopeType faceType = faceTypes[cf]; 368*412e9a14SMatthew G. Knepley const PetscInt *face = &faces[foff]; 36999836e0eSStefano Zampini PetscHashIJKLKey key; 37099836e0eSStefano Zampini PetscHashIter iter; 37199836e0eSStefano Zampini PetscBool missing; 37299836e0eSStefano Zampini 373*412e9a14SMatthew G. Knepley if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not support faces of size %D > 4", faceSize); 374*412e9a14SMatthew G. Knepley key.i = face[0]; 375*412e9a14SMatthew G. Knepley key.j = faceSize > 1 ? face[1] : PETSC_MAX_INT; 376*412e9a14SMatthew G. Knepley key.k = faceSize > 2 ? face[2] : PETSC_MAX_INT; 377*412e9a14SMatthew G. Knepley key.l = faceSize > 3 ? face[3] : PETSC_MAX_INT; 37899836e0eSStefano Zampini ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 37999836e0eSStefano Zampini ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 380*412e9a14SMatthew G. Knepley if (missing) {ierr = PetscHashIJKLIterSet(faceTable, iter, faceTypeStart[faceType]++);CHKERRQ(ierr);} 38199836e0eSStefano Zampini } 382*412e9a14SMatthew G. Knepley ierr = DMPlexRestoreRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);CHKERRQ(ierr); 38399836e0eSStefano Zampini } 384*412e9a14SMatthew G. Knepley for (ct = 1; ct < DM_NUM_POLYTOPES; ++ct) { 385*412e9a14SMatthew G. Knepley if (faceTypeStart[ct] != faceTypeStart[ct-1] + faceTypeNum[ct]) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent numbering for cell type %s, %D != %D + %D", DMPolytopeTypes[ct], faceTypeStart[ct], faceTypeStart[ct-1], faceTypeNum[ct]); 3869a5b13a2SMatthew G Knepley } 3879a5b13a2SMatthew G Knepley } 388*412e9a14SMatthew G. Knepley } 389*412e9a14SMatthew G. Knepley /* Add new points, always at the end of the numbering */ 390*412e9a14SMatthew G. Knepley ierr = DMPlexGetChart(dm, NULL, &Np);CHKERRQ(ierr); 391*412e9a14SMatthew G. Knepley ierr = DMPlexSetChart(idm, 0, Np + (fEnd - fStart));CHKERRQ(ierr); 3929a5b13a2SMatthew G Knepley /* Set cone sizes */ 393*412e9a14SMatthew G. Knepley /* Must create the celltype label here so that we do not automatically try to compute the types */ 394*412e9a14SMatthew G. Knepley ierr = DMCreateLabel(idm, "celltype");CHKERRQ(ierr); 395*412e9a14SMatthew G. Knepley ierr = DMPlexGetCellTypeLabel(idm, &ctLabel);CHKERRQ(ierr); 3969a5b13a2SMatthew G Knepley for (d = 0; d <= depth; ++d) { 397*412e9a14SMatthew G. Knepley DMPolytopeType ct; 398*412e9a14SMatthew G. Knepley PetscInt coneSize, pStart, pEnd, p; 3999f074e33SMatthew G Knepley 400*412e9a14SMatthew G. Knepley if (d == cellDepth) continue; 401*412e9a14SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 402*412e9a14SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 4039a5b13a2SMatthew G Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 4049a5b13a2SMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 405*412e9a14SMatthew G. Knepley ierr = DMPlexGetCellType(dm, p, &ct);CHKERRQ(ierr); 406*412e9a14SMatthew G. Knepley ierr = DMPlexSetCellType(idm, p, ct);CHKERRQ(ierr); 4079f074e33SMatthew G Knepley } 4089f074e33SMatthew G Knepley } 409*412e9a14SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 410*412e9a14SMatthew G. Knepley const PetscInt *cone, *faceSizes, *faces; 411*412e9a14SMatthew G. Knepley const DMPolytopeType *faceTypes; 412*412e9a14SMatthew G. Knepley DMPolytopeType ct; 413*412e9a14SMatthew G. Knepley PetscInt numFaces, cf, foff = 0; 414*412e9a14SMatthew G. Knepley 415*412e9a14SMatthew G. Knepley ierr = DMPlexGetCellType(dm, c, &ct);CHKERRQ(ierr); 416*412e9a14SMatthew G. Knepley ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 417*412e9a14SMatthew G. Knepley ierr = DMPlexGetRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);CHKERRQ(ierr); 418*412e9a14SMatthew G. Knepley ierr = DMPlexSetCellType(idm, c, ct);CHKERRQ(ierr); 419*412e9a14SMatthew G. Knepley ierr = DMPlexSetConeSize(idm, c, numFaces);CHKERRQ(ierr); 420*412e9a14SMatthew G. Knepley for (cf = 0; cf < numFaces; foff += faceSizes[cf], ++cf) { 421*412e9a14SMatthew G. Knepley const PetscInt faceSize = faceSizes[cf]; 422*412e9a14SMatthew G. Knepley const DMPolytopeType faceType = faceTypes[cf]; 423*412e9a14SMatthew G. Knepley const PetscInt *face = &faces[foff]; 424*412e9a14SMatthew G. Knepley PetscHashIJKLKey key; 425*412e9a14SMatthew G. Knepley PetscHashIter iter; 426*412e9a14SMatthew G. Knepley PetscBool missing; 427*412e9a14SMatthew G. Knepley PetscInt f; 428*412e9a14SMatthew G. Knepley 429*412e9a14SMatthew G. Knepley if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not support faces of size %D > 4", faceSize); 430*412e9a14SMatthew G. Knepley key.i = face[0]; 431*412e9a14SMatthew G. Knepley key.j = faceSize > 1 ? face[1] : PETSC_MAX_INT; 432*412e9a14SMatthew G. Knepley key.k = faceSize > 2 ? face[2] : PETSC_MAX_INT; 433*412e9a14SMatthew G. Knepley key.l = faceSize > 3 ? face[3] : PETSC_MAX_INT; 434*412e9a14SMatthew G. Knepley ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 435*412e9a14SMatthew G. Knepley ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 436*412e9a14SMatthew G. Knepley if (missing) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing face (cell %D, lf %D)", c, cf); 437*412e9a14SMatthew G. Knepley ierr = PetscHashIJKLIterGet(faceTable, iter, &f);CHKERRQ(ierr); 438*412e9a14SMatthew G. Knepley ierr = DMPlexSetConeSize(idm, f, faceSize);CHKERRQ(ierr); 439*412e9a14SMatthew G. Knepley ierr = DMPlexSetCellType(idm, f, faceType);CHKERRQ(ierr); 440*412e9a14SMatthew G. Knepley } 441*412e9a14SMatthew G. Knepley ierr = DMPlexRestoreRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);CHKERRQ(ierr); 4429f074e33SMatthew G Knepley } 4439f074e33SMatthew G Knepley ierr = DMSetUp(idm);CHKERRQ(ierr); 444*412e9a14SMatthew G. Knepley /* Initialize cones so we do not need the bash table to tell us that a cone has been set */ 445*412e9a14SMatthew G. Knepley { 446*412e9a14SMatthew G. Knepley PetscSection cs; 447*412e9a14SMatthew G. Knepley PetscInt *cones, csize; 4489a5b13a2SMatthew G Knepley 449*412e9a14SMatthew G. Knepley ierr = DMPlexGetConeSection(idm, &cs);CHKERRQ(ierr); 450*412e9a14SMatthew G. Knepley ierr = DMPlexGetCones(idm, &cones);CHKERRQ(ierr); 451*412e9a14SMatthew G. Knepley ierr = PetscSectionGetStorageSize(cs, &csize);CHKERRQ(ierr); 452*412e9a14SMatthew G. Knepley for (c = 0; c < csize; ++c) cones[c] = -1; 453*412e9a14SMatthew G. Knepley } 454*412e9a14SMatthew G. Knepley /* Set cones */ 455*412e9a14SMatthew G. Knepley for (d = 0; d <= depth; ++d) { 456*412e9a14SMatthew G. Knepley const PetscInt *cone; 457*412e9a14SMatthew G. Knepley PetscInt pStart, pEnd, p; 458*412e9a14SMatthew G. Knepley 459*412e9a14SMatthew G. Knepley if (d == cellDepth) continue; 460*412e9a14SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 461*412e9a14SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 4629a5b13a2SMatthew G Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 4639a5b13a2SMatthew G Knepley ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr); 4649a5b13a2SMatthew G Knepley ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr); 4659a5b13a2SMatthew G Knepley ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr); 4669f074e33SMatthew G Knepley } 4679a5b13a2SMatthew G Knepley } 468*412e9a14SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 469*412e9a14SMatthew G. Knepley const PetscInt *cone, *faceSizes, *faces; 470*412e9a14SMatthew G. Knepley const DMPolytopeType *faceTypes; 471ba2698f1SMatthew G. Knepley DMPolytopeType ct; 472*412e9a14SMatthew G. Knepley PetscInt numFaces, cf, foff = 0; 47399836e0eSStefano Zampini 474ba2698f1SMatthew G. Knepley ierr = DMPlexGetCellType(dm, c, &ct);CHKERRQ(ierr); 47599836e0eSStefano Zampini ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 476*412e9a14SMatthew G. Knepley ierr = DMPlexGetRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);CHKERRQ(ierr); 477*412e9a14SMatthew G. Knepley for (cf = 0; cf < numFaces; foff += faceSizes[cf], ++cf) { 478*412e9a14SMatthew G. Knepley DMPolytopeType faceType = faceTypes[cf]; 479*412e9a14SMatthew G. Knepley const PetscInt faceSize = faceSizes[cf]; 480*412e9a14SMatthew G. Knepley const PetscInt *face = &faces[foff]; 481*412e9a14SMatthew G. Knepley const PetscInt *fcone; 4829a5b13a2SMatthew G Knepley PetscHashIJKLKey key; 483e8f14785SLisandro Dalcin PetscHashIter iter; 484e8f14785SLisandro Dalcin PetscBool missing; 485*412e9a14SMatthew G. Knepley PetscInt f; 4869f074e33SMatthew G Knepley 487*412e9a14SMatthew G. Knepley if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not support faces of size %D > 4", faceSize); 488*412e9a14SMatthew G. Knepley key.i = face[0]; 489*412e9a14SMatthew G. Knepley key.j = faceSize > 1 ? face[1] : PETSC_MAX_INT; 490*412e9a14SMatthew G. Knepley key.k = faceSize > 2 ? face[2] : PETSC_MAX_INT; 491*412e9a14SMatthew G. Knepley key.l = faceSize > 3 ? face[3] : PETSC_MAX_INT; 49299836e0eSStefano Zampini ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 49399836e0eSStefano Zampini ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 49499836e0eSStefano Zampini ierr = PetscHashIJKLIterGet(faceTable, iter, &f);CHKERRQ(ierr); 49599836e0eSStefano Zampini ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 496*412e9a14SMatthew G. Knepley ierr = DMPlexGetCone(idm, f, &fcone);CHKERRQ(ierr); 497*412e9a14SMatthew G. Knepley if (fcone[0] < 0) {ierr = DMPlexSetCone(idm, f, face);CHKERRQ(ierr);} 498*412e9a14SMatthew G. Knepley /* TODO This should be unnecessary since we have autoamtic orientation */ 499*412e9a14SMatthew G. Knepley { 500a74221b0SStefano Zampini /* when matching hybrid faces in 3D, only few cases are possible. 501a74221b0SStefano Zampini Face traversal however can no longer follow the usual convention, this seems a serious issue to me */ 502a74221b0SStefano Zampini PetscInt tquad_map[4][4] = { {PETSC_MIN_INT, 0,PETSC_MIN_INT,PETSC_MIN_INT}, 503a74221b0SStefano Zampini { -1,PETSC_MIN_INT,PETSC_MIN_INT,PETSC_MIN_INT}, 504a74221b0SStefano Zampini {PETSC_MIN_INT,PETSC_MIN_INT,PETSC_MIN_INT, 1}, 505a74221b0SStefano Zampini {PETSC_MIN_INT,PETSC_MIN_INT, -2,PETSC_MIN_INT} }; 506*412e9a14SMatthew G. Knepley PetscInt i, i2, j; 507*412e9a14SMatthew G. Knepley const PetscInt *cone; 508*412e9a14SMatthew G. Knepley PetscInt coneSize, ornt; 509a74221b0SStefano Zampini 510*412e9a14SMatthew G. Knepley /* Orient face: Do not allow reverse orientation at the first vertex */ 511*412e9a14SMatthew G. Knepley ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 512*412e9a14SMatthew G. Knepley ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 513*412e9a14SMatthew 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); 514*412e9a14SMatthew G. Knepley /* - First find the initial vertex */ 515*412e9a14SMatthew G. Knepley for (i = 0; i < faceSize; ++i) if (face[0] == cone[i]) break; 516*412e9a14SMatthew G. Knepley /* If we want to compare tensor faces to regular faces, we have to flip them and take the else branch here */ 517*412e9a14SMatthew G. Knepley if (faceType == DM_POLYTOPE_SEG_PRISM_TENSOR) { 518a74221b0SStefano Zampini /* find the second vertex */ 519*412e9a14SMatthew G. Knepley for (i2 = 0; i2 < faceSize; ++i2) if (face[1] == cone[i2]) break; 520*412e9a14SMatthew G. Knepley switch (faceSize) { 521a74221b0SStefano Zampini case 2: 522a74221b0SStefano Zampini ornt = i ? -2 : 0; 523a74221b0SStefano Zampini break; 524a74221b0SStefano Zampini case 4: 525a74221b0SStefano Zampini ornt = tquad_map[i][i2]; 526a74221b0SStefano Zampini break; 527*412e9a14SMatthew G. Knepley default: SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unhandled face size %D for face %D in cell %D", faceSize, f, c); 528*412e9a14SMatthew G. Knepley } 529*412e9a14SMatthew G. Knepley } else { 530*412e9a14SMatthew G. Knepley /* Try forward comparison */ 531*412e9a14SMatthew G. Knepley for (j = 0; j < faceSize; ++j) if (face[j] != cone[(i+j)%faceSize]) break; 532*412e9a14SMatthew G. Knepley if (j == faceSize) { 533*412e9a14SMatthew G. Knepley if ((faceSize == 2) && (i == 1)) ornt = -2; 534*412e9a14SMatthew G. Knepley else ornt = i; 535*412e9a14SMatthew G. Knepley } else { 536*412e9a14SMatthew G. Knepley /* Try backward comparison */ 537*412e9a14SMatthew G. Knepley for (j = 0; j < faceSize; ++j) if (face[j] != cone[(i+faceSize-j)%faceSize]) break; 538*412e9a14SMatthew G. Knepley if (j == faceSize) { 539*412e9a14SMatthew G. Knepley if (i == 0) ornt = -faceSize; 540*412e9a14SMatthew G. Knepley else ornt = -i; 541*412e9a14SMatthew G. Knepley } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine orientation of face %D in cell %D", f, c); 542a74221b0SStefano Zampini } 543a74221b0SStefano Zampini } 54499836e0eSStefano Zampini ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 54599836e0eSStefano Zampini } 54699836e0eSStefano Zampini } 547*412e9a14SMatthew G. Knepley ierr = DMPlexRestoreRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);CHKERRQ(ierr); 54899836e0eSStefano Zampini } 5499a5b13a2SMatthew G Knepley ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 5509a5b13a2SMatthew G Knepley ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr); 5519a5b13a2SMatthew G Knepley ierr = DMPlexStratify(idm);CHKERRQ(ierr); 5529f074e33SMatthew G Knepley PetscFunctionReturn(0); 5539f074e33SMatthew G Knepley } 5549f074e33SMatthew G Knepley 555f80536cbSVaclav Hapla static PetscErrorCode SortRmineRremoteByRemote_Private(PetscSF sf, PetscInt *rmine1[], PetscInt *rremote1[]) 556f80536cbSVaclav Hapla { 557f80536cbSVaclav Hapla PetscInt nleaves; 558f80536cbSVaclav Hapla PetscInt nranks; 559a0d42dbcSVaclav Hapla const PetscMPIInt *ranks=NULL; 560a0d42dbcSVaclav Hapla const PetscInt *roffset=NULL, *rmine=NULL, *rremote=NULL; 561f80536cbSVaclav Hapla PetscInt n, o, r; 562f80536cbSVaclav Hapla PetscErrorCode ierr; 563f80536cbSVaclav Hapla 564f80536cbSVaclav Hapla PetscFunctionBegin; 565dec1416fSJunchao Zhang ierr = PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote);CHKERRQ(ierr); 566f80536cbSVaclav Hapla nleaves = roffset[nranks]; 567f80536cbSVaclav Hapla ierr = PetscMalloc2(nleaves, rmine1, nleaves, rremote1);CHKERRQ(ierr); 568f80536cbSVaclav Hapla for (r=0; r<nranks; r++) { 569f80536cbSVaclav Hapla /* simultaneously sort rank-wise portions of rmine & rremote by values in rremote 570f80536cbSVaclav Hapla - to unify order with the other side */ 571f80536cbSVaclav Hapla o = roffset[r]; 572f80536cbSVaclav Hapla n = roffset[r+1] - o; 573580bdb30SBarry Smith ierr = PetscArraycpy(&(*rmine1)[o], &rmine[o], n);CHKERRQ(ierr); 574580bdb30SBarry Smith ierr = PetscArraycpy(&(*rremote1)[o], &rremote[o], n);CHKERRQ(ierr); 575f80536cbSVaclav Hapla ierr = PetscSortIntWithArray(n, &(*rremote1)[o], &(*rmine1)[o]);CHKERRQ(ierr); 576f80536cbSVaclav Hapla } 577f80536cbSVaclav Hapla PetscFunctionReturn(0); 578f80536cbSVaclav Hapla } 579f80536cbSVaclav Hapla 58027d0e99aSVaclav Hapla PetscErrorCode DMPlexOrientInterface_Internal(DM dm) 5812e745bdaSMatthew G. Knepley { 58227d0e99aSVaclav Hapla /* Here we only compare first 2 points of the cone. Full cone size would lead to stronger self-checking. */ 58327d0e99aSVaclav Hapla PetscInt masterCone[2]; 584cae7fe92SVaclav Hapla PetscInt (*roots)[2], (*leaves)[2]; 5858a650c75SVaclav Hapla PetscMPIInt (*rootsRanks)[2], (*leavesRanks)[2]; 58627d0e99aSVaclav Hapla 58727d0e99aSVaclav Hapla PetscSF sf=NULL; 588a0d42dbcSVaclav Hapla const PetscInt *locals=NULL; 589a0d42dbcSVaclav Hapla const PetscSFNode *remotes=NULL; 5908a650c75SVaclav Hapla PetscInt nroots, nleaves, p, c; 591f80536cbSVaclav Hapla PetscInt nranks, n, o, r; 592a0d42dbcSVaclav Hapla const PetscMPIInt *ranks=NULL; 593a0d42dbcSVaclav Hapla const PetscInt *roffset=NULL; 594a0d42dbcSVaclav Hapla PetscInt *rmine1=NULL, *rremote1=NULL; /* rmine and rremote copies simultaneously sorted by rank and rremote */ 595a0d42dbcSVaclav Hapla const PetscInt *cone=NULL; 596adeface4SVaclav Hapla PetscInt coneSize, ind0; 5972e745bdaSMatthew G. Knepley MPI_Comm comm; 59827d0e99aSVaclav Hapla PetscMPIInt rank,size; 5992e745bdaSMatthew G. Knepley PetscInt debug = 0; 6002e745bdaSMatthew G. Knepley PetscErrorCode ierr; 6012e745bdaSMatthew G. Knepley 6022e745bdaSMatthew G. Knepley PetscFunctionBegin; 603df6a9fadSVaclav Hapla ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 6042e745bdaSMatthew G. Knepley ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &locals, &remotes);CHKERRQ(ierr); 6053ede9f65SMatthew G. Knepley if (nroots < 0) PetscFunctionReturn(0); 606f80536cbSVaclav Hapla ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 607dec1416fSJunchao Zhang ierr = PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, NULL, NULL);CHKERRQ(ierr); 608dc21a0bfSVaclav Hapla ierr = DMViewFromOptions(dm, NULL, "-before_fix_dm_view");CHKERRQ(ierr); 60927d0e99aSVaclav Hapla #if defined(PETSC_USE_DEBUG) 610dc21a0bfSVaclav Hapla ierr = DMPlexCheckPointSF(dm);CHKERRQ(ierr); 611dc21a0bfSVaclav Hapla #endif 612f80536cbSVaclav Hapla ierr = SortRmineRremoteByRemote_Private(sf, &rmine1, &rremote1);CHKERRQ(ierr); 6138a650c75SVaclav Hapla ierr = PetscMalloc4(nroots, &roots, nroots, &leaves, nroots, &rootsRanks, nroots, &leavesRanks);CHKERRQ(ierr); 6142e745bdaSMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 6152e745bdaSMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 61627d0e99aSVaclav Hapla ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 6179e24d8a0SVaclav Hapla for (p = 0; p < nroots; ++p) { 6189e24d8a0SVaclav Hapla ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 6199e24d8a0SVaclav Hapla ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 62027d0e99aSVaclav Hapla if (coneSize < 2) { 62127d0e99aSVaclav Hapla for (c = 0; c < 2; c++) { 62227d0e99aSVaclav Hapla roots[p][c] = -1; 62327d0e99aSVaclav Hapla rootsRanks[p][c] = -1; 62427d0e99aSVaclav Hapla } 62527d0e99aSVaclav Hapla continue; 62627d0e99aSVaclav Hapla } 6272e745bdaSMatthew G. Knepley /* Translate all points to root numbering */ 6288a650c75SVaclav Hapla for (c = 0; c < 2; c++) { 6298a650c75SVaclav Hapla ierr = PetscFindInt(cone[c], nleaves, locals, &ind0);CHKERRQ(ierr); 6308a650c75SVaclav Hapla if (ind0 < 0) { 6318a650c75SVaclav Hapla roots[p][c] = cone[c]; 6328a650c75SVaclav Hapla rootsRanks[p][c] = rank; 633c8148b97SVaclav Hapla } else { 6348a650c75SVaclav Hapla roots[p][c] = remotes[ind0].index; 6358a650c75SVaclav Hapla rootsRanks[p][c] = remotes[ind0].rank; 6368a650c75SVaclav Hapla } 6372e745bdaSMatthew G. Knepley } 6382e745bdaSMatthew G. Knepley } 6399e24d8a0SVaclav Hapla for (p = 0; p < nroots; ++p) { 6408ccb905dSVaclav Hapla for (c = 0; c < 2; c++) { 6418ccb905dSVaclav Hapla leaves[p][c] = -2; 6428ccb905dSVaclav Hapla leavesRanks[p][c] = -2; 6438ccb905dSVaclav Hapla } 644c8148b97SVaclav Hapla } 6452e745bdaSMatthew G. Knepley ierr = PetscSFBcastBegin(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr); 6468a650c75SVaclav Hapla ierr = PetscSFBcastBegin(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr); 6472e745bdaSMatthew G. Knepley ierr = PetscSFBcastEnd(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr); 6488a650c75SVaclav Hapla ierr = PetscSFBcastEnd(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr); 649c8148b97SVaclav Hapla if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);} 65027d0e99aSVaclav Hapla if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Referenced roots\n");CHKERRQ(ierr);} 6519e24d8a0SVaclav Hapla for (p = 0; p < nroots; ++p) { 6529e24d8a0SVaclav Hapla if (leaves[p][0] < 0) continue; 6539e24d8a0SVaclav Hapla ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 6549e24d8a0SVaclav Hapla ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 65527d0e99aSVaclav Hapla if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] %4D: cone=[%4D %4D] roots=[(%d,%4D) (%d,%4D)] leaves=[(%d,%4D) (%d,%4D)]", rank, p, cone[0], cone[1], rootsRanks[p][0], roots[p][0], rootsRanks[p][1], roots[p][1], leavesRanks[p][0], leaves[p][0], leavesRanks[p][1], leaves[p][1]);CHKERRQ(ierr);} 65682f5c0aeSVaclav Hapla if ((leaves[p][0] != roots[p][0]) || (leaves[p][1] != roots[p][1]) || (leavesRanks[p][0] != rootsRanks[p][0]) || (leavesRanks[p][1] != rootsRanks[p][1])) { 65727d0e99aSVaclav Hapla /* Translate these two leaves to my cone points; masterCone means desired order p's cone points */ 658f80536cbSVaclav Hapla for (c = 0; c < 2; c++) { 65927d0e99aSVaclav Hapla if (leavesRanks[p][c] == rank) { 66027d0e99aSVaclav Hapla /* A local leave is just taken as it is */ 66127d0e99aSVaclav Hapla masterCone[c] = leaves[p][c]; 66227d0e99aSVaclav Hapla continue; 66327d0e99aSVaclav Hapla } 664f80536cbSVaclav Hapla /* Find index of rank leavesRanks[p][c] among remote ranks */ 665f80536cbSVaclav Hapla /* No need for PetscMPIIntCast because these integers were originally cast from PetscMPIInt. */ 666f80536cbSVaclav Hapla ierr = PetscFindMPIInt((PetscMPIInt)leavesRanks[p][c], nranks, ranks, &r);CHKERRQ(ierr); 66727d0e99aSVaclav Hapla if (PetscUnlikely(r < 0)) SETERRQ7(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D cone[%D]=%D root (%d,%D) leave (%d,%D): leave rank not found among remote ranks",p,c,cone[c],rootsRanks[p][c],roots[p][c],leavesRanks[p][c],leaves[p][c]); 66827d0e99aSVaclav Hapla if (PetscUnlikely(ranks[r] < 0 || ranks[r] >= size)) SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_PLIB, "p=%D c=%D commsize=%d: ranks[%D] = %d makes no sense",p,c,size,r,ranks[r]); 669f80536cbSVaclav Hapla /* Find point leaves[p][c] among remote points aimed at rank leavesRanks[p][c] */ 670f80536cbSVaclav Hapla o = roffset[r]; 671f80536cbSVaclav Hapla n = roffset[r+1] - o; 672f80536cbSVaclav Hapla ierr = PetscFindInt(leaves[p][c], n, &rremote1[o], &ind0);CHKERRQ(ierr); 67327d0e99aSVaclav Hapla if (PetscUnlikely(ind0 < 0)) SETERRQ7(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D cone[%D]=%D root (%d,%D) leave (%d,%D): corresponding remote point not found - it seems there is missing connection in point SF!",p,c,cone[c],rootsRanks[p][c],roots[p][c],leavesRanks[p][c],leaves[p][c]); 674f80536cbSVaclav Hapla /* Get the corresponding local point */ 675f80536cbSVaclav Hapla masterCone[c] = rmine1[o+ind0];CHKERRQ(ierr); 676f80536cbSVaclav Hapla } 67727d0e99aSVaclav Hapla if (debug) {ierr = PetscSynchronizedPrintf(comm, " masterCone=[%4D %4D]\n", masterCone[0], masterCone[1]);CHKERRQ(ierr);} 67827d0e99aSVaclav Hapla /* Set the desired order of p's cone points and fix orientations accordingly */ 679f80536cbSVaclav Hapla /* Vaclav's note: Here we only compare first 2 points of the cone. Full cone size would lead to stronger self-checking. */ 680f80536cbSVaclav Hapla ierr = DMPlexOrientCell(dm, p, 2, masterCone);CHKERRQ(ierr); 68127d0e99aSVaclav Hapla } else if (debug) {ierr = PetscSynchronizedPrintf(comm, " ==\n");CHKERRQ(ierr);} 68223aaf07bSVaclav Hapla } 68327d0e99aSVaclav Hapla if (debug) { 68427d0e99aSVaclav Hapla ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 68527d0e99aSVaclav Hapla ierr = MPI_Barrier(comm);CHKERRQ(ierr); 6862e745bdaSMatthew G. Knepley } 687adeface4SVaclav Hapla ierr = DMViewFromOptions(dm, NULL, "-after_fix_dm_view");CHKERRQ(ierr); 6888a650c75SVaclav Hapla ierr = PetscFree4(roots, leaves, rootsRanks, leavesRanks);CHKERRQ(ierr); 6897c7bb832SVaclav Hapla ierr = PetscFree2(rmine1, rremote1);CHKERRQ(ierr); 6902e745bdaSMatthew G. Knepley PetscFunctionReturn(0); 6912e745bdaSMatthew G. Knepley } 6922e745bdaSMatthew G. Knepley 6932e72742cSMatthew 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[]) 6947bffde78SMichael Lange { 6952e72742cSMatthew G. Knepley PetscInt idx; 6962e72742cSMatthew G. Knepley PetscMPIInt rank; 6972e72742cSMatthew G. Knepley PetscBool flg; 6987bffde78SMichael Lange PetscErrorCode ierr; 6997bffde78SMichael Lange 7007bffde78SMichael Lange PetscFunctionBegin; 7012e72742cSMatthew G. Knepley ierr = PetscOptionsHasName(NULL, NULL, opt, &flg);CHKERRQ(ierr); 7022e72742cSMatthew G. Knepley if (!flg) PetscFunctionReturn(0); 7032e72742cSMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 7042e72742cSMatthew G. Knepley ierr = PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);CHKERRQ(ierr); 7052e72742cSMatthew 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);} 7062e72742cSMatthew G. Knepley ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 7072e72742cSMatthew G. Knepley PetscFunctionReturn(0); 7082e72742cSMatthew G. Knepley } 7092e72742cSMatthew G. Knepley 7102e72742cSMatthew G. Knepley static PetscErrorCode SFNodeArrayViewFromOptions(MPI_Comm comm, const char opt[], const char name[], const char idxname[], PetscInt n, const PetscSFNode a[]) 7112e72742cSMatthew G. Knepley { 7122e72742cSMatthew G. Knepley PetscInt idx; 7132e72742cSMatthew G. Knepley PetscMPIInt rank; 7142e72742cSMatthew G. Knepley PetscBool flg; 7152e72742cSMatthew G. Knepley PetscErrorCode ierr; 7162e72742cSMatthew G. Knepley 7172e72742cSMatthew G. Knepley PetscFunctionBegin; 7182e72742cSMatthew G. Knepley ierr = PetscOptionsHasName(NULL, NULL, opt, &flg);CHKERRQ(ierr); 7192e72742cSMatthew G. Knepley if (!flg) PetscFunctionReturn(0); 7202e72742cSMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 7212e72742cSMatthew G. Knepley ierr = PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);CHKERRQ(ierr); 7222e72742cSMatthew G. Knepley if (idxname) { 7232e72742cSMatthew 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);} 7242e72742cSMatthew G. Knepley } else { 7252e72742cSMatthew 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);} 7262e72742cSMatthew G. Knepley } 7272e72742cSMatthew G. Knepley ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 7282e72742cSMatthew G. Knepley PetscFunctionReturn(0); 7292e72742cSMatthew G. Knepley } 7302e72742cSMatthew G. Knepley 7313be36e83SMatthew G. Knepley static PetscErrorCode DMPlexMapToLocalPoint(DM dm, PetscHMapIJ remotehash, PetscSFNode remotePoint, PetscInt *localPoint) 7322e72742cSMatthew G. Knepley { 7333be36e83SMatthew G. Knepley PetscSF sf; 7343be36e83SMatthew G. Knepley const PetscInt *locals; 7353be36e83SMatthew G. Knepley PetscMPIInt rank; 7362e72742cSMatthew G. Knepley PetscErrorCode ierr; 7372e72742cSMatthew G. Knepley 7382e72742cSMatthew G. Knepley PetscFunctionBegin; 7393be36e83SMatthew G. Knepley ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 7403be36e83SMatthew G. Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 7413be36e83SMatthew G. Knepley ierr = PetscSFGetGraph(sf, NULL, NULL, &locals, NULL);CHKERRQ(ierr); 7422e72742cSMatthew G. Knepley if (remotePoint.rank == rank) { 7432e72742cSMatthew G. Knepley *localPoint = remotePoint.index; 7442e72742cSMatthew G. Knepley } else { 7452e72742cSMatthew G. Knepley PetscHashIJKey key; 7463be36e83SMatthew G. Knepley PetscInt l; 7472e72742cSMatthew G. Knepley 7482e72742cSMatthew G. Knepley key.i = remotePoint.index; 7492e72742cSMatthew G. Knepley key.j = remotePoint.rank; 7503be36e83SMatthew G. Knepley ierr = PetscHMapIJGet(remotehash, key, &l);CHKERRQ(ierr); 7513be36e83SMatthew G. Knepley if (l >= 0) { 7523be36e83SMatthew G. Knepley *localPoint = locals[l]; 7532e72742cSMatthew G. Knepley } else PetscFunctionReturn(1); 7542e72742cSMatthew G. Knepley } 7552e72742cSMatthew G. Knepley PetscFunctionReturn(0); 7562e72742cSMatthew G. Knepley } 7572e72742cSMatthew G. Knepley 7583be36e83SMatthew G. Knepley static PetscErrorCode DMPlexMapToGlobalPoint(DM dm, PetscInt localPoint, PetscSFNode *remotePoint) 7593be36e83SMatthew G. Knepley { 7603be36e83SMatthew G. Knepley PetscSF sf; 7613be36e83SMatthew G. Knepley const PetscInt *locals, *rootdegree; 7623be36e83SMatthew G. Knepley const PetscSFNode *remotes; 7633be36e83SMatthew G. Knepley PetscInt Nl, l; 7643be36e83SMatthew G. Knepley PetscMPIInt rank; 7653be36e83SMatthew G. Knepley PetscErrorCode ierr; 7663be36e83SMatthew G. Knepley 7673be36e83SMatthew G. Knepley PetscFunctionBegin; 7683be36e83SMatthew G. Knepley ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 7693be36e83SMatthew G. Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 7703be36e83SMatthew G. Knepley ierr = PetscSFGetGraph(sf, NULL, &Nl, &locals, &remotes);CHKERRQ(ierr); 7713be36e83SMatthew G. Knepley if (Nl < 0) goto owned; 7723be36e83SMatthew G. Knepley ierr = PetscSFComputeDegreeBegin(sf, &rootdegree);CHKERRQ(ierr); 7733be36e83SMatthew G. Knepley ierr = PetscSFComputeDegreeEnd(sf, &rootdegree);CHKERRQ(ierr); 7743be36e83SMatthew G. Knepley if (rootdegree[localPoint]) goto owned; 7753be36e83SMatthew G. Knepley ierr = PetscFindInt(localPoint, Nl, locals, &l);CHKERRQ(ierr); 7763be36e83SMatthew G. Knepley if (l < 0) PetscFunctionReturn(1); 7773be36e83SMatthew G. Knepley *remotePoint = remotes[l]; 7783be36e83SMatthew G. Knepley PetscFunctionReturn(0); 7793be36e83SMatthew G. Knepley owned: 7803be36e83SMatthew G. Knepley remotePoint->rank = rank; 7813be36e83SMatthew G. Knepley remotePoint->index = localPoint; 7823be36e83SMatthew G. Knepley PetscFunctionReturn(0); 7833be36e83SMatthew G. Knepley } 7843be36e83SMatthew G. Knepley 7853be36e83SMatthew G. Knepley 7863be36e83SMatthew G. Knepley static PetscErrorCode DMPlexPointIsShared(DM dm, PetscInt p, PetscBool *isShared) 7873be36e83SMatthew G. Knepley { 7883be36e83SMatthew G. Knepley PetscSF sf; 7893be36e83SMatthew G. Knepley const PetscInt *locals, *rootdegree; 7903be36e83SMatthew G. Knepley PetscInt Nl, idx; 7913be36e83SMatthew G. Knepley PetscErrorCode ierr; 7923be36e83SMatthew G. Knepley 7933be36e83SMatthew G. Knepley PetscFunctionBegin; 7943be36e83SMatthew G. Knepley *isShared = PETSC_FALSE; 7953be36e83SMatthew G. Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 7963be36e83SMatthew G. Knepley ierr = PetscSFGetGraph(sf, NULL, &Nl, &locals, NULL);CHKERRQ(ierr); 7973be36e83SMatthew G. Knepley if (Nl < 0) PetscFunctionReturn(0); 7983be36e83SMatthew G. Knepley ierr = PetscFindInt(p, Nl, locals, &idx);CHKERRQ(ierr); 7993be36e83SMatthew G. Knepley if (idx >= 0) {*isShared = PETSC_TRUE; PetscFunctionReturn(0);} 8003be36e83SMatthew G. Knepley ierr = PetscSFComputeDegreeBegin(sf, &rootdegree);CHKERRQ(ierr); 8013be36e83SMatthew G. Knepley ierr = PetscSFComputeDegreeEnd(sf, &rootdegree);CHKERRQ(ierr); 8023be36e83SMatthew G. Knepley if (rootdegree[p] > 0) *isShared = PETSC_TRUE; 8033be36e83SMatthew G. Knepley PetscFunctionReturn(0); 8043be36e83SMatthew G. Knepley } 8053be36e83SMatthew G. Knepley 8063be36e83SMatthew G. Knepley static PetscErrorCode DMPlexConeIsShared(DM dm, PetscInt p, PetscBool *isShared) 8073be36e83SMatthew G. Knepley { 8083be36e83SMatthew G. Knepley const PetscInt *cone; 8093be36e83SMatthew G. Knepley PetscInt coneSize, c; 8103be36e83SMatthew G. Knepley PetscBool cShared = PETSC_TRUE; 8113be36e83SMatthew G. Knepley PetscErrorCode ierr; 8123be36e83SMatthew G. Knepley 8133be36e83SMatthew G. Knepley PetscFunctionBegin; 8143be36e83SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 8153be36e83SMatthew G. Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 8163be36e83SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 8173be36e83SMatthew G. Knepley PetscBool pointShared; 8183be36e83SMatthew G. Knepley 8193be36e83SMatthew G. Knepley ierr = DMPlexPointIsShared(dm, cone[c], &pointShared);CHKERRQ(ierr); 8203be36e83SMatthew G. Knepley cShared = (PetscBool) (cShared && pointShared); 8213be36e83SMatthew G. Knepley } 8223be36e83SMatthew G. Knepley *isShared = coneSize ? cShared : PETSC_FALSE; 8233be36e83SMatthew G. Knepley PetscFunctionReturn(0); 8243be36e83SMatthew G. Knepley } 8253be36e83SMatthew G. Knepley 8263be36e83SMatthew G. Knepley static PetscErrorCode DMPlexGetConeMinimum(DM dm, PetscInt p, PetscSFNode *cpmin) 8273be36e83SMatthew G. Knepley { 8283be36e83SMatthew G. Knepley const PetscInt *cone; 8293be36e83SMatthew G. Knepley PetscInt coneSize, c; 8303be36e83SMatthew G. Knepley PetscSFNode cmin = {PETSC_MAX_INT, PETSC_MAX_INT}, missing = {-1, -1}; 8313be36e83SMatthew G. Knepley PetscErrorCode ierr; 8323be36e83SMatthew G. Knepley 8333be36e83SMatthew G. Knepley PetscFunctionBegin; 8343be36e83SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 8353be36e83SMatthew G. Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 8363be36e83SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 8373be36e83SMatthew G. Knepley PetscSFNode rcp; 8383be36e83SMatthew G. Knepley 8393be36e83SMatthew G. Knepley ierr = DMPlexMapToGlobalPoint(dm, cone[c], &rcp); 8403be36e83SMatthew G. Knepley if (ierr) { 8413be36e83SMatthew G. Knepley cmin = missing; 8423be36e83SMatthew G. Knepley } else { 8433be36e83SMatthew G. Knepley cmin = (rcp.rank < cmin.rank) || (rcp.rank == cmin.rank && rcp.index < cmin.index) ? rcp : cmin; 8443be36e83SMatthew G. Knepley } 8453be36e83SMatthew G. Knepley } 8463be36e83SMatthew G. Knepley *cpmin = coneSize ? cmin : missing; 8473be36e83SMatthew G. Knepley PetscFunctionReturn(0); 8483be36e83SMatthew G. Knepley } 8493be36e83SMatthew G. Knepley 8503be36e83SMatthew G. Knepley /* 8513be36e83SMatthew G. Knepley Each shared face has an entry in the candidates array: 8523be36e83SMatthew G. Knepley (-1, coneSize-1), {(global cone point)} 8533be36e83SMatthew G. Knepley where the set is missing the point p which we use as the key for the face 8543be36e83SMatthew G. Knepley */ 8553be36e83SMatthew G. Knepley static PetscErrorCode DMPlexAddSharedFace_Private(DM dm, PetscSection candidateSection, PetscSFNode candidates[], PetscHMapIJ faceHash, PetscInt p, PetscBool debug) 8563be36e83SMatthew G. Knepley { 8573be36e83SMatthew G. Knepley MPI_Comm comm; 8583be36e83SMatthew G. Knepley const PetscInt *support; 859cf4dc471SVaclav Hapla PetscInt supportSize, s, off = 0, idx = 0, overlap, cellHeight, height; 8603be36e83SMatthew G. Knepley PetscMPIInt rank; 8613be36e83SMatthew G. Knepley PetscErrorCode ierr; 8623be36e83SMatthew G. Knepley 8633be36e83SMatthew G. Knepley PetscFunctionBegin; 8643be36e83SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 8653be36e83SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 866cf4dc471SVaclav Hapla ierr = DMPlexGetOverlap(dm, &overlap);CHKERRQ(ierr); 867cf4dc471SVaclav Hapla ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 868cf4dc471SVaclav Hapla ierr = DMPlexGetPointHeight(dm, p, &height);CHKERRQ(ierr); 869cf4dc471SVaclav Hapla if (!overlap && height <= cellHeight+1) { 870cf4dc471SVaclav Hapla /* cells can't be shared for non-overlapping meshes */ 871cf4dc471SVaclav Hapla if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] Skipping face %D to avoid adding cell to hashmap since this is nonoverlapping mesh\n", rank, p);CHKERRQ(ierr);} 872cf4dc471SVaclav Hapla PetscFunctionReturn(0); 873cf4dc471SVaclav Hapla } 8743be36e83SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr); 8753be36e83SMatthew G. Knepley ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr); 8763be36e83SMatthew G. Knepley if (candidates) {ierr = PetscSectionGetOffset(candidateSection, p, &off);CHKERRQ(ierr);} 8773be36e83SMatthew G. Knepley for (s = 0; s < supportSize; ++s) { 8783be36e83SMatthew G. Knepley const PetscInt face = support[s]; 8793be36e83SMatthew G. Knepley const PetscInt *cone; 8803be36e83SMatthew G. Knepley PetscSFNode cpmin={-1,-1}, rp={-1,-1}; 8813be36e83SMatthew G. Knepley PetscInt coneSize, c, f; 8823be36e83SMatthew G. Knepley PetscBool isShared = PETSC_FALSE; 8833be36e83SMatthew G. Knepley PetscHashIJKey key; 8843be36e83SMatthew G. Knepley 8853be36e83SMatthew G. Knepley /* Only add point once */ 8863be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] Support face %D\n", rank, face);CHKERRQ(ierr);} 8873be36e83SMatthew G. Knepley key.i = p; 8883be36e83SMatthew G. Knepley key.j = face; 8893be36e83SMatthew G. Knepley ierr = PetscHMapIJGet(faceHash, key, &f);CHKERRQ(ierr); 8903be36e83SMatthew G. Knepley if (f >= 0) continue; 8913be36e83SMatthew G. Knepley ierr = DMPlexConeIsShared(dm, face, &isShared);CHKERRQ(ierr); 8923be36e83SMatthew G. Knepley ierr = DMPlexGetConeMinimum(dm, face, &cpmin);CHKERRQ(ierr); 8933be36e83SMatthew G. Knepley ierr = DMPlexMapToGlobalPoint(dm, p, &rp);CHKERRQ(ierr); 8943be36e83SMatthew G. Knepley if (debug) { 8953be36e83SMatthew G. Knepley ierr = PetscSynchronizedPrintf(comm, "[%d] Face point %D is shared: %d\n", rank, face, (int) isShared);CHKERRQ(ierr); 8963be36e83SMatthew G. Knepley ierr = PetscSynchronizedPrintf(comm, "[%d] Global point (%D, %D) Min Cone Point (%D, %D)\n", rank, rp.rank, rp.index, cpmin.rank, cpmin.index);CHKERRQ(ierr); 8973be36e83SMatthew G. Knepley } 8983be36e83SMatthew G. Knepley if (isShared && (rp.rank == cpmin.rank && rp.index == cpmin.index)) { 8993be36e83SMatthew G. Knepley ierr = PetscHMapIJSet(faceHash, key, p);CHKERRQ(ierr); 9003be36e83SMatthew G. Knepley if (candidates) { 9013be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] Adding shared face %D at idx %D\n[%d] ", rank, face, idx, rank);CHKERRQ(ierr);} 9023be36e83SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, face, &coneSize);CHKERRQ(ierr); 9033be36e83SMatthew G. Knepley ierr = DMPlexGetCone(dm, face, &cone);CHKERRQ(ierr); 9043be36e83SMatthew G. Knepley candidates[off+idx].rank = -1; 9053be36e83SMatthew G. Knepley candidates[off+idx++].index = coneSize-1; 9063be36e83SMatthew G. Knepley candidates[off+idx].rank = rank; 9073be36e83SMatthew G. Knepley candidates[off+idx++].index = face; 9083be36e83SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 9093be36e83SMatthew G. Knepley const PetscInt cp = cone[c]; 9103be36e83SMatthew G. Knepley 9113be36e83SMatthew G. Knepley if (cp == p) continue; 9123be36e83SMatthew G. Knepley ierr = DMPlexMapToGlobalPoint(dm, cp, &candidates[off+idx]);CHKERRQ(ierr); 9133be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, " (%D,%D)", candidates[off+idx].rank, candidates[off+idx].index);CHKERRQ(ierr);} 9143be36e83SMatthew G. Knepley ++idx; 9153be36e83SMatthew G. Knepley } 9163be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr);} 9173be36e83SMatthew G. Knepley } else { 9183be36e83SMatthew G. Knepley /* Add cone size to section */ 9193be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] Scheduling shared face %D\n", rank, face);CHKERRQ(ierr);} 9203be36e83SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, face, &coneSize);CHKERRQ(ierr); 9213be36e83SMatthew G. Knepley ierr = PetscHMapIJSet(faceHash, key, p);CHKERRQ(ierr); 9223be36e83SMatthew G. Knepley ierr = PetscSectionAddDof(candidateSection, p, coneSize+1);CHKERRQ(ierr); 9233be36e83SMatthew G. Knepley } 9243be36e83SMatthew G. Knepley } 9253be36e83SMatthew G. Knepley } 9263be36e83SMatthew G. Knepley PetscFunctionReturn(0); 9273be36e83SMatthew G. Knepley } 9283be36e83SMatthew G. Knepley 9292e72742cSMatthew G. Knepley /*@ 9302e72742cSMatthew G. Knepley DMPlexInterpolatePointSF - Insert interpolated points in the overlap into the PointSF in parallel, following local interpolation 9312e72742cSMatthew G. Knepley 932d083f849SBarry Smith Collective on dm 9332e72742cSMatthew G. Knepley 9342e72742cSMatthew G. Knepley Input Parameters: 9352e72742cSMatthew G. Knepley + dm - The interpolated DM 9362e72742cSMatthew G. Knepley - pointSF - The initial SF without interpolated points 9372e72742cSMatthew G. Knepley 9382e72742cSMatthew G. Knepley Output Parameter: 9392e72742cSMatthew G. Knepley . pointSF - The SF including interpolated points 9402e72742cSMatthew G. Knepley 941f0cfc026SVaclav Hapla Level: developer 9422e72742cSMatthew G. Knepley 9432e72742cSMatthew 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 9442e72742cSMatthew G. Knepley 9452e72742cSMatthew G. Knepley .seealso: DMPlexInterpolate(), DMPlexUninterpolate() 9462e72742cSMatthew G. Knepley @*/ 947e53487d3SMatthew G. Knepley PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF) 9482e72742cSMatthew G. Knepley { 9493be36e83SMatthew G. Knepley MPI_Comm comm; 9503be36e83SMatthew G. Knepley PetscHMapIJ remoteHash; 9513be36e83SMatthew G. Knepley PetscHMapI claimshash; 9523be36e83SMatthew G. Knepley PetscSection candidateSection, candidateRemoteSection, claimSection; 9533be36e83SMatthew G. Knepley PetscSFNode *candidates, *candidatesRemote, *claims; 9542e72742cSMatthew G. Knepley const PetscInt *localPoints, *rootdegree; 9552e72742cSMatthew G. Knepley const PetscSFNode *remotePoints; 9563be36e83SMatthew G. Knepley PetscInt ov, Nr, r, Nl, l; 9573be36e83SMatthew G. Knepley PetscInt candidatesSize, candidatesRemoteSize, claimsSize; 9583be36e83SMatthew G. Knepley PetscBool flg, debug = PETSC_FALSE; 959f0cfc026SVaclav Hapla PetscMPIInt rank; 9602e72742cSMatthew G. Knepley PetscErrorCode ierr; 9612e72742cSMatthew G. Knepley 9622e72742cSMatthew G. Knepley PetscFunctionBegin; 963f0cfc026SVaclav Hapla PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 9643be36e83SMatthew G. Knepley PetscValidHeaderSpecific(pointSF, PETSCSF_CLASSID, 3); 965f0cfc026SVaclav Hapla ierr = DMPlexIsDistributed(dm, &flg);CHKERRQ(ierr); 966f0cfc026SVaclav Hapla if (!flg) PetscFunctionReturn(0); 9673be36e83SMatthew G. Knepley /* Set initial SF so that lower level queries work */ 9683be36e83SMatthew G. Knepley ierr = DMSetPointSF(dm, pointSF);CHKERRQ(ierr); 9693be36e83SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 9703be36e83SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 9713be36e83SMatthew G. Knepley ierr = DMPlexGetOverlap(dm, &ov);CHKERRQ(ierr); 9723be36e83SMatthew G. Knepley if (ov) SETERRQ(comm, PETSC_ERR_SUP, "Interpolation of overlapped DMPlex not implemented yet"); 9733be36e83SMatthew G. Knepley ierr = PetscOptionsHasName(NULL, ((PetscObject) dm)->prefix, "-dmplex_interp_debug", &debug);CHKERRQ(ierr); 9742e72742cSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_interp_pre_view");CHKERRQ(ierr); 9752e72742cSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) pointSF, NULL, "-petscsf_interp_pre_view");CHKERRQ(ierr); 97625afeb17SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 9773be36e83SMatthew G. Knepley /* Step 0: Precalculations */ 9783be36e83SMatthew G. Knepley ierr = PetscSFGetGraph(pointSF, &Nr, &Nl, &localPoints, &remotePoints);CHKERRQ(ierr); 9793be36e83SMatthew G. Knepley if (Nr < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "This DMPlex is distributed but input PointSF has no graph set"); 9803be36e83SMatthew G. Knepley ierr = PetscHMapIJCreate(&remoteHash);CHKERRQ(ierr); 9813be36e83SMatthew G. Knepley for (l = 0; l < Nl; ++l) { 9823be36e83SMatthew G. Knepley PetscHashIJKey key; 9832e72742cSMatthew G. Knepley key.i = remotePoints[l].index; 9842e72742cSMatthew G. Knepley key.j = remotePoints[l].rank; 9853be36e83SMatthew G. Knepley ierr = PetscHMapIJSet(remoteHash, key, l);CHKERRQ(ierr); 9867bffde78SMichael Lange } 98766aa2a39SMatthew G. Knepley /* Compute root degree to identify shared points */ 9882e72742cSMatthew G. Knepley ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr); 9892e72742cSMatthew G. Knepley ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr); 9903be36e83SMatthew G. Knepley ierr = IntArrayViewFromOptions(comm, "-interp_root_degree_view", "Root degree", "point", "degree", Nr, rootdegree);CHKERRQ(ierr); 9913be36e83SMatthew G. Knepley /* 9923be36e83SMatthew G. Knepley 1) Loop over each leaf point $p$ at depth $d$ in the SF 9933be36e83SMatthew G. Knepley \item Get set $F(p)$ of faces $f$ in the support of $p$ for which 9943be36e83SMatthew G. Knepley \begin{itemize} 9953be36e83SMatthew G. Knepley \item all cone points of $f$ are shared 9963be36e83SMatthew G. Knepley \item $p$ is the cone point with smallest canonical number 9973be36e83SMatthew G. Knepley \end{itemize} 9983be36e83SMatthew G. Knepley \item Send $F(p)$ and the cone of each face to the active root point $r(p)$ 9993be36e83SMatthew G. Knepley \item At the root, if at least two faces with a given cone are present, including a local face, mark the face as shared \label{alg:rootStep} and choose the root face 10003be36e83SMatthew G. Knepley \item Send the root face from the root back to all leaf process 10013be36e83SMatthew G. Knepley \item Leaf processes add the shared face to the SF 10023be36e83SMatthew G. Knepley */ 10033be36e83SMatthew G. Knepley /* Step 1: Construct section+SFNode array 10043be36e83SMatthew G. Knepley The section has entries for all shared faces for which we have a leaf point in the cone 10053be36e83SMatthew G. Knepley The array holds candidate shared faces, each face is refered to by the leaf point */ 10063be36e83SMatthew G. Knepley ierr = PetscSectionCreate(comm, &candidateSection);CHKERRQ(ierr); 10073be36e83SMatthew G. Knepley ierr = PetscSectionSetChart(candidateSection, 0, Nr);CHKERRQ(ierr); 10087bffde78SMichael Lange { 10093be36e83SMatthew G. Knepley PetscHMapIJ faceHash; 10102e72742cSMatthew G. Knepley 10113be36e83SMatthew G. Knepley ierr = PetscHMapIJCreate(&faceHash);CHKERRQ(ierr); 10123be36e83SMatthew G. Knepley for (l = 0; l < Nl; ++l) { 10133be36e83SMatthew G. Knepley const PetscInt p = localPoints[l]; 10142e72742cSMatthew G. Knepley 10153be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] First pass leaf point %D\n", rank, p);CHKERRQ(ierr);} 10163be36e83SMatthew G. Knepley ierr = DMPlexAddSharedFace_Private(dm, candidateSection, NULL, faceHash, p, debug);CHKERRQ(ierr); 10172e72742cSMatthew G. Knepley } 10183be36e83SMatthew G. Knepley ierr = PetscHMapIJClear(faceHash);CHKERRQ(ierr); 10197bffde78SMichael Lange ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr); 10207bffde78SMichael Lange ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr); 10217bffde78SMichael Lange ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr); 10223be36e83SMatthew G. Knepley for (l = 0; l < Nl; ++l) { 10233be36e83SMatthew G. Knepley const PetscInt p = localPoints[l]; 10242e72742cSMatthew G. Knepley 10253be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] Second pass leaf point %D\n", rank, p);CHKERRQ(ierr);} 10263be36e83SMatthew G. Knepley ierr = DMPlexAddSharedFace_Private(dm, candidateSection, candidates, faceHash, p, debug);CHKERRQ(ierr); 10272e72742cSMatthew G. Knepley } 10283be36e83SMatthew G. Knepley ierr = PetscHMapIJDestroy(&faceHash);CHKERRQ(ierr); 10293be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);} 10307bffde78SMichael Lange } 10313be36e83SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) candidateSection, "Candidate Section");CHKERRQ(ierr); 10322e72742cSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) candidateSection, NULL, "-petscsection_interp_candidate_view");CHKERRQ(ierr); 10333be36e83SMatthew G. Knepley ierr = SFNodeArrayViewFromOptions(comm, "-petscsection_interp_candidate_view", "Candidates", NULL, candidatesSize, candidates);CHKERRQ(ierr); 10343be36e83SMatthew G. Knepley /* Step 2: Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */ 10352e72742cSMatthew G. Knepley /* Note that this section is indexed by offsets into leaves, not by point number */ 10367bffde78SMichael Lange { 10377bffde78SMichael Lange PetscSF sfMulti, sfInverse, sfCandidates; 10387bffde78SMichael Lange PetscInt *remoteOffsets; 10392e72742cSMatthew G. Knepley 10407bffde78SMichael Lange ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 10417bffde78SMichael Lange ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr); 10423be36e83SMatthew G. Knepley ierr = PetscSectionCreate(comm, &candidateRemoteSection);CHKERRQ(ierr); 10433be36e83SMatthew G. Knepley ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateRemoteSection);CHKERRQ(ierr); 10443be36e83SMatthew G. Knepley ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateRemoteSection, &sfCandidates);CHKERRQ(ierr); 10453be36e83SMatthew G. Knepley ierr = PetscSectionGetStorageSize(candidateRemoteSection, &candidatesRemoteSize);CHKERRQ(ierr); 10467bffde78SMichael Lange ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr); 10477bffde78SMichael Lange ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 10487bffde78SMichael Lange ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 10497bffde78SMichael Lange ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr); 10507bffde78SMichael Lange ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr); 10517bffde78SMichael Lange ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 10522e72742cSMatthew G. Knepley 10533be36e83SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) candidateRemoteSection, "Remote Candidate Section");CHKERRQ(ierr); 10543be36e83SMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) candidateRemoteSection, NULL, "-petscsection_interp_candidate_remote_view");CHKERRQ(ierr); 10553be36e83SMatthew G. Knepley ierr = SFNodeArrayViewFromOptions(comm, "-petscsection_interp_candidate_remote_view", "Remote Candidates", NULL, candidatesRemoteSize, candidatesRemote);CHKERRQ(ierr); 10567bffde78SMichael Lange } 10573be36e83SMatthew G. Knepley /* Step 3: At the root, if at least two faces with a given cone are present, including a local face, mark the face as shared and choose the root face */ 10587bffde78SMichael Lange { 10593be36e83SMatthew G. Knepley PetscHashIJKLRemote faceTable; 10603be36e83SMatthew G. Knepley PetscInt idx, idx2; 10613be36e83SMatthew G. Knepley 10623be36e83SMatthew G. Knepley ierr = PetscHashIJKLRemoteCreate(&faceTable);CHKERRQ(ierr); 10632e72742cSMatthew G. Knepley /* There is a section point for every leaf attached to a given root point */ 10643be36e83SMatthew G. Knepley for (r = 0, idx = 0, idx2 = 0; r < Nr; ++r) { 10652e72742cSMatthew G. Knepley PetscInt deg; 10663be36e83SMatthew G. Knepley 10672e72742cSMatthew G. Knepley for (deg = 0; deg < rootdegree[r]; ++deg, ++idx) { 10682e72742cSMatthew G. Knepley PetscInt offset, dof, d; 10692e72742cSMatthew G. Knepley 10703be36e83SMatthew G. Knepley ierr = PetscSectionGetDof(candidateRemoteSection, idx, &dof);CHKERRQ(ierr); 10713be36e83SMatthew G. Knepley ierr = PetscSectionGetOffset(candidateRemoteSection, idx, &offset);CHKERRQ(ierr); 10723be36e83SMatthew G. Knepley /* dof may include many faces from the remote process */ 10732e72742cSMatthew G. Knepley for (d = 0; d < dof; ++d) { 10743be36e83SMatthew G. Knepley const PetscInt hidx = offset+d; 10753be36e83SMatthew G. Knepley const PetscInt Np = candidatesRemote[hidx].index+1; 10763be36e83SMatthew G. Knepley const PetscSFNode rface = candidatesRemote[hidx+1]; 10773be36e83SMatthew G. Knepley const PetscSFNode *fcone = &candidatesRemote[hidx+2]; 10783be36e83SMatthew G. Knepley PetscSFNode fcp0; 10793be36e83SMatthew G. Knepley const PetscSFNode pmax = {PETSC_MAX_INT, PETSC_MAX_INT}; 10802e72742cSMatthew G. Knepley const PetscInt *join = NULL; 10813be36e83SMatthew G. Knepley PetscHashIJKLRemoteKey key; 10823be36e83SMatthew G. Knepley PetscHashIter iter; 10833be36e83SMatthew G. Knepley PetscBool missing; 10842e72742cSMatthew G. Knepley PetscInt points[1024], p, joinSize; 10852e72742cSMatthew G. Knepley 108666e92ce5SVaclav Hapla if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Checking face (%D, %D) at (%D, %D, %D) with cone size %D\n", rank, rface.rank, rface.index, r, idx, d, Np);CHKERRQ(ierr);} 108766e92ce5SVaclav Hapla if (Np > 4) SETERRQ6(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle face (%D, %D) at (%D, %D, %D) with %D cone points", rface.rank, rface.index, r, idx, d, Np); 10883be36e83SMatthew G. Knepley fcp0.rank = rank; 10893be36e83SMatthew G. Knepley fcp0.index = r; 10903be36e83SMatthew G. Knepley d += Np; 10913be36e83SMatthew G. Knepley /* Put remote face in hash table */ 10923be36e83SMatthew G. Knepley key.i = fcp0; 10933be36e83SMatthew G. Knepley key.j = fcone[0]; 10943be36e83SMatthew G. Knepley key.k = Np > 2 ? fcone[1] : pmax; 10953be36e83SMatthew G. Knepley key.l = Np > 3 ? fcone[2] : pmax; 10963be36e83SMatthew G. Knepley ierr = PetscSortSFNode(Np, (PetscSFNode *) &key);CHKERRQ(ierr); 10973be36e83SMatthew G. Knepley ierr = PetscHashIJKLRemotePut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 10983be36e83SMatthew G. Knepley if (missing) { 10993be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Setting remote face (%D, %D)\n", rank, rface.index, rface.rank);CHKERRQ(ierr);} 11003be36e83SMatthew G. Knepley ierr = PetscHashIJKLRemoteIterSet(faceTable, iter, rface);CHKERRQ(ierr); 11013be36e83SMatthew G. Knepley } else { 11023be36e83SMatthew G. Knepley PetscSFNode oface; 11033be36e83SMatthew G. Knepley 11043be36e83SMatthew G. Knepley ierr = PetscHashIJKLRemoteIterGet(faceTable, iter, &oface);CHKERRQ(ierr); 11053be36e83SMatthew G. Knepley if ((rface.rank < oface.rank) || (rface.rank == oface.rank && rface.index < oface.index)) { 11063be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Replacing with remote face (%D, %D)\n", rank, rface.index, rface.rank);CHKERRQ(ierr);} 11073be36e83SMatthew G. Knepley ierr = PetscHashIJKLRemoteIterSet(faceTable, iter, rface);CHKERRQ(ierr); 11083be36e83SMatthew G. Knepley } 11093be36e83SMatthew G. Knepley } 11103be36e83SMatthew G. Knepley /* Check for local face */ 11112e72742cSMatthew G. Knepley points[0] = r; 11123be36e83SMatthew G. Knepley for (p = 1; p < Np; ++p) { 11133be36e83SMatthew G. Knepley ierr = DMPlexMapToLocalPoint(dm, remoteHash, fcone[p-1], &points[p]); 11143be36e83SMatthew G. Knepley if (ierr) break; /* We got a point not in our overlap */ 11153be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Checking local candidate %D\n", rank, points[p]);CHKERRQ(ierr);} 11167bffde78SMichael Lange } 11172e72742cSMatthew G. Knepley if (ierr) continue; 11183be36e83SMatthew G. Knepley ierr = DMPlexGetJoin(dm, Np, points, &joinSize, &join);CHKERRQ(ierr); 11197bffde78SMichael Lange if (joinSize == 1) { 11203be36e83SMatthew G. Knepley PetscSFNode lface; 11213be36e83SMatthew G. Knepley PetscSFNode oface; 11223be36e83SMatthew G. Knepley 11233be36e83SMatthew G. Knepley /* Always replace with local face */ 11243be36e83SMatthew G. Knepley lface.rank = rank; 11253be36e83SMatthew G. Knepley lface.index = join[0]; 11263be36e83SMatthew G. Knepley ierr = PetscHashIJKLRemoteIterGet(faceTable, iter, &oface);CHKERRQ(ierr); 11273be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Replacing (%D, %D) with local face (%D, %D)\n", rank, oface.index, oface.rank, lface.index, lface.rank);CHKERRQ(ierr);} 11283be36e83SMatthew G. Knepley ierr = PetscHashIJKLRemoteIterSet(faceTable, iter, lface);CHKERRQ(ierr); 11297bffde78SMichael Lange } 11303be36e83SMatthew G. Knepley ierr = DMPlexRestoreJoin(dm, Np, points, &joinSize, &join);CHKERRQ(ierr); 11313be36e83SMatthew G. Knepley } 11323be36e83SMatthew G. Knepley } 11333be36e83SMatthew G. Knepley /* Put back faces for this root */ 11343be36e83SMatthew G. Knepley for (deg = 0; deg < rootdegree[r]; ++deg, ++idx2) { 11353be36e83SMatthew G. Knepley PetscInt offset, dof, d; 11363be36e83SMatthew G. Knepley 11373be36e83SMatthew G. Knepley ierr = PetscSectionGetDof(candidateRemoteSection, idx2, &dof);CHKERRQ(ierr); 11383be36e83SMatthew G. Knepley ierr = PetscSectionGetOffset(candidateRemoteSection, idx2, &offset);CHKERRQ(ierr); 11393be36e83SMatthew G. Knepley /* dof may include many faces from the remote process */ 11403be36e83SMatthew G. Knepley for (d = 0; d < dof; ++d) { 11413be36e83SMatthew G. Knepley const PetscInt hidx = offset+d; 11423be36e83SMatthew G. Knepley const PetscInt Np = candidatesRemote[hidx].index+1; 11433be36e83SMatthew G. Knepley const PetscSFNode *fcone = &candidatesRemote[hidx+2]; 11443be36e83SMatthew G. Knepley PetscSFNode fcp0; 11453be36e83SMatthew G. Knepley const PetscSFNode pmax = {PETSC_MAX_INT, PETSC_MAX_INT}; 11463be36e83SMatthew G. Knepley PetscHashIJKLRemoteKey key; 11473be36e83SMatthew G. Knepley PetscHashIter iter; 11483be36e83SMatthew G. Knepley PetscBool missing; 11493be36e83SMatthew G. Knepley 11503be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Entering face at (%D, %D)\n", rank, r, idx);CHKERRQ(ierr);} 11513be36e83SMatthew G. Knepley if (Np > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle faces with %D cone points", Np); 11523be36e83SMatthew G. Knepley fcp0.rank = rank; 11533be36e83SMatthew G. Knepley fcp0.index = r; 11543be36e83SMatthew G. Knepley d += Np; 11553be36e83SMatthew G. Knepley /* Find remote face in hash table */ 11563be36e83SMatthew G. Knepley key.i = fcp0; 11573be36e83SMatthew G. Knepley key.j = fcone[0]; 11583be36e83SMatthew G. Knepley key.k = Np > 2 ? fcone[1] : pmax; 11593be36e83SMatthew G. Knepley key.l = Np > 3 ? fcone[2] : pmax; 11603be36e83SMatthew G. Knepley ierr = PetscSortSFNode(Np, (PetscSFNode *) &key);CHKERRQ(ierr); 11613be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] key (%D, %D) (%D, %D) (%D, %D) (%D, %D)\n", rank, key.i.rank, key.i.index, key.j.rank, key.j.index, key.k.rank, key.k.index, key.l.rank, key.l.index);CHKERRQ(ierr);} 11623be36e83SMatthew G. Knepley ierr = PetscHashIJKLRemotePut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 11633be36e83SMatthew G. Knepley if (missing) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Root %D Idx %D ought to have an assoicated face", r, idx2); 11643be36e83SMatthew G. Knepley else {ierr = PetscHashIJKLRemoteIterGet(faceTable, iter, &candidatesRemote[hidx]);CHKERRQ(ierr);} 11657bffde78SMichael Lange } 11667bffde78SMichael Lange } 11677bffde78SMichael Lange } 11682e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);} 11693be36e83SMatthew G. Knepley ierr = PetscHashIJKLRemoteDestroy(&faceTable);CHKERRQ(ierr); 11707bffde78SMichael Lange } 11713be36e83SMatthew G. Knepley /* Step 4: Push back owned faces */ 11727bffde78SMichael Lange { 11737bffde78SMichael Lange PetscSF sfMulti, sfClaims, sfPointNew; 11747bffde78SMichael Lange PetscSFNode *remotePointsNew; 11752e72742cSMatthew G. Knepley PetscInt *remoteOffsets, *localPointsNew; 11763be36e83SMatthew G. Knepley PetscInt pStart, pEnd, r, NlNew, p; 11772e72742cSMatthew G. Knepley 11783be36e83SMatthew G. Knepley /* 4) Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */ 11797bffde78SMichael Lange ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 11803be36e83SMatthew G. Knepley ierr = PetscSectionCreate(comm, &claimSection);CHKERRQ(ierr); 11813be36e83SMatthew G. Knepley ierr = PetscSFDistributeSection(sfMulti, candidateRemoteSection, &remoteOffsets, claimSection);CHKERRQ(ierr); 11823be36e83SMatthew G. Knepley ierr = PetscSFCreateSectionSF(sfMulti, candidateRemoteSection, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr); 11832e72742cSMatthew G. Knepley ierr = PetscSectionGetStorageSize(claimSection, &claimsSize);CHKERRQ(ierr); 11842e72742cSMatthew G. Knepley ierr = PetscMalloc1(claimsSize, &claims);CHKERRQ(ierr); 11853be36e83SMatthew G. Knepley for (p = 0; p < claimsSize; ++p) claims[p].rank = -1; 11867bffde78SMichael Lange ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 11877bffde78SMichael Lange ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 11887bffde78SMichael Lange ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr); 11897bffde78SMichael Lange ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 11903be36e83SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) claimSection, "Claim Section");CHKERRQ(ierr); 11912e72742cSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) claimSection, NULL, "-petscsection_interp_claim_view");CHKERRQ(ierr); 11923be36e83SMatthew G. Knepley ierr = SFNodeArrayViewFromOptions(comm, "-petscsection_interp_claim_view", "Claims", NULL, claimsSize, claims);CHKERRQ(ierr); 11933be36e83SMatthew G. Knepley /* Step 5) Walk the original section of local supports and add an SF entry for each updated item */ 11943be36e83SMatthew G. Knepley /* TODO I should not have to do a join here since I already put the face and its cone in the candidate section */ 1195e8f14785SLisandro Dalcin ierr = PetscHMapICreate(&claimshash);CHKERRQ(ierr); 11963be36e83SMatthew G. Knepley for (r = 0; r < Nr; ++r) { 11973be36e83SMatthew G. Knepley PetscInt dof, off, d; 11982e72742cSMatthew G. Knepley 11993be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] Checking root for claims %D\n", rank, r);CHKERRQ(ierr);} 12003be36e83SMatthew G. Knepley ierr = PetscSectionGetDof(candidateSection, r, &dof);CHKERRQ(ierr); 12013be36e83SMatthew G. Knepley ierr = PetscSectionGetOffset(candidateSection, r, &off);CHKERRQ(ierr); 12022e72742cSMatthew G. Knepley for (d = 0; d < dof;) { 12033be36e83SMatthew G. Knepley if (claims[off+d].rank >= 0) { 12043be36e83SMatthew G. Knepley const PetscInt faceInd = off+d; 12053be36e83SMatthew G. Knepley const PetscInt Np = candidates[off+d].index; 12062e72742cSMatthew G. Knepley const PetscInt *join = NULL; 12072e72742cSMatthew G. Knepley PetscInt joinSize, points[1024], c; 12082e72742cSMatthew G. Knepley 12093be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] Found claim for remote point (%D, %D)\n", rank, claims[faceInd].rank, claims[faceInd].index);CHKERRQ(ierr);} 12103be36e83SMatthew G. Knepley points[0] = r; 12113be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] point %D\n", rank, points[0]);CHKERRQ(ierr);} 12123be36e83SMatthew G. Knepley for (c = 0, d += 2; c < Np; ++c, ++d) { 12133be36e83SMatthew G. Knepley ierr = DMPlexMapToLocalPoint(dm, remoteHash, candidates[off+d], &points[c+1]);CHKERRQ(ierr); 12143be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] point %D\n", rank, points[c+1]);CHKERRQ(ierr);} 12152e72742cSMatthew G. Knepley } 12163be36e83SMatthew G. Knepley ierr = DMPlexGetJoin(dm, Np+1, points, &joinSize, &join);CHKERRQ(ierr); 12172e72742cSMatthew G. Knepley if (joinSize == 1) { 12183be36e83SMatthew G. Knepley if (claims[faceInd].rank == rank) { 12193be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] Ignoring local face %D for non-remote partner\n", rank, join[0]);CHKERRQ(ierr);} 12203be36e83SMatthew G. Knepley } else { 12213be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] Found local face %D\n", rank, join[0]);CHKERRQ(ierr);} 12222e72742cSMatthew G. Knepley ierr = PetscHMapISet(claimshash, join[0], faceInd);CHKERRQ(ierr); 12232e72742cSMatthew G. Knepley } 12243be36e83SMatthew G. Knepley } else { 12253be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] Failed to find face\n", rank);CHKERRQ(ierr);} 12263be36e83SMatthew G. Knepley } 12273be36e83SMatthew G. Knepley ierr = DMPlexRestoreJoin(dm, Np+1, points, &joinSize, &join);CHKERRQ(ierr); 12283be36e83SMatthew G. Knepley } else { 12293be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] No claim for point %D\n", rank, r);CHKERRQ(ierr);} 12303be36e83SMatthew G. Knepley d += claims[off+d].index+1; 12317bffde78SMichael Lange } 12327bffde78SMichael Lange } 12333be36e83SMatthew G. Knepley } 12343be36e83SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);} 12353be36e83SMatthew G. Knepley /* Step 6) Create new pointSF from hashed claims */ 12363be36e83SMatthew G. Knepley ierr = PetscHMapIGetSize(claimshash, &NlNew);CHKERRQ(ierr); 12377bffde78SMichael Lange ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 12383be36e83SMatthew G. Knepley ierr = PetscMalloc1(Nl + NlNew, &localPointsNew);CHKERRQ(ierr); 12393be36e83SMatthew G. Knepley ierr = PetscMalloc1(Nl + NlNew, &remotePointsNew);CHKERRQ(ierr); 12403be36e83SMatthew G. Knepley for (l = 0; l < Nl; ++l) { 12413be36e83SMatthew G. Knepley localPointsNew[l] = localPoints[l]; 12423be36e83SMatthew G. Knepley remotePointsNew[l].index = remotePoints[l].index; 12433be36e83SMatthew G. Knepley remotePointsNew[l].rank = remotePoints[l].rank; 12447bffde78SMichael Lange } 12453be36e83SMatthew G. Knepley p = Nl; 1246e8f14785SLisandro Dalcin ierr = PetscHMapIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr); 12473be36e83SMatthew G. Knepley /* We sort new points, and assume they are numbered after all existing points */ 12483be36e83SMatthew G. Knepley ierr = PetscSortInt(NlNew, &localPointsNew[Nl]);CHKERRQ(ierr); 12493be36e83SMatthew G. Knepley for (p = Nl; p < Nl + NlNew; ++p) { 12503be36e83SMatthew G. Knepley PetscInt off; 12513be36e83SMatthew G. Knepley ierr = PetscHMapIGet(claimshash, localPointsNew[p], &off);CHKERRQ(ierr); 12523be36e83SMatthew G. Knepley if (claims[off].rank < 0 || claims[off].index < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid claim for local point %D, (%D, %D)", localPointsNew[p], claims[off].rank, claims[off].index); 12533be36e83SMatthew G. Knepley remotePointsNew[p] = claims[off]; 12547bffde78SMichael Lange } 12553be36e83SMatthew G. Knepley ierr = PetscSFCreate(comm, &sfPointNew);CHKERRQ(ierr); 12563be36e83SMatthew G. Knepley ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, Nl+NlNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); 12573be36e83SMatthew G. Knepley ierr = PetscSFSetUp(sfPointNew);CHKERRQ(ierr); 12587bffde78SMichael Lange ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr); 12593be36e83SMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) sfPointNew, NULL, "-petscsf_interp_view");CHKERRQ(ierr); 12607bffde78SMichael Lange ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr); 1261e8f14785SLisandro Dalcin ierr = PetscHMapIDestroy(&claimshash);CHKERRQ(ierr); 12627bffde78SMichael Lange } 12633be36e83SMatthew G. Knepley ierr = PetscHMapIJDestroy(&remoteHash);CHKERRQ(ierr); 12647bffde78SMichael Lange ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr); 12653be36e83SMatthew G. Knepley ierr = PetscSectionDestroy(&candidateRemoteSection);CHKERRQ(ierr); 12667bffde78SMichael Lange ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr); 12677bffde78SMichael Lange ierr = PetscFree(candidates);CHKERRQ(ierr); 12687bffde78SMichael Lange ierr = PetscFree(candidatesRemote);CHKERRQ(ierr); 12697bffde78SMichael Lange ierr = PetscFree(claims);CHKERRQ(ierr); 127025afeb17SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 12717bffde78SMichael Lange PetscFunctionReturn(0); 12727bffde78SMichael Lange } 12737bffde78SMichael Lange 1274248eb259SVaclav Hapla /*@ 127580330477SMatthew G. Knepley DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc. 127680330477SMatthew G. Knepley 1277d083f849SBarry Smith Collective on dm 127880330477SMatthew G. Knepley 1279e56d480eSMatthew G. Knepley Input Parameters: 1280e56d480eSMatthew G. Knepley + dm - The DMPlex object with only cells and vertices 128110a67516SMatthew G. Knepley - dmInt - The interpolated DM 128280330477SMatthew G. Knepley 128380330477SMatthew G. Knepley Output Parameter: 12844e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object 128580330477SMatthew G. Knepley 128680330477SMatthew G. Knepley Level: intermediate 128780330477SMatthew G. Knepley 128895452b02SPatrick Sanan Notes: 128995452b02SPatrick Sanan It does not copy over the coordinates. 129043eeeb2dSStefano Zampini 12919ade3219SVaclav Hapla Developer Notes: 12929ade3219SVaclav Hapla It sets plex->interpolated = DMPLEX_INTERPOLATED_FULL. 12939ade3219SVaclav Hapla 129443eeeb2dSStefano Zampini .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates() 129580330477SMatthew G. Knepley @*/ 12969f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt) 12979f074e33SMatthew G Knepley { 1298a7748839SVaclav Hapla DMPlexInterpolatedFlag interpolated; 12999a5b13a2SMatthew G Knepley DM idm, odm = dm; 13007bffde78SMichael Lange PetscSF sfPoint; 13012e1b13c2SMatthew G. Knepley PetscInt depth, dim, d; 130210a67516SMatthew G. Knepley const char *name; 1303b325530aSVaclav Hapla PetscBool flg=PETSC_TRUE; 13049f074e33SMatthew G Knepley PetscErrorCode ierr; 13059f074e33SMatthew G Knepley 13069f074e33SMatthew G Knepley PetscFunctionBegin; 130710a67516SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 130810a67516SMatthew G. Knepley PetscValidPointer(dmInt, 2); 1309a72f3261SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 13102e1b13c2SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1311c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1312827c4036SVaclav Hapla ierr = DMPlexIsInterpolated(dm, &interpolated);CHKERRQ(ierr); 1313827c4036SVaclav Hapla if (interpolated == DMPLEX_INTERPOLATED_PARTIAL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not for partially interpolated meshes"); 1314827c4036SVaclav Hapla if (interpolated == DMPLEX_INTERPOLATED_FULL) { 131576b791aaSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 131676b791aaSMatthew G. Knepley idm = dm; 1317b21b8912SMatthew G. Knepley } else { 13189a5b13a2SMatthew G Knepley for (d = 1; d < dim; ++d) { 13199a5b13a2SMatthew G Knepley /* Create interpolated mesh */ 132010a67516SMatthew G. Knepley ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr); 13219a5b13a2SMatthew G Knepley ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr); 1322c73cfb54SMatthew G. Knepley ierr = DMSetDimension(idm, dim);CHKERRQ(ierr); 13237bffde78SMichael Lange if (depth > 0) { 13247bffde78SMichael Lange ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr); 13257bffde78SMichael Lange ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr); 132694488200SVaclav Hapla { 13273be36e83SMatthew G. Knepley /* TODO: We need to systematically fix cases of distributed Plexes with no graph set */ 132894488200SVaclav Hapla PetscInt nroots; 132994488200SVaclav Hapla ierr = PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 133094488200SVaclav Hapla if (nroots >= 0) {ierr = DMPlexInterpolatePointSF(idm, sfPoint);CHKERRQ(ierr);} 133194488200SVaclav Hapla } 13327bffde78SMichael Lange } 13339a5b13a2SMatthew G Knepley if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 13349a5b13a2SMatthew G Knepley odm = idm; 13359f074e33SMatthew G Knepley } 133610a67516SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr); 133710a67516SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) idm, name);CHKERRQ(ierr); 133810a67516SMatthew G. Knepley ierr = DMPlexCopyCoordinates(dm, idm);CHKERRQ(ierr); 13395d80c0bfSVaclav Hapla ierr = DMCopyLabels(dm, idm, PETSC_COPY_VALUES, PETSC_FALSE);CHKERRQ(ierr); 1340b325530aSVaclav Hapla ierr = PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_interpolate_orient_interfaces", &flg, NULL);CHKERRQ(ierr); 134127d0e99aSVaclav Hapla if (flg) {ierr = DMPlexOrientInterface_Internal(idm);CHKERRQ(ierr);} 134284699f70SSatish Balay } 134343eeeb2dSStefano Zampini { 134443eeeb2dSStefano Zampini PetscBool isper; 134543eeeb2dSStefano Zampini const PetscReal *maxCell, *L; 134643eeeb2dSStefano Zampini const DMBoundaryType *bd; 134743eeeb2dSStefano Zampini 134843eeeb2dSStefano Zampini ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr); 134943eeeb2dSStefano Zampini ierr = DMSetPeriodicity(idm,isper,maxCell,L,bd);CHKERRQ(ierr); 135043eeeb2dSStefano Zampini } 1351827c4036SVaclav Hapla /* This function makes the mesh fully interpolated on all ranks */ 1352827c4036SVaclav Hapla { 1353d6e9e4f7SVaclav Hapla DM_Plex *plex = (DM_Plex *) idm->data; 1354827c4036SVaclav Hapla plex->interpolated = plex->interpolatedCollective = DMPLEX_INTERPOLATED_FULL; 1355827c4036SVaclav Hapla } 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 1364d083f849SBarry Smith Collective on dmA 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 137665f90189SJed Brown .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection() 137780330477SMatthew G. Knepley @*/ 137807b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB) 137907b0a1fcSMatthew G Knepley { 138007b0a1fcSMatthew G Knepley Vec coordinatesA, coordinatesB; 138143eeeb2dSStefano Zampini VecType vtype; 138207b0a1fcSMatthew G Knepley PetscSection coordSectionA, coordSectionB; 138307b0a1fcSMatthew G Knepley PetscScalar *coordsA, *coordsB; 13840bedd6beSMatthew G. Knepley PetscInt spaceDim, Nf, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d; 138590a8aa44SJed Brown PetscInt cStartA, cEndA, cStartB, cEndB, cS, cE, cdim; 138643eeeb2dSStefano Zampini PetscBool lc = PETSC_FALSE; 138707b0a1fcSMatthew G Knepley PetscErrorCode ierr; 138807b0a1fcSMatthew G Knepley 138907b0a1fcSMatthew G Knepley PetscFunctionBegin; 139043eeeb2dSStefano Zampini PetscValidHeaderSpecific(dmA, DM_CLASSID, 1); 139143eeeb2dSStefano Zampini PetscValidHeaderSpecific(dmB, DM_CLASSID, 2); 139276b791aaSMatthew G. Knepley if (dmA == dmB) PetscFunctionReturn(0); 139390a8aa44SJed Brown ierr = DMGetCoordinateDim(dmA, &cdim);CHKERRQ(ierr); 139490a8aa44SJed Brown ierr = DMSetCoordinateDim(dmB, cdim);CHKERRQ(ierr); 139507b0a1fcSMatthew G Knepley ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr); 139607b0a1fcSMatthew G Knepley ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr); 139707b0a1fcSMatthew 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); 139843eeeb2dSStefano Zampini ierr = DMPlexGetHeightStratum(dmA, 0, &cStartA, &cEndA);CHKERRQ(ierr); 139943eeeb2dSStefano Zampini ierr = DMPlexGetHeightStratum(dmB, 0, &cStartB, &cEndB);CHKERRQ(ierr); 140069d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr); 140169d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr); 1402972bc18aSToby Isaac if (coordSectionA == coordSectionB) PetscFunctionReturn(0); 14030bedd6beSMatthew G. Knepley ierr = PetscSectionGetNumFields(coordSectionA, &Nf);CHKERRQ(ierr); 14040bedd6beSMatthew G. Knepley if (!Nf) PetscFunctionReturn(0); 14050bedd6beSMatthew G. Knepley if (Nf > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate fields must be 1, not %D", Nf); 1406df26b574SMatthew G. Knepley if (!coordSectionB) { 1407df26b574SMatthew G. Knepley PetscInt dim; 1408df26b574SMatthew G. Knepley 1409df26b574SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr); 1410df26b574SMatthew G. Knepley ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr); 1411df26b574SMatthew G. Knepley ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr); 1412df26b574SMatthew G. Knepley ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr); 1413df26b574SMatthew G. Knepley } 141407b0a1fcSMatthew G Knepley ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr); 141507b0a1fcSMatthew G Knepley ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr); 141607b0a1fcSMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr); 141743eeeb2dSStefano Zampini ierr = PetscSectionGetChart(coordSectionA, &cS, &cE);CHKERRQ(ierr); 141843eeeb2dSStefano Zampini if (cStartA <= cS && cS < cEndA) { /* localized coordinates */ 1419367003a6SStefano 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); 142043eeeb2dSStefano Zampini cS = cS - cStartA + cStartB; 142143eeeb2dSStefano Zampini cE = vEndB; 142243eeeb2dSStefano Zampini lc = PETSC_TRUE; 142343eeeb2dSStefano Zampini } else { 142443eeeb2dSStefano Zampini cS = vStartB; 142543eeeb2dSStefano Zampini cE = vEndB; 142643eeeb2dSStefano Zampini } 142743eeeb2dSStefano Zampini ierr = PetscSectionSetChart(coordSectionB, cS, cE);CHKERRQ(ierr); 142807b0a1fcSMatthew G Knepley for (v = vStartB; v < vEndB; ++v) { 142907b0a1fcSMatthew G Knepley ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr); 143007b0a1fcSMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr); 143107b0a1fcSMatthew G Knepley } 143243eeeb2dSStefano Zampini if (lc) { /* localized coordinates */ 143343eeeb2dSStefano Zampini PetscInt c; 143443eeeb2dSStefano Zampini 143543eeeb2dSStefano Zampini for (c = cS-cStartB; c < cEndB-cStartB; c++) { 143643eeeb2dSStefano Zampini PetscInt dof; 143743eeeb2dSStefano Zampini 143843eeeb2dSStefano Zampini ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr); 143943eeeb2dSStefano Zampini ierr = PetscSectionSetDof(coordSectionB, c + cStartB, dof);CHKERRQ(ierr); 144043eeeb2dSStefano Zampini ierr = PetscSectionSetFieldDof(coordSectionB, c + cStartB, 0, dof);CHKERRQ(ierr); 144143eeeb2dSStefano Zampini } 144243eeeb2dSStefano Zampini } 144307b0a1fcSMatthew G Knepley ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr); 144407b0a1fcSMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr); 144507b0a1fcSMatthew G Knepley ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr); 14468b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinatesB);CHKERRQ(ierr); 144707b0a1fcSMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr); 144807b0a1fcSMatthew G Knepley ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr); 144943eeeb2dSStefano Zampini ierr = VecGetBlockSize(coordinatesA, &d);CHKERRQ(ierr); 145043eeeb2dSStefano Zampini ierr = VecSetBlockSize(coordinatesB, d);CHKERRQ(ierr); 145143eeeb2dSStefano Zampini ierr = VecGetType(coordinatesA, &vtype);CHKERRQ(ierr); 145243eeeb2dSStefano Zampini ierr = VecSetType(coordinatesB, vtype);CHKERRQ(ierr); 145307b0a1fcSMatthew G Knepley ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr); 145407b0a1fcSMatthew G Knepley ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr); 145507b0a1fcSMatthew G Knepley for (v = 0; v < vEndB-vStartB; ++v) { 145643eeeb2dSStefano Zampini PetscInt offA, offB; 145743eeeb2dSStefano Zampini 145843eeeb2dSStefano Zampini ierr = PetscSectionGetOffset(coordSectionA, v + vStartA, &offA);CHKERRQ(ierr); 145943eeeb2dSStefano Zampini ierr = PetscSectionGetOffset(coordSectionB, v + vStartB, &offB);CHKERRQ(ierr); 146007b0a1fcSMatthew G Knepley for (d = 0; d < spaceDim; ++d) { 146143eeeb2dSStefano Zampini coordsB[offB+d] = coordsA[offA+d]; 146243eeeb2dSStefano Zampini } 146343eeeb2dSStefano Zampini } 146443eeeb2dSStefano Zampini if (lc) { /* localized coordinates */ 146543eeeb2dSStefano Zampini PetscInt c; 146643eeeb2dSStefano Zampini 146743eeeb2dSStefano Zampini for (c = cS-cStartB; c < cEndB-cStartB; c++) { 146843eeeb2dSStefano Zampini PetscInt dof, offA, offB; 146943eeeb2dSStefano Zampini 147043eeeb2dSStefano Zampini ierr = PetscSectionGetOffset(coordSectionA, c + cStartA, &offA);CHKERRQ(ierr); 147143eeeb2dSStefano Zampini ierr = PetscSectionGetOffset(coordSectionB, c + cStartB, &offB);CHKERRQ(ierr); 147243eeeb2dSStefano Zampini ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr); 1473580bdb30SBarry Smith ierr = PetscArraycpy(coordsB + offB,coordsA + offA,dof);CHKERRQ(ierr); 147407b0a1fcSMatthew G Knepley } 147507b0a1fcSMatthew G Knepley } 147607b0a1fcSMatthew G Knepley ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr); 147707b0a1fcSMatthew G Knepley ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr); 147807b0a1fcSMatthew G Knepley ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr); 147907b0a1fcSMatthew G Knepley ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr); 148007b0a1fcSMatthew G Knepley PetscFunctionReturn(0); 148107b0a1fcSMatthew G Knepley } 14825c386225SMatthew G. Knepley 14834e3744c5SMatthew G. Knepley /*@ 14844e3744c5SMatthew G. Knepley DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh 14854e3744c5SMatthew G. Knepley 1486d083f849SBarry Smith Collective on dm 14874e3744c5SMatthew G. Knepley 14884e3744c5SMatthew G. Knepley Input Parameter: 14894e3744c5SMatthew G. Knepley . dm - The complete DMPlex object 14904e3744c5SMatthew G. Knepley 14914e3744c5SMatthew G. Knepley Output Parameter: 14924e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices 14934e3744c5SMatthew G. Knepley 14944e3744c5SMatthew G. Knepley Level: intermediate 14954e3744c5SMatthew G. Knepley 149695452b02SPatrick Sanan Notes: 149795452b02SPatrick Sanan It does not copy over the coordinates. 149843eeeb2dSStefano Zampini 14999ade3219SVaclav Hapla Developer Notes: 15009ade3219SVaclav Hapla It sets plex->interpolated = DMPLEX_INTERPOLATED_NONE. 15019ade3219SVaclav Hapla 150243eeeb2dSStefano Zampini .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates() 15034e3744c5SMatthew G. Knepley @*/ 15044e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint) 15054e3744c5SMatthew G. Knepley { 1506827c4036SVaclav Hapla DMPlexInterpolatedFlag interpolated; 15074e3744c5SMatthew G. Knepley DM udm; 1508*412e9a14SMatthew G. Knepley PetscInt dim, vStart, vEnd, cStart, cEnd, c, maxConeSize = 0, *cone; 15094e3744c5SMatthew G. Knepley PetscErrorCode ierr; 15104e3744c5SMatthew G. Knepley 15114e3744c5SMatthew G. Knepley PetscFunctionBegin; 151243eeeb2dSStefano Zampini PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 151343eeeb2dSStefano Zampini PetscValidPointer(dmUnint, 2); 1514c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1515827c4036SVaclav Hapla ierr = DMPlexIsInterpolated(dm, &interpolated);CHKERRQ(ierr); 1516827c4036SVaclav Hapla if (interpolated == DMPLEX_INTERPOLATED_PARTIAL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not for partially interpolated meshes"); 1517827c4036SVaclav Hapla if (interpolated == DMPLEX_INTERPOLATED_NONE || dim <= 1) { 1518827c4036SVaclav Hapla /* in case dim <= 1 just keep the DMPLEX_INTERPOLATED_FULL flag */ 15194e3744c5SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 1520595d4abbSMatthew G. Knepley *dmUnint = dm; 1521595d4abbSMatthew G. Knepley PetscFunctionReturn(0); 15224e3744c5SMatthew G. Knepley } 1523595d4abbSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1524595d4abbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 15254e3744c5SMatthew G. Knepley ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr); 15264e3744c5SMatthew G. Knepley ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr); 1527c73cfb54SMatthew G. Knepley ierr = DMSetDimension(udm, dim);CHKERRQ(ierr); 15284e3744c5SMatthew G. Knepley ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr); 15294e3744c5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1530595d4abbSMatthew G. Knepley PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 15314e3744c5SMatthew G. Knepley 15324e3744c5SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 15334e3744c5SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 15344e3744c5SMatthew G. Knepley const PetscInt p = closure[cl]; 15354e3744c5SMatthew G. Knepley 15364e3744c5SMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) ++coneSize; 15374e3744c5SMatthew G. Knepley } 15384e3744c5SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 15394e3744c5SMatthew G. Knepley ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr); 1540595d4abbSMatthew G. Knepley maxConeSize = PetscMax(maxConeSize, coneSize); 15414e3744c5SMatthew G. Knepley } 15424e3744c5SMatthew G. Knepley ierr = DMSetUp(udm);CHKERRQ(ierr); 1543785e854fSJed Brown ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr); 15444e3744c5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1545595d4abbSMatthew G. Knepley PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 15464e3744c5SMatthew G. Knepley 15474e3744c5SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 15484e3744c5SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 15494e3744c5SMatthew G. Knepley const PetscInt p = closure[cl]; 15504e3744c5SMatthew G. Knepley 15514e3744c5SMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p; 15524e3744c5SMatthew G. Knepley } 15534e3744c5SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 15544e3744c5SMatthew G. Knepley ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr); 15554e3744c5SMatthew G. Knepley } 15564e3744c5SMatthew G. Knepley ierr = PetscFree(cone);CHKERRQ(ierr); 15574e3744c5SMatthew G. Knepley ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr); 15584e3744c5SMatthew G. Knepley ierr = DMPlexStratify(udm);CHKERRQ(ierr); 15595c7de58aSMatthew G. Knepley /* Reduce SF */ 15605c7de58aSMatthew G. Knepley { 15615c7de58aSMatthew G. Knepley PetscSF sfPoint, sfPointUn; 15625c7de58aSMatthew G. Knepley const PetscSFNode *remotePoints; 15635c7de58aSMatthew G. Knepley const PetscInt *localPoints; 15645c7de58aSMatthew G. Knepley PetscSFNode *remotePointsUn; 15655c7de58aSMatthew G. Knepley PetscInt *localPointsUn; 15665c7de58aSMatthew G. Knepley PetscInt vEnd, numRoots, numLeaves, l; 15675c7de58aSMatthew G. Knepley PetscInt numLeavesUn = 0, n = 0; 15685c7de58aSMatthew G. Knepley PetscErrorCode ierr; 15695c7de58aSMatthew G. Knepley 15705c7de58aSMatthew G. Knepley /* Get original SF information */ 15715c7de58aSMatthew G. Knepley ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 15725c7de58aSMatthew G. Knepley ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr); 15735c7de58aSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr); 15745c7de58aSMatthew G. Knepley ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 15755c7de58aSMatthew G. Knepley /* Allocate space for cells and vertices */ 15765c7de58aSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++; 15775c7de58aSMatthew G. Knepley /* Fill in leaves */ 15785c7de58aSMatthew G. Knepley if (vEnd >= 0) { 15795c7de58aSMatthew G. Knepley ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr); 15805c7de58aSMatthew G. Knepley ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr); 15815c7de58aSMatthew G. Knepley for (l = 0; l < numLeaves; l++) { 15825c7de58aSMatthew G. Knepley if (localPoints[l] < vEnd) { 15835c7de58aSMatthew G. Knepley localPointsUn[n] = localPoints[l]; 15845c7de58aSMatthew G. Knepley remotePointsUn[n].rank = remotePoints[l].rank; 15855c7de58aSMatthew G. Knepley remotePointsUn[n].index = remotePoints[l].index; 15865c7de58aSMatthew G. Knepley ++n; 15875c7de58aSMatthew G. Knepley } 15885c7de58aSMatthew G. Knepley } 15895c7de58aSMatthew G. Knepley if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn); 15905c7de58aSMatthew G. Knepley ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr); 15915c7de58aSMatthew G. Knepley } 15925c7de58aSMatthew G. Knepley } 159343eeeb2dSStefano Zampini { 159443eeeb2dSStefano Zampini PetscBool isper; 159543eeeb2dSStefano Zampini const PetscReal *maxCell, *L; 159643eeeb2dSStefano Zampini const DMBoundaryType *bd; 159743eeeb2dSStefano Zampini 159843eeeb2dSStefano Zampini ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr); 159943eeeb2dSStefano Zampini ierr = DMSetPeriodicity(udm,isper,maxCell,L,bd);CHKERRQ(ierr); 160043eeeb2dSStefano Zampini } 1601827c4036SVaclav Hapla /* This function makes the mesh fully uninterpolated on all ranks */ 1602827c4036SVaclav Hapla { 1603d6e9e4f7SVaclav Hapla DM_Plex *plex = (DM_Plex *) udm->data; 1604827c4036SVaclav Hapla plex->interpolated = plex->interpolatedCollective = DMPLEX_INTERPOLATED_NONE; 1605827c4036SVaclav Hapla } 16064e3744c5SMatthew G. Knepley *dmUnint = udm; 16074e3744c5SMatthew G. Knepley PetscFunctionReturn(0); 16084e3744c5SMatthew G. Knepley } 1609a7748839SVaclav Hapla 1610a7748839SVaclav Hapla static PetscErrorCode DMPlexIsInterpolated_Internal(DM dm, DMPlexInterpolatedFlag *interpolated) 1611a7748839SVaclav Hapla { 1612a7748839SVaclav Hapla PetscInt coneSize, depth, dim, h, p, pStart, pEnd; 1613a7748839SVaclav Hapla MPI_Comm comm; 1614a7748839SVaclav Hapla PetscErrorCode ierr; 1615a7748839SVaclav Hapla 1616a7748839SVaclav Hapla PetscFunctionBegin; 1617a7748839SVaclav Hapla ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 1618a7748839SVaclav Hapla ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1619a7748839SVaclav Hapla ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1620a7748839SVaclav Hapla 1621a7748839SVaclav Hapla if (depth == dim) { 1622a7748839SVaclav Hapla *interpolated = DMPLEX_INTERPOLATED_FULL; 1623a7748839SVaclav Hapla if (!dim) goto finish; 1624a7748839SVaclav Hapla 1625a7748839SVaclav Hapla /* Check points at height = dim are vertices (have no cones) */ 1626a7748839SVaclav Hapla ierr = DMPlexGetHeightStratum(dm, dim, &pStart, &pEnd);CHKERRQ(ierr); 1627a7748839SVaclav Hapla for (p=pStart; p<pEnd; p++) { 1628a7748839SVaclav Hapla ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 1629a7748839SVaclav Hapla if (coneSize) { 1630a7748839SVaclav Hapla *interpolated = DMPLEX_INTERPOLATED_PARTIAL; 1631a7748839SVaclav Hapla goto finish; 1632a7748839SVaclav Hapla } 1633a7748839SVaclav Hapla } 1634a7748839SVaclav Hapla 1635a7748839SVaclav Hapla /* Check points at height < dim have cones */ 1636a7748839SVaclav Hapla for (h=0; h<dim; h++) { 1637a7748839SVaclav Hapla ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr); 1638a7748839SVaclav Hapla for (p=pStart; p<pEnd; p++) { 1639a7748839SVaclav Hapla ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 1640a7748839SVaclav Hapla if (!coneSize) { 1641a7748839SVaclav Hapla *interpolated = DMPLEX_INTERPOLATED_PARTIAL; 1642a7748839SVaclav Hapla goto finish; 1643a7748839SVaclav Hapla } 1644a7748839SVaclav Hapla } 1645a7748839SVaclav Hapla } 1646a7748839SVaclav Hapla } else if (depth == 1) { 1647a7748839SVaclav Hapla *interpolated = DMPLEX_INTERPOLATED_NONE; 1648a7748839SVaclav Hapla } else { 1649a7748839SVaclav Hapla *interpolated = DMPLEX_INTERPOLATED_PARTIAL; 1650a7748839SVaclav Hapla } 1651a7748839SVaclav Hapla finish: 1652a7748839SVaclav Hapla PetscFunctionReturn(0); 1653a7748839SVaclav Hapla } 1654a7748839SVaclav Hapla 1655a7748839SVaclav Hapla /*@ 16569ade3219SVaclav Hapla DMPlexIsInterpolated - Find out to what extent the DMPlex is topologically interpolated. 1657a7748839SVaclav Hapla 1658a7748839SVaclav Hapla Not Collective 1659a7748839SVaclav Hapla 1660a7748839SVaclav Hapla Input Parameter: 1661a7748839SVaclav Hapla . dm - The DM object 1662a7748839SVaclav Hapla 1663a7748839SVaclav Hapla Output Parameter: 1664a7748839SVaclav Hapla . interpolated - Flag whether the DM is interpolated 1665a7748839SVaclav Hapla 1666a7748839SVaclav Hapla Level: intermediate 1667a7748839SVaclav Hapla 1668a7748839SVaclav Hapla Notes: 16699ade3219SVaclav Hapla Unlike DMPlexIsInterpolatedCollective(), this is NOT collective 16709ade3219SVaclav Hapla so the results can be different on different ranks in special cases. 1671a7748839SVaclav Hapla However, DMPlexInterpolate() guarantees the result is the same on all. 16729ade3219SVaclav Hapla 1673a7748839SVaclav Hapla Unlike DMPlexIsInterpolatedCollective(), this cannot return DMPLEX_INTERPOLATED_MIXED. 1674a7748839SVaclav Hapla 16759ade3219SVaclav Hapla Developer Notes: 16769ade3219SVaclav Hapla Initially, plex->interpolated = DMPLEX_INTERPOLATED_INVALID. 16779ade3219SVaclav Hapla 16789ade3219SVaclav Hapla If plex->interpolated == DMPLEX_INTERPOLATED_INVALID, DMPlexIsInterpolated_Internal() is called. 16799ade3219SVaclav Hapla It checks the actual topology and sets plex->interpolated on each rank separately to one of 16809ade3219SVaclav Hapla DMPLEX_INTERPOLATED_NONE, DMPLEX_INTERPOLATED_PARTIAL or DMPLEX_INTERPOLATED_FULL. 16819ade3219SVaclav Hapla 16829ade3219SVaclav Hapla If plex->interpolated != DMPLEX_INTERPOLATED_INVALID, this function just returns plex->interpolated. 16839ade3219SVaclav Hapla 16849ade3219SVaclav Hapla DMPlexInterpolate() sets plex->interpolated = DMPLEX_INTERPOLATED_FULL, 16859ade3219SVaclav Hapla and DMPlexUninterpolate() sets plex->interpolated = DMPLEX_INTERPOLATED_NONE. 16869ade3219SVaclav Hapla 1687a7748839SVaclav Hapla .seealso: DMPlexInterpolate(), DMPlexIsInterpolatedCollective() 1688a7748839SVaclav Hapla @*/ 1689a7748839SVaclav Hapla PetscErrorCode DMPlexIsInterpolated(DM dm, DMPlexInterpolatedFlag *interpolated) 1690a7748839SVaclav Hapla { 1691a7748839SVaclav Hapla DM_Plex *plex = (DM_Plex *) dm->data; 1692a7748839SVaclav Hapla PetscErrorCode ierr; 1693a7748839SVaclav Hapla 1694a7748839SVaclav Hapla PetscFunctionBegin; 1695a7748839SVaclav Hapla PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1696a7748839SVaclav Hapla PetscValidPointer(interpolated,2); 1697a7748839SVaclav Hapla if (plex->interpolated < 0) { 1698a7748839SVaclav Hapla ierr = DMPlexIsInterpolated_Internal(dm, &plex->interpolated);CHKERRQ(ierr); 1699a7748839SVaclav Hapla } else { 1700a7748839SVaclav Hapla #if defined (PETSC_USE_DEBUG) 1701a7748839SVaclav Hapla DMPlexInterpolatedFlag flg; 1702a7748839SVaclav Hapla 1703a7748839SVaclav Hapla ierr = DMPlexIsInterpolated_Internal(dm, &flg);CHKERRQ(ierr); 17047fc06600SVaclav Hapla if (flg != plex->interpolated) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Stashed DMPlexInterpolatedFlag %s is inconsistent with current %s", DMPlexInterpolatedFlags[plex->interpolated], DMPlexInterpolatedFlags[flg]); 1705a7748839SVaclav Hapla #endif 1706a7748839SVaclav Hapla } 1707a7748839SVaclav Hapla *interpolated = plex->interpolated; 1708a7748839SVaclav Hapla PetscFunctionReturn(0); 1709a7748839SVaclav Hapla } 1710a7748839SVaclav Hapla 1711a7748839SVaclav Hapla /*@ 17129ade3219SVaclav Hapla DMPlexIsInterpolatedCollective - Find out to what extent the DMPlex is topologically interpolated (in collective manner). 1713a7748839SVaclav Hapla 17142666ff3cSVaclav Hapla Collective 1715a7748839SVaclav Hapla 1716a7748839SVaclav Hapla Input Parameter: 1717a7748839SVaclav Hapla . dm - The DM object 1718a7748839SVaclav Hapla 1719a7748839SVaclav Hapla Output Parameter: 1720a7748839SVaclav Hapla . interpolated - Flag whether the DM is interpolated 1721a7748839SVaclav Hapla 1722a7748839SVaclav Hapla Level: intermediate 1723a7748839SVaclav Hapla 1724a7748839SVaclav Hapla Notes: 17259ade3219SVaclav Hapla Unlike DMPlexIsInterpolated(), this is collective so the results are guaranteed to be the same on all ranks. 17269ade3219SVaclav Hapla 17279ade3219SVaclav Hapla This function will return DMPLEX_INTERPOLATED_MIXED if the results of DMPlexIsInterpolated() are different on different ranks. 17289ade3219SVaclav Hapla 17299ade3219SVaclav Hapla Developer Notes: 17309ade3219SVaclav Hapla Initially, plex->interpolatedCollective = DMPLEX_INTERPOLATED_INVALID. 17319ade3219SVaclav Hapla 17329ade3219SVaclav Hapla If plex->interpolatedCollective == DMPLEX_INTERPOLATED_INVALID, this function calls DMPlexIsInterpolated() which sets plex->interpolated. 17339ade3219SVaclav Hapla MPI_Allreduce() is then called and collectively consistent flag plex->interpolatedCollective is set and returned; 17349ade3219SVaclav Hapla if plex->interpolated varies on different ranks, plex->interpolatedCollective = DMPLEX_INTERPOLATED_MIXED, 17359ade3219SVaclav Hapla otherwise sets plex->interpolatedCollective = plex->interpolated. 17369ade3219SVaclav Hapla 17379ade3219SVaclav Hapla If plex->interpolatedCollective != DMPLEX_INTERPOLATED_INVALID, this function just returns plex->interpolatedCollective. 1738a7748839SVaclav Hapla 1739a7748839SVaclav Hapla .seealso: DMPlexInterpolate(), DMPlexIsInterpolated() 1740a7748839SVaclav Hapla @*/ 1741a7748839SVaclav Hapla PetscErrorCode DMPlexIsInterpolatedCollective(DM dm, DMPlexInterpolatedFlag *interpolated) 1742a7748839SVaclav Hapla { 1743a7748839SVaclav Hapla DM_Plex *plex = (DM_Plex *) dm->data; 1744a7748839SVaclav Hapla PetscBool debug=PETSC_FALSE; 1745a7748839SVaclav Hapla PetscErrorCode ierr; 1746a7748839SVaclav Hapla 1747a7748839SVaclav Hapla PetscFunctionBegin; 1748a7748839SVaclav Hapla PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1749a7748839SVaclav Hapla PetscValidPointer(interpolated,2); 1750a7748839SVaclav Hapla ierr = PetscOptionsGetBool(((PetscObject) dm)->options, ((PetscObject) dm)->prefix, "-dm_plex_is_interpolated_collective_debug", &debug, NULL);CHKERRQ(ierr); 1751a7748839SVaclav Hapla if (plex->interpolatedCollective < 0) { 1752a7748839SVaclav Hapla DMPlexInterpolatedFlag min, max; 1753a7748839SVaclav Hapla MPI_Comm comm; 1754a7748839SVaclav Hapla 1755a7748839SVaclav Hapla ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 1756a7748839SVaclav Hapla ierr = DMPlexIsInterpolated(dm, &plex->interpolatedCollective);CHKERRQ(ierr); 1757a7748839SVaclav Hapla ierr = MPI_Allreduce(&plex->interpolatedCollective, &min, 1, MPIU_ENUM, MPI_MIN, comm);CHKERRQ(ierr); 1758a7748839SVaclav Hapla ierr = MPI_Allreduce(&plex->interpolatedCollective, &max, 1, MPIU_ENUM, MPI_MAX, comm);CHKERRQ(ierr); 1759a7748839SVaclav Hapla if (min != max) plex->interpolatedCollective = DMPLEX_INTERPOLATED_MIXED; 1760a7748839SVaclav Hapla if (debug) { 1761a7748839SVaclav Hapla PetscMPIInt rank; 1762a7748839SVaclav Hapla 1763a7748839SVaclav Hapla ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1764a7748839SVaclav Hapla ierr = PetscSynchronizedPrintf(comm, "[%d] interpolated=%s interpolatedCollective=%s\n", rank, DMPlexInterpolatedFlags[plex->interpolated], DMPlexInterpolatedFlags[plex->interpolatedCollective]);CHKERRQ(ierr); 1765a7748839SVaclav Hapla ierr = PetscSynchronizedFlush(comm, PETSC_STDOUT);CHKERRQ(ierr); 1766a7748839SVaclav Hapla } 1767a7748839SVaclav Hapla } 1768a7748839SVaclav Hapla *interpolated = plex->interpolatedCollective; 1769a7748839SVaclav Hapla PetscFunctionReturn(0); 1770a7748839SVaclav Hapla } 1771