1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2e8f14785SLisandro Dalcin #include <petsc/private/hashmapi.h> 3e8f14785SLisandro Dalcin #include <petsc/private/hashmapij.h> 4e8f14785SLisandro Dalcin 5e8f14785SLisandro Dalcin /* HashIJKL */ 6e8f14785SLisandro Dalcin 7e8f14785SLisandro Dalcin #include <petsc/private/hashmap.h> 8e8f14785SLisandro Dalcin 9e8f14785SLisandro Dalcin typedef struct _PetscHashIJKLKey { PetscInt i, j, k, l; } PetscHashIJKLKey; 10e8f14785SLisandro Dalcin 11e8f14785SLisandro Dalcin #define PetscHashIJKLKeyHash(key) \ 12e8f14785SLisandro Dalcin PetscHashCombine(PetscHashCombine(PetscHashInt((key).i),PetscHashInt((key).j)), \ 13e8f14785SLisandro Dalcin PetscHashCombine(PetscHashInt((key).k),PetscHashInt((key).l))) 14e8f14785SLisandro Dalcin 15e8f14785SLisandro Dalcin #define PetscHashIJKLKeyEqual(k1,k2) \ 16e8f14785SLisandro Dalcin (((k1).i==(k2).i) ? ((k1).j==(k2).j) ? ((k1).k==(k2).k) ? ((k1).l==(k2).l) : 0 : 0 : 0) 17e8f14785SLisandro Dalcin 18e8f14785SLisandro Dalcin PETSC_HASH_MAP(HashIJKL, PetscHashIJKLKey, PetscInt, PetscHashIJKLKeyHash, PetscHashIJKLKeyEqual, -1) 19e8f14785SLisandro Dalcin 209f074e33SMatthew G Knepley 219f074e33SMatthew G Knepley /* 229a5b13a2SMatthew G Knepley DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell 23cd921712SStefano Zampini This assumes that the mesh is not interpolated from the depth of point p to the vertices 249f074e33SMatthew G Knepley */ 25439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 269f074e33SMatthew G Knepley { 279f074e33SMatthew G Knepley const PetscInt *cone = NULL; 28cd921712SStefano Zampini PetscInt coneSize; 299f074e33SMatthew G Knepley PetscErrorCode ierr; 309f074e33SMatthew G Knepley 319f074e33SMatthew G Knepley PetscFunctionBegin; 329f074e33SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 339f074e33SMatthew G Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 349f074e33SMatthew G Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 35439ece16SMatthew G. Knepley ierr = DMPlexGetRawFaces_Internal(dm, dim, coneSize, cone, numFaces, faceSize, faces);CHKERRQ(ierr); 36439ece16SMatthew G. Knepley PetscFunctionReturn(0); 37439ece16SMatthew G. Knepley } 38439ece16SMatthew G. Knepley 39439ece16SMatthew G. Knepley /* 40439ece16SMatthew G. Knepley DMPlexRestoreFaces_Internal - Restores the array 41439ece16SMatthew G. Knepley */ 42439ece16SMatthew G. Knepley PetscErrorCode DMPlexRestoreFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 43439ece16SMatthew G. Knepley { 44439ece16SMatthew G. Knepley PetscErrorCode ierr; 45439ece16SMatthew G. Knepley 46439ece16SMatthew G. Knepley PetscFunctionBegin; 47cd921712SStefano Zampini if (faces) { ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);CHKERRQ(ierr); } 48439ece16SMatthew G. Knepley PetscFunctionReturn(0); 49439ece16SMatthew G. Knepley } 50439ece16SMatthew G. Knepley 51439ece16SMatthew G. Knepley /* 52439ece16SMatthew G. Knepley DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone 53439ece16SMatthew G. Knepley */ 54439ece16SMatthew G. Knepley PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 55439ece16SMatthew G. Knepley { 56439ece16SMatthew G. Knepley PetscInt *facesTmp; 57439ece16SMatthew G. Knepley PetscInt maxConeSize, maxSupportSize; 58439ece16SMatthew G. Knepley PetscErrorCode ierr; 59439ece16SMatthew G. Knepley 60439ece16SMatthew G. Knepley PetscFunctionBegin; 61439ece16SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 62cd921712SStefano Zampini if (faces && coneSize) PetscValidIntPointer(cone,4); 63439ece16SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 6469291d52SBarry Smith if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);CHKERRQ(ierr);} 659f074e33SMatthew G Knepley switch (dim) { 66c49d9212SMatthew G. Knepley case 1: 67c49d9212SMatthew G. Knepley switch (coneSize) { 68c49d9212SMatthew G. Knepley case 2: 69c49d9212SMatthew G. Knepley if (faces) { 70c49d9212SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 71c49d9212SMatthew G. Knepley *faces = facesTmp; 72c49d9212SMatthew G. Knepley } 73c49d9212SMatthew G. Knepley if (numFaces) *numFaces = 2; 74c49d9212SMatthew G. Knepley if (faceSize) *faceSize = 1; 75c49d9212SMatthew G. Knepley break; 76c49d9212SMatthew G. Knepley default: 7799836e0eSStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 78c49d9212SMatthew G. Knepley } 79c49d9212SMatthew G. Knepley break; 809f074e33SMatthew G Knepley case 2: 819f074e33SMatthew G Knepley switch (coneSize) { 829f074e33SMatthew G Knepley case 3: 839a5b13a2SMatthew G Knepley if (faces) { 849a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 859a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 869a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 879a5b13a2SMatthew G Knepley *faces = facesTmp; 889a5b13a2SMatthew G Knepley } 899a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 3; 909a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 919f074e33SMatthew G Knepley break; 929f074e33SMatthew G Knepley case 4: 939a5b13a2SMatthew G Knepley /* Vertices follow right hand rule */ 949a5b13a2SMatthew G Knepley if (faces) { 959a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 969a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 979a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[3]; 989a5b13a2SMatthew G Knepley facesTmp[6] = cone[3]; facesTmp[7] = cone[0]; 999a5b13a2SMatthew G Knepley *faces = facesTmp; 1009a5b13a2SMatthew G Knepley } 1019a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 4; 1029a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 1039f074e33SMatthew G Knepley break; 1049f074e33SMatthew G Knepley default: 10599836e0eSStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 1069f074e33SMatthew G Knepley } 1079f074e33SMatthew G Knepley break; 1089f074e33SMatthew G Knepley case 3: 1099f074e33SMatthew G Knepley switch (coneSize) { 1109f074e33SMatthew G Knepley case 3: 1119a5b13a2SMatthew G Knepley if (faces) { 1129a5b13a2SMatthew G Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 1139a5b13a2SMatthew G Knepley facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 1149a5b13a2SMatthew G Knepley facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 1159a5b13a2SMatthew G Knepley *faces = facesTmp; 1169a5b13a2SMatthew G Knepley } 1179a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 3; 1189a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 2; 1199f074e33SMatthew G Knepley break; 1209f074e33SMatthew G Knepley case 4: 1212e1b13c2SMatthew G. Knepley /* Vertices of first face follow right hand rule and normal points away from last vertex */ 1229a5b13a2SMatthew G Knepley if (faces) { 1232e1b13c2SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; 1242e1b13c2SMatthew G. Knepley facesTmp[3] = cone[0]; facesTmp[4] = cone[3]; facesTmp[5] = cone[1]; 1252e1b13c2SMatthew G. Knepley facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[3]; 1262e1b13c2SMatthew G. Knepley facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3]; 1279a5b13a2SMatthew G Knepley *faces = facesTmp; 1289a5b13a2SMatthew G Knepley } 1299a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 4; 1309a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 3; 1319a5b13a2SMatthew G Knepley break; 1329a5b13a2SMatthew G Knepley case 8: 133e368ef20SMatthew G. Knepley /* 7--------6 134e368ef20SMatthew G. Knepley /| /| 135e368ef20SMatthew G. Knepley / | / | 136e368ef20SMatthew G. Knepley 4--------5 | 137e368ef20SMatthew G. Knepley | | | | 138e368ef20SMatthew G. Knepley | | | | 139e368ef20SMatthew G. Knepley | 1--------2 140e368ef20SMatthew G. Knepley | / | / 141e368ef20SMatthew G. Knepley |/ |/ 142e368ef20SMatthew G. Knepley 0--------3 143e368ef20SMatthew G. Knepley */ 1449a5b13a2SMatthew G Knepley if (faces) { 14551cfd6a4SMatthew G. Knepley facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; /* Bottom */ 14651cfd6a4SMatthew G. Knepley facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; /* Top */ 14751cfd6a4SMatthew G. Knepley facesTmp[8] = cone[0]; facesTmp[9] = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; /* Front */ 14851cfd6a4SMatthew G. Knepley facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; /* Back */ 14951cfd6a4SMatthew G. Knepley facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; /* Right */ 15051cfd6a4SMatthew G. Knepley facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; /* Left */ 1519a5b13a2SMatthew G Knepley *faces = facesTmp; 1529a5b13a2SMatthew G Knepley } 1539a5b13a2SMatthew G Knepley if (numFaces) *numFaces = 6; 1549a5b13a2SMatthew G Knepley if (faceSize) *faceSize = 4; 1559f074e33SMatthew G Knepley break; 1569f074e33SMatthew G Knepley default: 15799836e0eSStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 1589f074e33SMatthew G Knepley } 1599f074e33SMatthew G Knepley break; 1609f074e33SMatthew G Knepley default: 16199836e0eSStefano Zampini SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim); 1629f074e33SMatthew G Knepley } 1639f074e33SMatthew G Knepley PetscFunctionReturn(0); 1649f074e33SMatthew G Knepley } 1659f074e33SMatthew G Knepley 16699836e0eSStefano Zampini /* 16799836e0eSStefano Zampini DMPlexGetRawFacesHybrid_Internal - Gets groups of vertices that correspond to faces for the given cone using hybrid ordering (prisms) 16899836e0eSStefano Zampini */ 16999836e0eSStefano Zampini static PetscErrorCode DMPlexGetRawFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[]) 17099836e0eSStefano Zampini { 17199836e0eSStefano Zampini PetscInt *facesTmp; 17299836e0eSStefano Zampini PetscInt maxConeSize, maxSupportSize; 17399836e0eSStefano Zampini PetscErrorCode ierr; 17499836e0eSStefano Zampini 17599836e0eSStefano Zampini PetscFunctionBegin; 17699836e0eSStefano Zampini PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 17799836e0eSStefano Zampini if (faces && coneSize) PetscValidIntPointer(cone,4); 17899836e0eSStefano Zampini ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 17999836e0eSStefano Zampini if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);CHKERRQ(ierr);} 18099836e0eSStefano Zampini switch (dim) { 18199836e0eSStefano Zampini case 1: 18299836e0eSStefano Zampini switch (coneSize) { 18399836e0eSStefano Zampini case 2: 18499836e0eSStefano Zampini if (faces) { 18599836e0eSStefano Zampini facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 18699836e0eSStefano Zampini *faces = facesTmp; 18799836e0eSStefano Zampini } 18899836e0eSStefano Zampini if (numFaces) *numFaces = 2; 18999836e0eSStefano Zampini if (numFacesNotH) *numFacesNotH = 2; 19099836e0eSStefano Zampini if (faceSize) *faceSize = 1; 19199836e0eSStefano Zampini break; 19299836e0eSStefano Zampini default: 19399836e0eSStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 19499836e0eSStefano Zampini } 19599836e0eSStefano Zampini break; 19699836e0eSStefano Zampini case 2: 19799836e0eSStefano Zampini switch (coneSize) { 19899836e0eSStefano Zampini case 4: 19999836e0eSStefano Zampini if (faces) { 20099836e0eSStefano Zampini facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 20199836e0eSStefano Zampini facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; 20299836e0eSStefano Zampini facesTmp[4] = cone[0]; facesTmp[5] = cone[2]; 20399836e0eSStefano Zampini facesTmp[6] = cone[1]; facesTmp[7] = cone[3]; 20499836e0eSStefano Zampini *faces = facesTmp; 20599836e0eSStefano Zampini } 20699836e0eSStefano Zampini if (numFaces) *numFaces = 4; 20799836e0eSStefano Zampini if (numFacesNotH) *numFacesNotH = 2; 20899836e0eSStefano Zampini if (faceSize) *faceSize = 2; 20999836e0eSStefano Zampini break; 21099836e0eSStefano Zampini default: 21199836e0eSStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 21299836e0eSStefano Zampini } 21399836e0eSStefano Zampini break; 21499836e0eSStefano Zampini case 3: 21599836e0eSStefano Zampini switch (coneSize) { 21699836e0eSStefano Zampini case 6: /* triangular prism */ 21799836e0eSStefano Zampini if (faces) { 21899836e0eSStefano Zampini facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = -1; /* Bottom */ 21999836e0eSStefano Zampini facesTmp[4] = cone[3]; facesTmp[5] = cone[4]; facesTmp[6] = cone[5]; facesTmp[7] = -1; /* Top */ 22099836e0eSStefano Zampini facesTmp[8] = cone[0]; facesTmp[9] = cone[1]; facesTmp[10] = cone[3]; facesTmp[11] = cone[4]; /* Back left */ 22199836e0eSStefano Zampini facesTmp[12] = cone[1]; facesTmp[13] = cone[2]; facesTmp[14] = cone[4]; facesTmp[15] = cone[5]; /* Back right */ 22299836e0eSStefano Zampini facesTmp[16] = cone[2]; facesTmp[17] = cone[0]; facesTmp[18] = cone[5]; facesTmp[19] = cone[3]; /* Front */ 22399836e0eSStefano Zampini *faces = facesTmp; 22499836e0eSStefano Zampini } 22599836e0eSStefano Zampini if (numFaces) *numFaces = 5; 22699836e0eSStefano Zampini if (numFacesNotH) *numFacesNotH = 2; 22799836e0eSStefano Zampini if (faceSize) *faceSize = -4; 22899836e0eSStefano Zampini break; 22999836e0eSStefano Zampini default: 23099836e0eSStefano Zampini SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 23199836e0eSStefano Zampini } 23299836e0eSStefano Zampini break; 23399836e0eSStefano Zampini default: 23499836e0eSStefano Zampini SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim); 23599836e0eSStefano Zampini } 23699836e0eSStefano Zampini PetscFunctionReturn(0); 23799836e0eSStefano Zampini } 23899836e0eSStefano Zampini 23999836e0eSStefano Zampini static PetscErrorCode DMPlexRestoreRawFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[]) 24099836e0eSStefano Zampini { 24199836e0eSStefano Zampini PetscErrorCode ierr; 24299836e0eSStefano Zampini 24399836e0eSStefano Zampini PetscFunctionBegin; 24499836e0eSStefano Zampini if (faces) { ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);CHKERRQ(ierr); } 24599836e0eSStefano Zampini PetscFunctionReturn(0); 24699836e0eSStefano Zampini } 24799836e0eSStefano Zampini 24899836e0eSStefano Zampini static PetscErrorCode DMPlexGetFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[]) 24999836e0eSStefano Zampini { 25099836e0eSStefano Zampini const PetscInt *cone = NULL; 25199836e0eSStefano Zampini PetscInt coneSize; 25299836e0eSStefano Zampini PetscErrorCode ierr; 25399836e0eSStefano Zampini 25499836e0eSStefano Zampini PetscFunctionBegin; 25599836e0eSStefano Zampini PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 25699836e0eSStefano Zampini ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 25799836e0eSStefano Zampini ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 25899836e0eSStefano Zampini ierr = DMPlexGetRawFacesHybrid_Internal(dm, dim, coneSize, cone, numFaces, numFacesNotH, faceSize, faces);CHKERRQ(ierr); 25999836e0eSStefano Zampini PetscFunctionReturn(0); 26099836e0eSStefano Zampini } 26199836e0eSStefano Zampini 2629a5b13a2SMatthew G Knepley /* This interpolates faces for cells at some stratum */ 2639a5b13a2SMatthew G Knepley static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm) 2649f074e33SMatthew G Knepley { 265d4efc6f1SMatthew G. Knepley DMLabel subpointMap; 2669a5b13a2SMatthew G Knepley PetscHashIJKL faceTable; 2679a5b13a2SMatthew G Knepley PetscInt *pStart, *pEnd; 2689a5b13a2SMatthew G Knepley PetscInt cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d; 269e9fa77a1SMatthew G. Knepley PetscInt coneSizeH = 0, faceSizeAllH = 0, faceSizeAllT = 0, numCellFacesH = 0, faceT = 0, faceH, pMax = -1, dim, outerloop; 27099836e0eSStefano Zampini PetscInt cMax, fMax, eMax, vMax; 2719f074e33SMatthew G Knepley PetscErrorCode ierr; 2729f074e33SMatthew G Knepley 2739f074e33SMatthew G Knepley PetscFunctionBegin; 274c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &cellDim);CHKERRQ(ierr); 275d4efc6f1SMatthew G. Knepley /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */ 276d4efc6f1SMatthew G. Knepley ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 277d4efc6f1SMatthew G. Knepley if (subpointMap) ++cellDim; 2789a5b13a2SMatthew G Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 2799a5b13a2SMatthew G Knepley ++depth; 2809a5b13a2SMatthew G Knepley ++cellDepth; 2819a5b13a2SMatthew G Knepley cellDim -= depth - cellDepth; 282dcca6d9dSJed Brown ierr = PetscMalloc2(depth+1,&pStart,depth+1,&pEnd);CHKERRQ(ierr); 2839a5b13a2SMatthew G Knepley for (d = depth-1; d >= faceDepth; --d) { 2849a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr); 2859f074e33SMatthew G Knepley } 2869a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr); 2879a5b13a2SMatthew G Knepley pEnd[faceDepth] = pStart[faceDepth]; 2889a5b13a2SMatthew G Knepley for (d = faceDepth-1; d >= 0; --d) { 2899a5b13a2SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 2909f074e33SMatthew G Knepley } 29199836e0eSStefano Zampini cMax = fMax = eMax = vMax = PETSC_DETERMINE; 29299836e0eSStefano Zampini ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 29399836e0eSStefano Zampini if (cellDim == dim) { 29499836e0eSStefano Zampini ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 29599836e0eSStefano Zampini pMax = cMax; 29699836e0eSStefano Zampini } else if (cellDim == dim -1) { 29799836e0eSStefano Zampini ierr = DMPlexGetHybridBounds(dm, &cMax, &fMax, NULL, NULL);CHKERRQ(ierr); 29899836e0eSStefano Zampini pMax = fMax; 29999836e0eSStefano Zampini } 30099836e0eSStefano Zampini pMax = pMax < 0 ? pEnd[cellDepth] : pMax; 30199836e0eSStefano Zampini if (pMax < pEnd[cellDepth]) { 30299836e0eSStefano Zampini const PetscInt *cellFaces, *cone; 30399836e0eSStefano Zampini PetscInt numCellFacesT, faceSize, cf; 3049f074e33SMatthew G Knepley 305e9fa77a1SMatthew G. Knepley /* First get normal cell face size (we now allow hybrid cells to meet normal cells on either hybrid or normal faces */ 306*5ebdaad1SMatthew G. Knepley if (pStart[cellDepth] < pMax) {ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);} 307e9fa77a1SMatthew G. Knepley 30899836e0eSStefano Zampini ierr = DMPlexGetConeSize(dm, pMax, &coneSizeH);CHKERRQ(ierr); 30999836e0eSStefano Zampini ierr = DMPlexGetCone(dm, pMax, &cone);CHKERRQ(ierr); 31099836e0eSStefano Zampini ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSizeH, cone, &numCellFacesH, &numCellFacesT, &faceSize, &cellFaces);CHKERRQ(ierr); 31199836e0eSStefano Zampini if (faceSize < 0) { 31299836e0eSStefano Zampini PetscInt *sizes, minv, maxv; 31399836e0eSStefano Zampini 31499836e0eSStefano Zampini /* count vertices of hybrid and non-hybrid faces */ 31599836e0eSStefano Zampini ierr = PetscCalloc1(numCellFacesH, &sizes);CHKERRQ(ierr); 31699836e0eSStefano Zampini for (cf = 0; cf < numCellFacesT; ++cf) { /* These are the non-hybrid faces */ 31799836e0eSStefano Zampini const PetscInt *cellFace = &cellFaces[-cf*faceSize]; 31899836e0eSStefano Zampini PetscInt f; 31999836e0eSStefano Zampini 32099836e0eSStefano Zampini for (f = 0; f < -faceSize; ++f) sizes[cf] += (cellFace[f] >= 0 ? 1 : 0); 32199836e0eSStefano Zampini } 32299836e0eSStefano Zampini ierr = PetscSortInt(numCellFacesT, sizes);CHKERRQ(ierr); 32399836e0eSStefano Zampini minv = sizes[0]; 32499836e0eSStefano Zampini maxv = sizes[PetscMax(numCellFacesT-1, 0)]; 32599836e0eSStefano Zampini if (minv != maxv) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Different number of vertices for non-hybrid face %D != %D", minv, maxv); 326e9fa77a1SMatthew G. Knepley faceSizeAllT = minv; 32799836e0eSStefano Zampini ierr = PetscMemzero(sizes, numCellFacesH*sizeof(PetscInt));CHKERRQ(ierr); 32899836e0eSStefano Zampini for (cf = numCellFacesT; cf < numCellFacesH; ++cf) { /* These are the hybrid faces */ 32999836e0eSStefano Zampini const PetscInt *cellFace = &cellFaces[-cf*faceSize]; 33099836e0eSStefano Zampini PetscInt f; 33199836e0eSStefano Zampini 33299836e0eSStefano Zampini for (f = 0; f < -faceSize; ++f) sizes[cf-numCellFacesT] += (cellFace[f] >= 0 ? 1 : 0); 33399836e0eSStefano Zampini } 33499836e0eSStefano Zampini ierr = PetscSortInt(numCellFacesH - numCellFacesT, sizes);CHKERRQ(ierr); 33599836e0eSStefano Zampini minv = sizes[0]; 33699836e0eSStefano Zampini maxv = sizes[PetscMax(numCellFacesH - numCellFacesT-1, 0)]; 33799836e0eSStefano Zampini ierr = PetscFree(sizes);CHKERRQ(ierr); 33899836e0eSStefano Zampini if (minv != maxv) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Different number of vertices for hybrid face %D != %D", minv, maxv); 33999836e0eSStefano Zampini faceSizeAllH = minv; 340*5ebdaad1SMatthew G. Knepley if (!faceSizeAll) faceSizeAll = faceSizeAllT; 34199836e0eSStefano Zampini } else { /* the size of the faces in hybrid cells is the same */ 342e9fa77a1SMatthew G. Knepley faceSizeAll = faceSizeAllH = faceSizeAllT = faceSize; 34399836e0eSStefano Zampini } 34499836e0eSStefano Zampini ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, cone, &numCellFacesH, &numCellFacesT, &faceSize, &cellFaces);CHKERRQ(ierr); 34599836e0eSStefano Zampini } else if (pEnd[cellDepth] > pStart[cellDepth]) { 34699836e0eSStefano Zampini ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr); 34799836e0eSStefano Zampini } 34899836e0eSStefano Zampini if (faceSizeAll > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 34999836e0eSStefano Zampini 350b135d7daSStefano Zampini /* With hybrid grids, we first iterate on hybrid cells and start numbering the non-hybrid faces 35199836e0eSStefano Zampini Then, faces for non-hybrid cells are numbered. 35299836e0eSStefano Zampini This is to guarantee consistent orientations (all 0) of all the points in the cone of the hybrid cells */ 35399836e0eSStefano Zampini ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 35499836e0eSStefano Zampini for (outerloop = 0, face = pStart[faceDepth]; outerloop < 2; outerloop++) { 35599836e0eSStefano Zampini PetscInt start, end; 35699836e0eSStefano Zampini 35799836e0eSStefano Zampini start = outerloop == 0 ? pMax : pStart[cellDepth]; 35899836e0eSStefano Zampini end = outerloop == 0 ? pEnd[cellDepth] : pMax; 35999836e0eSStefano Zampini for (c = start; c < end; ++c) { 36099836e0eSStefano Zampini const PetscInt *cellFaces; 361e9fa77a1SMatthew G. Knepley PetscInt numCellFaces, faceSize, faceSizeInc, faceSizeCheck, cf; 36299836e0eSStefano Zampini 36399836e0eSStefano Zampini if (c < pMax) { 3649a5b13a2SMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 3659a5b13a2SMatthew G Knepley if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 366e9fa77a1SMatthew G. Knepley faceSizeCheck = faceSizeAll; 36799836e0eSStefano Zampini } else { /* Hybrid cell */ 36899836e0eSStefano Zampini const PetscInt *cone; 36999836e0eSStefano Zampini PetscInt numCellFacesN, coneSize; 37099836e0eSStefano Zampini 37199836e0eSStefano Zampini ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 37299836e0eSStefano Zampini if (coneSize != coneSizeH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid coneSize %D != %D", coneSize, coneSizeH); 37399836e0eSStefano Zampini ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 37499836e0eSStefano Zampini ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr); 37599836e0eSStefano Zampini if (numCellFaces != numCellFacesH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected numCellFaces %D != %D for hybrid cell %D", numCellFaces, numCellFacesH, c); 37699836e0eSStefano Zampini faceSize = PetscMax(faceSize, -faceSize); 37799836e0eSStefano Zampini if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSize); 37899836e0eSStefano Zampini numCellFaces = numCellFacesN; /* process only non-hybrid faces */ 379e9fa77a1SMatthew G. Knepley faceSizeCheck = faceSizeAllT; 38099836e0eSStefano Zampini } 38199836e0eSStefano Zampini faceSizeInc = faceSize; 3829f074e33SMatthew G Knepley for (cf = 0; cf < numCellFaces; ++cf) { 38399836e0eSStefano Zampini const PetscInt *cellFace = &cellFaces[cf*faceSizeInc]; 38499836e0eSStefano Zampini PetscInt faceSizeH = faceSize; 3859a5b13a2SMatthew G Knepley PetscHashIJKLKey key; 386e8f14785SLisandro Dalcin PetscHashIter iter; 387e8f14785SLisandro Dalcin PetscBool missing; 3889f074e33SMatthew G Knepley 38999836e0eSStefano Zampini if (faceSizeInc == 2) { 3909a5b13a2SMatthew G Knepley key.i = PetscMin(cellFace[0], cellFace[1]); 3919a5b13a2SMatthew G Knepley key.j = PetscMax(cellFace[0], cellFace[1]); 39299836e0eSStefano Zampini key.k = PETSC_MAX_INT; 39399836e0eSStefano Zampini key.l = PETSC_MAX_INT; 3949a5b13a2SMatthew G Knepley } else { 39599836e0eSStefano Zampini key.i = cellFace[0]; 39699836e0eSStefano Zampini key.j = cellFace[1]; 39799836e0eSStefano Zampini key.k = cellFace[2]; 39899836e0eSStefano Zampini key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT; 399302440fdSBarry Smith ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 4009f074e33SMatthew G Knepley } 40199836e0eSStefano Zampini /* this check is redundant for non-hybrid meshes */ 402e9fa77a1SMatthew G. Knepley if (faceSizeH != faceSizeCheck) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected number of vertices for face %D of point %D -> %D != %D", cf, c, faceSizeH, faceSizeCheck); 403e8f14785SLisandro Dalcin ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 404e9fa77a1SMatthew G. Knepley if (missing) { 405e9fa77a1SMatthew G. Knepley ierr = PetscHashIJKLIterSet(faceTable, iter, face++);CHKERRQ(ierr); 406e9fa77a1SMatthew G. Knepley if (c >= pMax) ++faceT; 407e9fa77a1SMatthew G. Knepley } 4089a5b13a2SMatthew G Knepley } 40999836e0eSStefano Zampini if (c < pMax) { 410439ece16SMatthew G. Knepley ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 41199836e0eSStefano Zampini } else { 41299836e0eSStefano Zampini ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, NULL, NULL, NULL, NULL, &cellFaces);CHKERRQ(ierr); 41399836e0eSStefano Zampini } 41499836e0eSStefano Zampini } 41599836e0eSStefano Zampini } 41699836e0eSStefano Zampini pEnd[faceDepth] = face; 41799836e0eSStefano Zampini 41899836e0eSStefano Zampini /* Second pass for hybrid meshes: number hybrid faces */ 41999836e0eSStefano Zampini for (c = pMax; c < pEnd[cellDepth]; ++c) { 42099836e0eSStefano Zampini const PetscInt *cellFaces, *cone; 42199836e0eSStefano Zampini PetscInt numCellFaces, numCellFacesN, faceSize, cf, coneSize; 42299836e0eSStefano Zampini 42399836e0eSStefano Zampini ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 42499836e0eSStefano Zampini ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 42599836e0eSStefano Zampini ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr); 42699836e0eSStefano Zampini if (numCellFaces != numCellFacesH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid numCellFaces %D != %D", numCellFaces, numCellFacesH); 42799836e0eSStefano Zampini faceSize = PetscMax(faceSize, -faceSize); 42899836e0eSStefano Zampini for (cf = numCellFacesN; cf < numCellFaces; ++cf) { /* These are the hybrid faces */ 42999836e0eSStefano Zampini const PetscInt *cellFace = &cellFaces[cf*faceSize]; 43099836e0eSStefano Zampini PetscHashIJKLKey key; 43199836e0eSStefano Zampini PetscHashIter iter; 43299836e0eSStefano Zampini PetscBool missing; 43399836e0eSStefano Zampini PetscInt faceSizeH = faceSize; 43499836e0eSStefano Zampini 43599836e0eSStefano Zampini if (faceSize == 2) { 43699836e0eSStefano Zampini key.i = PetscMin(cellFace[0], cellFace[1]); 43799836e0eSStefano Zampini key.j = PetscMax(cellFace[0], cellFace[1]); 43899836e0eSStefano Zampini key.k = PETSC_MAX_INT; 43999836e0eSStefano Zampini key.l = PETSC_MAX_INT; 44099836e0eSStefano Zampini } else { 44199836e0eSStefano Zampini key.i = cellFace[0]; 44299836e0eSStefano Zampini key.j = cellFace[1]; 44399836e0eSStefano Zampini key.k = cellFace[2]; 44499836e0eSStefano Zampini key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT; 44599836e0eSStefano Zampini ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 44699836e0eSStefano Zampini } 44799836e0eSStefano Zampini if (faceSizeH != faceSizeAllH) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected number of vertices for hybrid face %D of point %D -> %D != %D", cf, c, faceSizeH, faceSizeAllH); 44899836e0eSStefano Zampini ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 44999836e0eSStefano Zampini if (missing) {ierr = PetscHashIJKLIterSet(faceTable, iter, face++);CHKERRQ(ierr);} 45099836e0eSStefano Zampini } 45199836e0eSStefano Zampini ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr); 45299836e0eSStefano Zampini } 45399836e0eSStefano Zampini faceH = face - pEnd[faceDepth]; 45499836e0eSStefano Zampini if (faceH) { 45599836e0eSStefano Zampini if (fMax == PETSC_DETERMINE) fMax = pEnd[faceDepth]; 45699836e0eSStefano Zampini else if (eMax == PETSC_DETERMINE) eMax = pEnd[faceDepth]; 45799836e0eSStefano Zampini else SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of unassigned hybrid facets %D for cellDim %D and dimension %D", faceH, cellDim, dim); 4589a5b13a2SMatthew G Knepley } 4599a5b13a2SMatthew G Knepley pEnd[faceDepth] = face; 4609a5b13a2SMatthew G Knepley ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 4619a5b13a2SMatthew G Knepley /* Count new points */ 4629a5b13a2SMatthew G Knepley for (d = 0; d <= depth; ++d) { 4639a5b13a2SMatthew G Knepley numPoints += pEnd[d]-pStart[d]; 4649a5b13a2SMatthew G Knepley } 4659a5b13a2SMatthew G Knepley ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr); 4669a5b13a2SMatthew G Knepley /* Set cone sizes */ 4679a5b13a2SMatthew G Knepley for (d = 0; d <= depth; ++d) { 4689a5b13a2SMatthew G Knepley PetscInt coneSize, p; 4699f074e33SMatthew G Knepley 4709a5b13a2SMatthew G Knepley if (d == faceDepth) { 471e9fa77a1SMatthew G. Knepley /* Now we have two cases: */ 472e9fa77a1SMatthew G. Knepley if (faceSizeAll == faceSizeAllT) { 4739a5b13a2SMatthew G Knepley /* I see no way to do this if we admit faces of different shapes */ 47499836e0eSStefano Zampini for (p = pStart[d]; p < pEnd[d]-faceH; ++p) { 4759a5b13a2SMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr); 4769a5b13a2SMatthew G Knepley } 47799836e0eSStefano Zampini for (p = pEnd[d]-faceH; p < pEnd[d]; ++p) { 47899836e0eSStefano Zampini ierr = DMPlexSetConeSize(idm, p, faceSizeAllH);CHKERRQ(ierr); 47999836e0eSStefano Zampini } 480e9fa77a1SMatthew G. Knepley } else if (faceSizeAll == faceSizeAllH) { 481e9fa77a1SMatthew G. Knepley for (p = pStart[d]; p < pStart[d]+faceT; ++p) { 482e9fa77a1SMatthew G. Knepley ierr = DMPlexSetConeSize(idm, p, faceSizeAllT);CHKERRQ(ierr); 483e9fa77a1SMatthew G. Knepley } 484e9fa77a1SMatthew G. Knepley for (p = pStart[d]+faceT; p < pEnd[d]-faceH; ++p) { 485e9fa77a1SMatthew G. Knepley ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr); 486e9fa77a1SMatthew G. Knepley } 487e9fa77a1SMatthew G. Knepley for (p = pEnd[d]-faceH; p < pEnd[d]; ++p) { 488e9fa77a1SMatthew G. Knepley ierr = DMPlexSetConeSize(idm, p, faceSizeAllH);CHKERRQ(ierr); 489e9fa77a1SMatthew G. Knepley } 490e9fa77a1SMatthew G. Knepley } else SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent faces sizes N: %D T: %D H: %D", faceSizeAll, faceSizeAllT, faceSizeAllH); 491a014044eSMatthew G Knepley } else if (d == cellDepth) { 492a014044eSMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 493a014044eSMatthew G Knepley /* Number of cell faces may be different from number of cell vertices*/ 49499836e0eSStefano Zampini if (p < pMax) { 495a014044eSMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr); 49699836e0eSStefano Zampini } else { 49799836e0eSStefano Zampini ierr = DMPlexGetFacesHybrid_Internal(dm, cellDim, p, &coneSize, NULL, NULL, NULL);CHKERRQ(ierr); 49899836e0eSStefano Zampini } 499a014044eSMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 500a014044eSMatthew G Knepley } 5019a5b13a2SMatthew G Knepley } else { 5029a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 5039a5b13a2SMatthew G Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 5049a5b13a2SMatthew G Knepley ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 5059f074e33SMatthew G Knepley } 5069f074e33SMatthew G Knepley } 5079f074e33SMatthew G Knepley } 5089f074e33SMatthew G Knepley ierr = DMSetUp(idm);CHKERRQ(ierr); 5099f074e33SMatthew G Knepley /* Get face cones from subsets of cell vertices */ 51099836e0eSStefano Zampini if (faceSizeAll > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 5119a5b13a2SMatthew G Knepley ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 5129a5b13a2SMatthew G Knepley for (d = depth; d > cellDepth; --d) { 5139a5b13a2SMatthew G Knepley const PetscInt *cone; 5149a5b13a2SMatthew G Knepley PetscInt p; 5159a5b13a2SMatthew G Knepley 5169a5b13a2SMatthew G Knepley for (p = pStart[d]; p < pEnd[d]; ++p) { 5179a5b13a2SMatthew G Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 5189a5b13a2SMatthew G Knepley ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr); 5199a5b13a2SMatthew G Knepley ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr); 5209a5b13a2SMatthew G Knepley ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr); 5219f074e33SMatthew G Knepley } 5229a5b13a2SMatthew G Knepley } 52399836e0eSStefano Zampini for (outerloop = 0, face = pStart[faceDepth]; outerloop < 2; outerloop++) { 52499836e0eSStefano Zampini PetscInt start, end; 5259f074e33SMatthew G Knepley 52699836e0eSStefano Zampini start = outerloop == 0 ? pMax : pStart[cellDepth]; 52799836e0eSStefano Zampini end = outerloop == 0 ? pEnd[cellDepth] : pMax; 52899836e0eSStefano Zampini for (c = start; c < end; ++c) { 52999836e0eSStefano Zampini const PetscInt *cellFaces; 53099836e0eSStefano Zampini PetscInt numCellFaces, faceSize, faceSizeInc, cf; 53199836e0eSStefano Zampini 53299836e0eSStefano Zampini if (c < pMax) { 5339a5b13a2SMatthew G Knepley ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 5349a5b13a2SMatthew G Knepley if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 53599836e0eSStefano Zampini } else { 53699836e0eSStefano Zampini const PetscInt *cone; 53799836e0eSStefano Zampini PetscInt numCellFacesN, coneSize; 53899836e0eSStefano Zampini 53999836e0eSStefano Zampini ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 54099836e0eSStefano Zampini if (coneSize != coneSizeH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid coneSize %D != %D", coneSize, coneSizeH); 54199836e0eSStefano Zampini ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 54299836e0eSStefano Zampini ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr); 54399836e0eSStefano Zampini if (numCellFaces != numCellFacesH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected numCellFaces %D != %D for hybrid cell %D", numCellFaces, numCellFacesH, c); 54499836e0eSStefano Zampini faceSize = PetscMax(faceSize, -faceSize); 54599836e0eSStefano Zampini if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSize); 54699836e0eSStefano Zampini numCellFaces = numCellFacesN; /* process only non-hybrid faces */ 54799836e0eSStefano Zampini } 54899836e0eSStefano Zampini faceSizeInc = faceSize; 5499f074e33SMatthew G Knepley for (cf = 0; cf < numCellFaces; ++cf) { 55099836e0eSStefano Zampini const PetscInt *cellFace = &cellFaces[cf*faceSizeInc]; 5519a5b13a2SMatthew G Knepley PetscHashIJKLKey key; 552e8f14785SLisandro Dalcin PetscHashIter iter; 553e8f14785SLisandro Dalcin PetscBool missing; 5549f074e33SMatthew G Knepley 55599836e0eSStefano Zampini if (faceSizeInc == 2) { 5569a5b13a2SMatthew G Knepley key.i = PetscMin(cellFace[0], cellFace[1]); 5579a5b13a2SMatthew G Knepley key.j = PetscMax(cellFace[0], cellFace[1]); 55899836e0eSStefano Zampini key.k = PETSC_MAX_INT; 55999836e0eSStefano Zampini key.l = PETSC_MAX_INT; 5609a5b13a2SMatthew G Knepley } else { 561e8f14785SLisandro Dalcin key.i = cellFace[0]; 562e8f14785SLisandro Dalcin key.j = cellFace[1]; 563e8f14785SLisandro Dalcin key.k = cellFace[2]; 56499836e0eSStefano Zampini key.l = faceSizeInc > 3 ? (cellFace[3] < 0 ? faceSize = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT; 56599836e0eSStefano Zampini ierr = PetscSortInt(faceSizeInc, (PetscInt *) &key);CHKERRQ(ierr); 5669f074e33SMatthew G Knepley } 567e8f14785SLisandro Dalcin ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 568735a0255SMatthew G. Knepley if (missing) { 5699a5b13a2SMatthew G Knepley ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr); 570e8f14785SLisandro Dalcin ierr = PetscHashIJKLIterSet(faceTable, iter, face);CHKERRQ(ierr); 571735a0255SMatthew G. Knepley ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr); 5729a5b13a2SMatthew G Knepley } else { 5739a5b13a2SMatthew G Knepley const PetscInt *cone; 574735a0255SMatthew G. Knepley PetscInt coneSize, ornt, i, j, f; 5759f074e33SMatthew G Knepley 576e8f14785SLisandro Dalcin ierr = PetscHashIJKLIterGet(faceTable, iter, &f);CHKERRQ(ierr); 5779a5b13a2SMatthew G Knepley ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 5782e1b13c2SMatthew G. Knepley /* Orient face: Do not allow reverse orientation at the first vertex */ 5799f074e33SMatthew G Knepley ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 5809f074e33SMatthew G Knepley ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 5819a5b13a2SMatthew 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); 5829a5b13a2SMatthew G Knepley /* - First find the initial vertex */ 5839a5b13a2SMatthew G Knepley for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break; 5849a5b13a2SMatthew G Knepley /* - Try forward comparison */ 5859a5b13a2SMatthew G Knepley for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break; 5869a5b13a2SMatthew G Knepley if (j == faceSize) { 5879a5b13a2SMatthew G Knepley if ((faceSize == 2) && (i == 1)) ornt = -2; 5889a5b13a2SMatthew G Knepley else ornt = i; 5899a5b13a2SMatthew G Knepley } else { 5909a5b13a2SMatthew G Knepley /* - Try backward comparison */ 5919a5b13a2SMatthew G Knepley for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break; 5922e1b13c2SMatthew G. Knepley if (j == faceSize) { 5932e1b13c2SMatthew G. Knepley if (i == 0) ornt = -faceSize; 594dc1a705cSMatthew G. Knepley else ornt = -i; 595e9fa77a1SMatthew G. Knepley } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine orientation of face %D in cell %D", f, c); 5969a5b13a2SMatthew G Knepley } 5979a5b13a2SMatthew G Knepley ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 5989f074e33SMatthew G Knepley } 5999f074e33SMatthew G Knepley } 60099836e0eSStefano Zampini if (c < pMax) { 601439ece16SMatthew G. Knepley ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 60299836e0eSStefano Zampini } else { 60399836e0eSStefano Zampini ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, NULL, NULL, NULL, NULL, &cellFaces);CHKERRQ(ierr); 6049f074e33SMatthew G Knepley } 60599836e0eSStefano Zampini } 60699836e0eSStefano Zampini } 60799836e0eSStefano Zampini /* Second pass for hybrid meshes: orient hybrid faces */ 60899836e0eSStefano Zampini for (c = pMax; c < pEnd[cellDepth]; ++c) { 60999836e0eSStefano Zampini const PetscInt *cellFaces, *cone; 61099836e0eSStefano Zampini PetscInt numCellFaces, numCellFacesN, faceSize, cf, coneSize; 61199836e0eSStefano Zampini 61299836e0eSStefano Zampini ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 61399836e0eSStefano Zampini ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 61499836e0eSStefano Zampini ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr); 61599836e0eSStefano Zampini if (numCellFaces != numCellFacesH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid numCellFaces %D != %D", numCellFaces, numCellFacesH); 61699836e0eSStefano Zampini faceSize = PetscMax(faceSize, -faceSize); 61799836e0eSStefano Zampini for (cf = numCellFacesN; cf < numCellFaces; ++cf) { /* These are the hybrid faces */ 61899836e0eSStefano Zampini const PetscInt *cellFace = &cellFaces[cf*faceSize]; 61999836e0eSStefano Zampini PetscHashIJKLKey key; 62099836e0eSStefano Zampini PetscHashIter iter; 62199836e0eSStefano Zampini PetscBool missing; 62299836e0eSStefano Zampini PetscInt faceSizeH = faceSize; 62399836e0eSStefano Zampini 62499836e0eSStefano Zampini if (faceSize == 2) { 62599836e0eSStefano Zampini key.i = PetscMin(cellFace[0], cellFace[1]); 62699836e0eSStefano Zampini key.j = PetscMax(cellFace[0], cellFace[1]); 62799836e0eSStefano Zampini key.k = PETSC_MAX_INT; 62899836e0eSStefano Zampini key.l = PETSC_MAX_INT; 62999836e0eSStefano Zampini } else { 63099836e0eSStefano Zampini key.i = cellFace[0]; 63199836e0eSStefano Zampini key.j = cellFace[1]; 63299836e0eSStefano Zampini key.k = cellFace[2]; 63399836e0eSStefano Zampini key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT; 63499836e0eSStefano Zampini ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 63599836e0eSStefano Zampini } 63699836e0eSStefano Zampini if (faceSizeH != faceSizeAllH) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected number of vertices for hybrid face %D of point %D -> %D != %D", cf, c, faceSizeH, faceSizeAllH); 63799836e0eSStefano Zampini ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 63899836e0eSStefano Zampini if (missing) { 63999836e0eSStefano Zampini ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr); 64099836e0eSStefano Zampini ierr = PetscHashIJKLIterSet(faceTable, iter, face);CHKERRQ(ierr); 64199836e0eSStefano Zampini ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr); 64299836e0eSStefano Zampini } else { 643e9fa77a1SMatthew G. Knepley PetscInt fv[4] = {0, 1, 2, 3}; 64499836e0eSStefano Zampini const PetscInt *cone; 64599836e0eSStefano Zampini PetscInt coneSize, ornt, i, j, f; 64699836e0eSStefano Zampini 64799836e0eSStefano Zampini ierr = PetscHashIJKLIterGet(faceTable, iter, &f);CHKERRQ(ierr); 64899836e0eSStefano Zampini ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 64999836e0eSStefano Zampini /* Orient face: Do not allow reverse orientation at the first vertex */ 65099836e0eSStefano Zampini ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 65199836e0eSStefano Zampini ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 65299836e0eSStefano Zampini if (coneSize != faceSizeH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of face vertices %D for face %D should be %D", coneSize, f, faceSizeH); 653e9fa77a1SMatthew G. Knepley /* Hybrid faces are stored as tensor products of edges, so to compare them to normal faces, we have to flip */ 654e9fa77a1SMatthew G. Knepley if (faceSize == 4 && c >= pMax && faceSizeAll != faceSizeAllT && f < pEnd[faceDepth] - faceH) {fv[2] = 3; fv[3] = 2;} 65599836e0eSStefano Zampini /* - First find the initial vertex */ 656e9fa77a1SMatthew G. Knepley for (i = 0; i < faceSizeH; ++i) if (cellFace[fv[0]] == cone[i]) break; 65799836e0eSStefano Zampini /* - Try forward comparison */ 658e9fa77a1SMatthew G. Knepley for (j = 0; j < faceSizeH; ++j) if (cellFace[fv[j]] != cone[(i+j)%faceSizeH]) break; 65999836e0eSStefano Zampini if (j == faceSizeH) { 66099836e0eSStefano Zampini if ((faceSizeH == 2) && (i == 1)) ornt = -2; 66199836e0eSStefano Zampini else ornt = i; 66299836e0eSStefano Zampini } else { 66399836e0eSStefano Zampini /* - Try backward comparison */ 664e9fa77a1SMatthew G. Knepley for (j = 0; j < faceSizeH; ++j) if (cellFace[fv[j]] != cone[(i+faceSizeH-j)%faceSizeH]) break; 66599836e0eSStefano Zampini if (j == faceSizeH) { 66699836e0eSStefano Zampini if (i == 0) ornt = -faceSizeH; 66799836e0eSStefano Zampini else ornt = -i; 668e9fa77a1SMatthew G. Knepley } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine orientation of face %D in cell %D", f, c); 66999836e0eSStefano Zampini } 67099836e0eSStefano Zampini ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 67199836e0eSStefano Zampini } 67299836e0eSStefano Zampini } 67399836e0eSStefano Zampini ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr); 67499836e0eSStefano Zampini } 67599836e0eSStefano Zampini if (face != pEnd[faceDepth]) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of faces %D should be %D", face-pStart[faceDepth], pEnd[faceDepth]-pStart[faceDepth]); 676c907b753SJed Brown ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 6779a5b13a2SMatthew G Knepley ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 6786551a8c7SMatthew G. Knepley ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 67999836e0eSStefano Zampini ierr = DMPlexSetHybridBounds(idm, cMax, fMax, eMax, vMax);CHKERRQ(ierr); 6809a5b13a2SMatthew G Knepley ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr); 6819a5b13a2SMatthew G Knepley ierr = DMPlexStratify(idm);CHKERRQ(ierr); 6829f074e33SMatthew G Knepley PetscFunctionReturn(0); 6839f074e33SMatthew G Knepley } 6849f074e33SMatthew G Knepley 685fa01033eSVaclav Hapla PetscErrorCode DMPlexOrientCell(DM dm, PetscInt p, PetscInt masterConeSize, const PetscInt masterCone[]) 686fa01033eSVaclav Hapla { 68767aa0e3dSVaclav Hapla PetscInt coneSize; 68867aa0e3dSVaclav Hapla PetscInt start1=0; 68967aa0e3dSVaclav Hapla PetscBool reverse1=PETSC_FALSE; 690a0d42dbcSVaclav Hapla const PetscInt *cone=NULL; 691b8e38d09SVaclav Hapla PetscErrorCode ierr; 692b8e38d09SVaclav Hapla 693b8e38d09SVaclav Hapla PetscFunctionBegin; 694b8e38d09SVaclav Hapla ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 695b8e38d09SVaclav Hapla if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */ 696b8e38d09SVaclav Hapla ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 697b8e38d09SVaclav Hapla ierr = DMPlexFixFaceOrientations_Orient_Private(coneSize, masterConeSize, masterCone, cone, &start1, &reverse1);CHKERRQ(ierr); 698b8e38d09SVaclav Hapla #if defined(PETSC_USE_DEBUG) 699b8e38d09SVaclav Hapla if (PetscUnlikely(cone[start1] != masterCone[0])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "The algorithm above is wrong as cone[%d] = %d != %d = masterCone[0]", start1, cone[start1], masterCone[0]); 700b8e38d09SVaclav Hapla #endif 701b8e38d09SVaclav Hapla ierr = DMPlexOrientCell_Internal(dm, p, start1, reverse1);CHKERRQ(ierr); 702b8e38d09SVaclav Hapla #if defined(PETSC_USE_DEBUG) 70367aa0e3dSVaclav Hapla { 70467aa0e3dSVaclav Hapla PetscInt c; 705b8e38d09SVaclav Hapla ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 706b8e38d09SVaclav Hapla for (c = 0; c < 2; c++) { 707b8e38d09SVaclav Hapla if (PetscUnlikely(cone[c] != masterCone[c])) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "The algorithm above is wrong as cone[%d] = %d != %d = masterCone[%d]", c, cone[c], masterCone[c], c); 708b8e38d09SVaclav Hapla } 70967aa0e3dSVaclav Hapla } 710b8e38d09SVaclav Hapla #endif 711b8e38d09SVaclav Hapla PetscFunctionReturn(0); 712b8e38d09SVaclav Hapla } 713b8e38d09SVaclav Hapla 7140fe46fbbSVaclav Hapla PetscErrorCode DMPlexOrientCell_Internal(DM dm, PetscInt p, PetscInt start1, PetscBool reverse1) 715b8e38d09SVaclav Hapla { 716fa01033eSVaclav Hapla PetscInt i, j, k, maxConeSize, coneSize, coneConeSize, supportSize, supportConeSize; 717b8e38d09SVaclav Hapla PetscInt start0, start; 718b8e38d09SVaclav Hapla PetscBool reverse0, reverse; 719fa01033eSVaclav Hapla PetscInt newornt; 720a0d42dbcSVaclav Hapla const PetscInt *cone=NULL, *support=NULL, *supportCone=NULL, *ornts=NULL; 721a0d42dbcSVaclav Hapla PetscInt *newcone=NULL, *newornts=NULL; 722fa01033eSVaclav Hapla PetscErrorCode ierr; 723fa01033eSVaclav Hapla 724fa01033eSVaclav Hapla PetscFunctionBegin; 725b8e38d09SVaclav Hapla if (!start1 && !reverse1) PetscFunctionReturn(0); 726fa01033eSVaclav Hapla ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 727fa01033eSVaclav Hapla if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */ 728fa01033eSVaclav Hapla ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 729b8e38d09SVaclav Hapla ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 730fa01033eSVaclav Hapla /* permute p's cone and orientations */ 731fa01033eSVaclav Hapla ierr = DMPlexGetConeOrientation(dm, p, &ornts);CHKERRQ(ierr); 732fa01033eSVaclav Hapla ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr); 733fa01033eSVaclav Hapla ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr); 734fa01033eSVaclav Hapla ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, cone, start1, reverse1, newcone);CHKERRQ(ierr); 735fa01033eSVaclav Hapla ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, ornts, start1, reverse1, newornts);CHKERRQ(ierr); 736fa01033eSVaclav Hapla /* if direction of p (face) is flipped, flip also p's cone points (edges) */ 737fa01033eSVaclav Hapla if (reverse1) { 738fa01033eSVaclav Hapla for (i=0; i<coneSize; i++) { 739fa01033eSVaclav Hapla ierr = DMPlexGetConeSize(dm, cone[i], &coneConeSize);CHKERRQ(ierr); 740fa01033eSVaclav Hapla ierr = DMPlexFixFaceOrientations_Translate_Private(newornts[i], &start0, &reverse0);CHKERRQ(ierr); 741fa01033eSVaclav Hapla ierr = DMPlexFixFaceOrientations_Combine_Private(coneConeSize, start0, reverse0, 1, PETSC_FALSE, &start, &reverse);CHKERRQ(ierr); 742fa01033eSVaclav Hapla ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneConeSize, start, reverse, &newornts[i]);CHKERRQ(ierr); 743fa01033eSVaclav Hapla } 744fa01033eSVaclav Hapla } 745fa01033eSVaclav Hapla ierr = DMPlexSetConeOrientation(dm, p, newornts);CHKERRQ(ierr); 746fa01033eSVaclav Hapla /* fix oriention of p within cones of p's support points */ 747fa01033eSVaclav Hapla ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr); 748fa01033eSVaclav Hapla ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr); 749fa01033eSVaclav Hapla for (j=0; j<supportSize; j++) { 750fa01033eSVaclav Hapla ierr = DMPlexGetCone(dm, support[j], &supportCone);CHKERRQ(ierr); 751fa01033eSVaclav Hapla ierr = DMPlexGetConeSize(dm, support[j], &supportConeSize);CHKERRQ(ierr); 752fa01033eSVaclav Hapla ierr = DMPlexGetConeOrientation(dm, support[j], &ornts);CHKERRQ(ierr); 753fa01033eSVaclav Hapla for (k=0; k<supportConeSize; k++) { 754fa01033eSVaclav Hapla if (supportCone[k] != p) continue; 755fa01033eSVaclav Hapla ierr = DMPlexFixFaceOrientations_Translate_Private(ornts[k], &start0, &reverse0);CHKERRQ(ierr); 756fa01033eSVaclav Hapla ierr = DMPlexFixFaceOrientations_Combine_Private(coneSize, start0, reverse0, start1, reverse1, &start, &reverse);CHKERRQ(ierr); 757fa01033eSVaclav Hapla ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneSize, start, reverse, &newornt);CHKERRQ(ierr); 758fa01033eSVaclav Hapla ierr = DMPlexInsertConeOrientation(dm, support[j], k, newornt);CHKERRQ(ierr); 759fa01033eSVaclav Hapla } 760fa01033eSVaclav Hapla } 761fa01033eSVaclav Hapla /* rewrite cone */ 762fa01033eSVaclav Hapla ierr = DMPlexSetCone(dm, p, newcone);CHKERRQ(ierr); 763fa01033eSVaclav Hapla ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr); 764fa01033eSVaclav Hapla ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr); 765fa01033eSVaclav Hapla PetscFunctionReturn(0); 766fa01033eSVaclav Hapla } 767fa01033eSVaclav Hapla 768f80536cbSVaclav Hapla static PetscErrorCode SortRmineRremoteByRemote_Private(PetscSF sf, PetscInt *rmine1[], PetscInt *rremote1[]) 769f80536cbSVaclav Hapla { 770f80536cbSVaclav Hapla PetscInt nleaves; 771f80536cbSVaclav Hapla PetscInt nranks; 772a0d42dbcSVaclav Hapla const PetscMPIInt *ranks=NULL; 773a0d42dbcSVaclav Hapla const PetscInt *roffset=NULL, *rmine=NULL, *rremote=NULL; 774f80536cbSVaclav Hapla PetscInt n, o, r; 775f80536cbSVaclav Hapla PetscErrorCode ierr; 776f80536cbSVaclav Hapla 777f80536cbSVaclav Hapla PetscFunctionBegin; 778f80536cbSVaclav Hapla ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote);CHKERRQ(ierr); 779f80536cbSVaclav Hapla nleaves = roffset[nranks]; 780f80536cbSVaclav Hapla ierr = PetscMalloc2(nleaves, rmine1, nleaves, rremote1);CHKERRQ(ierr); 781f80536cbSVaclav Hapla for (r=0; r<nranks; r++) { 782f80536cbSVaclav Hapla /* simultaneously sort rank-wise portions of rmine & rremote by values in rremote 783f80536cbSVaclav Hapla - to unify order with the other side */ 784f80536cbSVaclav Hapla o = roffset[r]; 785f80536cbSVaclav Hapla n = roffset[r+1] - o; 786f80536cbSVaclav Hapla ierr = PetscMemcpy(&(*rmine1)[o], &rmine[o], n*sizeof(PetscInt));CHKERRQ(ierr); 787f80536cbSVaclav Hapla ierr = PetscMemcpy(&(*rremote1)[o], &rremote[o], n*sizeof(PetscInt));CHKERRQ(ierr); 788f80536cbSVaclav Hapla ierr = PetscSortIntWithArray(n, &(*rremote1)[o], &(*rmine1)[o]);CHKERRQ(ierr); 789f80536cbSVaclav Hapla } 790f80536cbSVaclav Hapla PetscFunctionReturn(0); 791f80536cbSVaclav Hapla } 792f80536cbSVaclav Hapla 793bba29ab8SVaclav Hapla PetscErrorCode DMPlexOrientInterface(DM dm) 7942e745bdaSMatthew G. Knepley { 795a0d42dbcSVaclav Hapla PetscSF sf=NULL; 796cae7fe92SVaclav Hapla PetscInt (*roots)[2], (*leaves)[2]; 7978a650c75SVaclav Hapla PetscMPIInt (*rootsRanks)[2], (*leavesRanks)[2]; 798a0d42dbcSVaclav Hapla const PetscInt *locals=NULL; 799a0d42dbcSVaclav Hapla const PetscSFNode *remotes=NULL; 8008a650c75SVaclav Hapla PetscInt nroots, nleaves, p, c; 801f80536cbSVaclav Hapla PetscInt nranks, n, o, r; 802a0d42dbcSVaclav Hapla const PetscMPIInt *ranks=NULL; 803a0d42dbcSVaclav Hapla const PetscInt *roffset=NULL; 804a0d42dbcSVaclav Hapla PetscInt *rmine1=NULL, *rremote1=NULL; /* rmine and rremote copies simultaneously sorted by rank and rremote */ 805a0d42dbcSVaclav Hapla const PetscInt *cone=NULL; 806adeface4SVaclav Hapla PetscInt coneSize, ind0; 8072e745bdaSMatthew G. Knepley MPI_Comm comm; 8082e745bdaSMatthew G. Knepley PetscMPIInt rank; 8092e745bdaSMatthew G. Knepley PetscInt debug = 0; 8102e745bdaSMatthew G. Knepley PetscErrorCode ierr; 8112e745bdaSMatthew G. Knepley 8122e745bdaSMatthew G. Knepley PetscFunctionBegin; 813df6a9fadSVaclav Hapla ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 8142e745bdaSMatthew G. Knepley ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &locals, &remotes);CHKERRQ(ierr); 8153ede9f65SMatthew G. Knepley if (nroots < 0) PetscFunctionReturn(0); 816f80536cbSVaclav Hapla ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 817f80536cbSVaclav Hapla ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, NULL, NULL);CHKERRQ(ierr); 818dc21a0bfSVaclav Hapla #if defined(PETSC_USE_DEBUG) 819dc21a0bfSVaclav Hapla ierr = DMViewFromOptions(dm, NULL, "-before_fix_dm_view");CHKERRQ(ierr); 820dc21a0bfSVaclav Hapla ierr = DMPlexCheckPointSF(dm);CHKERRQ(ierr); 821dc21a0bfSVaclav Hapla #endif 822f80536cbSVaclav Hapla ierr = SortRmineRremoteByRemote_Private(sf, &rmine1, &rremote1);CHKERRQ(ierr); 8238a650c75SVaclav Hapla ierr = PetscMalloc4(nroots, &roots, nroots, &leaves, nroots, &rootsRanks, nroots, &leavesRanks);CHKERRQ(ierr); 8242e745bdaSMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 8252e745bdaSMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 8262e745bdaSMatthew G. Knepley if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Roots\n");CHKERRQ(ierr);} 8279e24d8a0SVaclav Hapla for (p = 0; p < nroots; ++p) { 8289e24d8a0SVaclav Hapla ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 8299e24d8a0SVaclav Hapla ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 8302e745bdaSMatthew G. Knepley /* Translate all points to root numbering */ 8318a650c75SVaclav Hapla for (c = 0; c < 2; c++) { 8328a650c75SVaclav Hapla if (coneSize > 1) { 8338a650c75SVaclav Hapla ierr = PetscFindInt(cone[c], nleaves, locals, &ind0);CHKERRQ(ierr); 8348a650c75SVaclav Hapla if (ind0 < 0) { 8358a650c75SVaclav Hapla roots[p][c] = cone[c]; 8368a650c75SVaclav Hapla rootsRanks[p][c] = rank; 837c8148b97SVaclav Hapla } else { 8388a650c75SVaclav Hapla roots[p][c] = remotes[ind0].index; 8398a650c75SVaclav Hapla rootsRanks[p][c] = remotes[ind0].rank; 8408a650c75SVaclav Hapla } 8418a650c75SVaclav Hapla } else { 8428a650c75SVaclav Hapla roots[p][c] = -1; 8438a650c75SVaclav Hapla rootsRanks[p][c] = -1; 8442e745bdaSMatthew G. Knepley } 8452e745bdaSMatthew G. Knepley } 8468a650c75SVaclav Hapla } 8478a650c75SVaclav Hapla if (debug) { 8488a650c75SVaclav Hapla for (p = 0; p < nroots; ++p) { 8498a650c75SVaclav Hapla ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 8508a650c75SVaclav Hapla ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 8518a650c75SVaclav Hapla if (coneSize > 1) { 8528a650c75SVaclav Hapla ierr = PetscSynchronizedPrintf(comm, "[%d] %D: cone=[%D %D] roots=[%D %D] rootsRanks=[%D %D]\n", rank, p, cone[0], cone[1], roots[p][0], roots[p][1], rootsRanks[p][0], rootsRanks[p][1]);CHKERRQ(ierr); 8538a650c75SVaclav Hapla } 8548a650c75SVaclav Hapla } 8558a650c75SVaclav Hapla ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 8568a650c75SVaclav Hapla } 8579e24d8a0SVaclav Hapla for (p = 0; p < nroots; ++p) { 8588ccb905dSVaclav Hapla for (c = 0; c < 2; c++) { 8598ccb905dSVaclav Hapla leaves[p][c] = -2; 8608ccb905dSVaclav Hapla leavesRanks[p][c] = -2; 8618ccb905dSVaclav Hapla } 862c8148b97SVaclav Hapla } 8632e745bdaSMatthew G. Knepley ierr = PetscSFBcastBegin(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr); 8648a650c75SVaclav Hapla ierr = PetscSFBcastBegin(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr); 8652e745bdaSMatthew G. Knepley ierr = PetscSFBcastEnd(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr); 8668a650c75SVaclav Hapla ierr = PetscSFBcastEnd(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr); 867c8148b97SVaclav Hapla if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);} 868c8148b97SVaclav Hapla if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Referred leaves\n");CHKERRQ(ierr);} 8699e24d8a0SVaclav Hapla for (p = 0; p < nroots; ++p) { 8709e24d8a0SVaclav Hapla if (leaves[p][0] < 0) continue; 8719e24d8a0SVaclav Hapla ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 8729e24d8a0SVaclav Hapla ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 873ea211068SVaclav Hapla if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] %D: cone=[%D %D] leaves=[%D %D] roots=[%D %D] leavesRanks=[%D %D] rootsRanks=[%D %D]\n", rank, p, cone[0], cone[1], leaves[p][0], leaves[p][1], roots[p][0], roots[p][1], leavesRanks[p][0], leavesRanks[p][1], rootsRanks[p][0], rootsRanks[p][1]);CHKERRQ(ierr);} 874ea211068SVaclav Hapla if ((leaves[p][0] != roots[p][0]) || (leaves[p][1] != roots[p][1]) || (leavesRanks[p][0] != rootsRanks[p][0]) || (leavesRanks[p][0] != rootsRanks[p][0])) { 875f80536cbSVaclav Hapla PetscInt masterCone[2]; 876f80536cbSVaclav Hapla /* Translate these two cone points back to leave numbering */ 877f80536cbSVaclav Hapla for (c = 0; c < 2; c++) { 878def636c8SVaclav Hapla if (leavesRanks[p][c] == rank) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "this should never happen - remote rank of point %D is the same rank",leavesRanks[p][c]); 879f80536cbSVaclav Hapla /* Find index of rank leavesRanks[p][c] among remote ranks */ 880f80536cbSVaclav Hapla /* No need for PetscMPIIntCast because these integers were originally cast from PetscMPIInt. */ 881f80536cbSVaclav Hapla ierr = PetscFindMPIInt((PetscMPIInt)leavesRanks[p][c], nranks, ranks, &r);CHKERRQ(ierr); 882f80536cbSVaclav Hapla if (PetscUnlikely(r < 0)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "this should never happen - rank %D not found among remote ranks",leavesRanks[p][c]); 883f80536cbSVaclav Hapla /* Find point leaves[p][c] among remote points aimed at rank leavesRanks[p][c] */ 884f80536cbSVaclav Hapla o = roffset[r]; 885f80536cbSVaclav Hapla n = roffset[r+1] - o; 886f80536cbSVaclav Hapla ierr = PetscFindInt(leaves[p][c], n, &rremote1[o], &ind0);CHKERRQ(ierr); 887adeface4SVaclav Hapla if (PetscUnlikely(ind0 < 0)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "No cone point of %D is connected to (%D, %D) - it seems there is missing connection in point SF!",p,ranks[r],leaves[p][c]); 888f80536cbSVaclav Hapla /* Get the corresponding local point */ 889f80536cbSVaclav Hapla masterCone[c] = rmine1[o+ind0];CHKERRQ(ierr); 890f80536cbSVaclav Hapla } 891f80536cbSVaclav Hapla if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] %D: masterCone=[%D %D]\n", rank, p, masterCone[0], masterCone[1]);CHKERRQ(ierr);} 892f80536cbSVaclav Hapla /* Vaclav's note: Here we only compare first 2 points of the cone. Full cone size would lead to stronger self-checking. */ 893f80536cbSVaclav Hapla ierr = DMPlexOrientCell(dm, p, 2, masterCone);CHKERRQ(ierr); 89423aaf07bSVaclav Hapla } 8952e745bdaSMatthew G. Knepley } 896adeface4SVaclav Hapla #if defined(PETSC_USE_DEBUG) 897adeface4SVaclav Hapla ierr = DMViewFromOptions(dm, NULL, "-after_fix_dm_view");CHKERRQ(ierr); 898adeface4SVaclav Hapla for (r = 0; r < nleaves; ++r) { 899adeface4SVaclav Hapla p = locals[r]; 900adeface4SVaclav Hapla ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 901adeface4SVaclav Hapla if (!coneSize) continue; 902adeface4SVaclav Hapla ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 903adeface4SVaclav Hapla for (c = 0; c < 2; c++) { 904adeface4SVaclav Hapla ierr = PetscFindInt(cone[c], nleaves, locals, &ind0);CHKERRQ(ierr); 905adeface4SVaclav Hapla if (ind0 < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point SF contains %D but is missing its cone point cone[%D] = %D!", p, c, cone[c]); 906adeface4SVaclav Hapla if (leaves[p][c] != remotes[ind0].index || leavesRanks[p][c] != remotes[ind0].rank) { 907adeface4SVaclav Hapla if (leavesRanks[p][c] == rank) { 908adeface4SVaclav Hapla PetscInt ind1; 909adeface4SVaclav Hapla ierr = PetscFindInt(leaves[p][c], nleaves, locals, &ind1);CHKERRQ(ierr); 910adeface4SVaclav Hapla if (ind1 < 0) { 911adeface4SVaclav Hapla SETERRQ8(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D = locals[%d]: cone[%D]=%D --> (%D, %D) differs from the enforced (%D, %D). The latter was not even found among the local SF points - it is probably broken!", p, r, c, cone[c], remotes[ind0].rank, remotes[ind0].index, leavesRanks[p][c], leaves[p][c]); 912adeface4SVaclav Hapla } else { 913adeface4SVaclav Hapla SETERRQ9(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D = locals[%d]: cone[%D]=%D --> (%D, %D) differs from the enforced %D --> (%D, %D). Is the algorithm above or the point SF broken?", p, r, c, cone[c], remotes[ind0].rank, remotes[ind0].index, leaves[p][c], remotes[ind1].rank, remotes[ind1].index); 914adeface4SVaclav Hapla } 915adeface4SVaclav Hapla } else { 916adeface4SVaclav Hapla SETERRQ8(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D = locals[%d]: cone[%D]=%D --> (%D, %D) differs from the enforced (%D, %D). Is the algorithm above or the point SF broken?", p, r, c, cone[c], remotes[ind0].rank, remotes[ind0].index, leavesRanks[p][c], leaves[p][c]); 917adeface4SVaclav Hapla } 918adeface4SVaclav Hapla } 919adeface4SVaclav Hapla } 920adeface4SVaclav Hapla } 921adeface4SVaclav Hapla #endif 9222e745bdaSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);} 9238a650c75SVaclav Hapla ierr = PetscFree4(roots, leaves, rootsRanks, leavesRanks);CHKERRQ(ierr); 9247c7bb832SVaclav Hapla ierr = PetscFree2(rmine1, rremote1);CHKERRQ(ierr); 9252e745bdaSMatthew G. Knepley PetscFunctionReturn(0); 9262e745bdaSMatthew G. Knepley } 9272e745bdaSMatthew G. Knepley 9282e72742cSMatthew 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[]) 9297bffde78SMichael Lange { 9302e72742cSMatthew G. Knepley PetscInt idx; 9312e72742cSMatthew G. Knepley PetscMPIInt rank; 9322e72742cSMatthew G. Knepley PetscBool flg; 9337bffde78SMichael Lange PetscErrorCode ierr; 9347bffde78SMichael Lange 9357bffde78SMichael Lange PetscFunctionBegin; 9362e72742cSMatthew G. Knepley ierr = PetscOptionsHasName(NULL, NULL, opt, &flg);CHKERRQ(ierr); 9372e72742cSMatthew G. Knepley if (!flg) PetscFunctionReturn(0); 9382e72742cSMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 9392e72742cSMatthew G. Knepley ierr = PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);CHKERRQ(ierr); 9402e72742cSMatthew 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);} 9412e72742cSMatthew G. Knepley ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 9422e72742cSMatthew G. Knepley PetscFunctionReturn(0); 9432e72742cSMatthew G. Knepley } 9442e72742cSMatthew G. Knepley 9452e72742cSMatthew G. Knepley static PetscErrorCode SFNodeArrayViewFromOptions(MPI_Comm comm, const char opt[], const char name[], const char idxname[], PetscInt n, const PetscSFNode a[]) 9462e72742cSMatthew G. Knepley { 9472e72742cSMatthew G. Knepley PetscInt idx; 9482e72742cSMatthew G. Knepley PetscMPIInt rank; 9492e72742cSMatthew G. Knepley PetscBool flg; 9502e72742cSMatthew G. Knepley PetscErrorCode ierr; 9512e72742cSMatthew G. Knepley 9522e72742cSMatthew G. Knepley PetscFunctionBegin; 9532e72742cSMatthew G. Knepley ierr = PetscOptionsHasName(NULL, NULL, opt, &flg);CHKERRQ(ierr); 9542e72742cSMatthew G. Knepley if (!flg) PetscFunctionReturn(0); 9552e72742cSMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 9562e72742cSMatthew G. Knepley ierr = PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);CHKERRQ(ierr); 9572e72742cSMatthew G. Knepley if (idxname) { 9582e72742cSMatthew 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);} 9592e72742cSMatthew G. Knepley } else { 9602e72742cSMatthew 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);} 9612e72742cSMatthew G. Knepley } 9622e72742cSMatthew G. Knepley ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 9632e72742cSMatthew G. Knepley PetscFunctionReturn(0); 9642e72742cSMatthew G. Knepley } 9652e72742cSMatthew G. Knepley 9662e72742cSMatthew G. Knepley static PetscErrorCode DMPlexMapToLocalPoint(PetscHMapIJ roothash, const PetscInt localPoints[], PetscMPIInt rank, PetscSFNode remotePoint, PetscInt *localPoint) 9672e72742cSMatthew G. Knepley { 9682e72742cSMatthew G. Knepley PetscErrorCode ierr; 9692e72742cSMatthew G. Knepley 9702e72742cSMatthew G. Knepley PetscFunctionBegin; 9712e72742cSMatthew G. Knepley if (remotePoint.rank == rank) { 9722e72742cSMatthew G. Knepley *localPoint = remotePoint.index; 9732e72742cSMatthew G. Knepley } else { 9742e72742cSMatthew G. Knepley PetscHashIJKey key; 9752e72742cSMatthew G. Knepley PetscInt root; 9762e72742cSMatthew G. Knepley 9772e72742cSMatthew G. Knepley key.i = remotePoint.index; 9782e72742cSMatthew G. Knepley key.j = remotePoint.rank; 9792e72742cSMatthew G. Knepley ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr); 9802e72742cSMatthew G. Knepley if (root >= 0) { 9812e72742cSMatthew G. Knepley *localPoint = localPoints[root]; 9822e72742cSMatthew G. Knepley } else PetscFunctionReturn(1); 9832e72742cSMatthew G. Knepley } 9842e72742cSMatthew G. Knepley PetscFunctionReturn(0); 9852e72742cSMatthew G. Knepley } 9862e72742cSMatthew G. Knepley 9872e72742cSMatthew G. Knepley /*@ 9882e72742cSMatthew G. Knepley DMPlexInterpolatePointSF - Insert interpolated points in the overlap into the PointSF in parallel, following local interpolation 9892e72742cSMatthew G. Knepley 9902e72742cSMatthew G. Knepley Collective on DM 9912e72742cSMatthew G. Knepley 9922e72742cSMatthew G. Knepley Input Parameters: 9932e72742cSMatthew G. Knepley + dm - The interpolated DM 9942e72742cSMatthew G. Knepley - pointSF - The initial SF without interpolated points 9952e72742cSMatthew G. Knepley 9962e72742cSMatthew G. Knepley Output Parameter: 9972e72742cSMatthew G. Knepley . pointSF - The SF including interpolated points 9982e72742cSMatthew G. Knepley 9992e72742cSMatthew G. Knepley Level: intermediate 10002e72742cSMatthew G. Knepley 10012e72742cSMatthew 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 10022e72742cSMatthew G. Knepley 10032e72742cSMatthew G. Knepley .keywords: mesh 10042e72742cSMatthew G. Knepley .seealso: DMPlexInterpolate(), DMPlexUninterpolate() 10052e72742cSMatthew G. Knepley @*/ 1006e53487d3SMatthew G. Knepley PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF) 10072e72742cSMatthew G. Knepley { 10082e72742cSMatthew G. Knepley /* 10092e72742cSMatthew G. Knepley Okay, the algorithm is: 10102e72742cSMatthew G. Knepley - Take each point in the overlap (root) 10112e72742cSMatthew G. Knepley - Look at the neighboring points in the overlap (candidates) 10122e72742cSMatthew G. Knepley - Send these candidate points to neighbors 10132e72742cSMatthew G. Knepley - Neighbor checks for edge between root and candidate 10142e72742cSMatthew G. Knepley - If edge is found, it replaces candidate point with edge point 10152e72742cSMatthew G. Knepley - Send back the overwritten candidates (claims) 10162e72742cSMatthew G. Knepley - Original guy checks for edges, different from original candidate, and gets its own edge 10172e72742cSMatthew G. Knepley - This pair is put into SF 10182e72742cSMatthew G. Knepley 10192e72742cSMatthew G. Knepley We need a new algorithm that tolerates groups larger than 2. 10202e72742cSMatthew G. Knepley - Take each point in the overlap (root) 10212e72742cSMatthew G. Knepley - Find all collections of points in the overlap which make faces (do early join) 10222e72742cSMatthew G. Knepley - Send collections as candidates (add size as first number) 102366aa2a39SMatthew G. Knepley - Make sure to send collection to all owners of all overlap points in collection 10242e72742cSMatthew G. Knepley - Neighbor check for face in collections 10252e72742cSMatthew G. Knepley - If face is found, it replaces candidate point with face point 10262e72742cSMatthew G. Knepley - Send back the overwritten candidates (claims) 10272e72742cSMatthew G. Knepley - Original guy checks for faces, different from original candidate, and gets its own face 10282e72742cSMatthew G. Knepley - This pair is put into SF 10292e72742cSMatthew G. Knepley */ 10302e72742cSMatthew G. Knepley PetscHMapI leafhash; 10312e72742cSMatthew G. Knepley PetscHMapIJ roothash; 10322e72742cSMatthew G. Knepley const PetscInt *localPoints, *rootdegree; 10332e72742cSMatthew G. Knepley const PetscSFNode *remotePoints; 10342e72742cSMatthew G. Knepley PetscSFNode *candidates, *candidatesRemote, *claims; 10352e72742cSMatthew G. Knepley PetscSection candidateSection, candidateSectionRemote, claimSection; 10362e72742cSMatthew G. Knepley PetscInt numLeaves, l, numRoots, r, candidatesSize, candidatesRemoteSize; 10372e72742cSMatthew G. Knepley PetscMPIInt size, rank; 10382e72742cSMatthew G. Knepley PetscHashIJKey key; 10392e72742cSMatthew G. Knepley PetscBool debug = PETSC_FALSE; 10402e72742cSMatthew G. Knepley PetscErrorCode ierr; 10412e72742cSMatthew G. Knepley 10422e72742cSMatthew G. Knepley PetscFunctionBegin; 10432e72742cSMatthew G. Knepley ierr = PetscOptionsHasName(NULL, ((PetscObject) dm)->prefix, "-dmplex_interp_debug", &debug);CHKERRQ(ierr); 10449852e123SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr); 10457bffde78SMichael Lange ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 10467bffde78SMichael Lange ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 10479852e123SBarry Smith if (size < 2 || numRoots < 0) PetscFunctionReturn(0); 10482e72742cSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_interp_pre_view");CHKERRQ(ierr); 10492e72742cSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) pointSF, NULL, "-petscsf_interp_pre_view");CHKERRQ(ierr); 105025afeb17SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 10517bffde78SMichael Lange /* Build hashes of points in the SF for efficient lookup */ 1052e8f14785SLisandro Dalcin ierr = PetscHMapICreate(&leafhash);CHKERRQ(ierr); 1053e8f14785SLisandro Dalcin ierr = PetscHMapIJCreate(&roothash);CHKERRQ(ierr); 10542e72742cSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) { 10552e72742cSMatthew G. Knepley key.i = remotePoints[l].index; 10562e72742cSMatthew G. Knepley key.j = remotePoints[l].rank; 10572e72742cSMatthew G. Knepley ierr = PetscHMapISet(leafhash, localPoints[l], l);CHKERRQ(ierr); 10582e72742cSMatthew G. Knepley ierr = PetscHMapIJSet(roothash, key, l);CHKERRQ(ierr); 10597bffde78SMichael Lange } 106066aa2a39SMatthew G. Knepley /* Compute root degree to identify shared points */ 10612e72742cSMatthew G. Knepley ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr); 10622e72742cSMatthew G. Knepley ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr); 10632e72742cSMatthew G. Knepley ierr = IntArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-interp_root_degree_view", "Root degree", "point", "degree", numRoots, rootdegree);CHKERRQ(ierr); 106466aa2a39SMatthew G. Knepley /* Build a section / SFNode array of candidate points (face bd points) in the cone(support(leaf)), 106566aa2a39SMatthew G. Knepley where each candidate is defined by a set of remote points (roots) for the other points that define the face. */ 10667bffde78SMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr); 10677bffde78SMichael Lange ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr); 10687bffde78SMichael Lange { 106966aa2a39SMatthew G. Knepley PetscHMapIJ facehash; 10702e72742cSMatthew G. Knepley 107166aa2a39SMatthew G. Knepley ierr = PetscHMapIJCreate(&facehash);CHKERRQ(ierr); 10722e72742cSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) { 10732e72742cSMatthew G. Knepley const PetscInt localPoint = localPoints[l]; 10742e72742cSMatthew G. Knepley const PetscInt *support; 10752e72742cSMatthew G. Knepley PetscInt supportSize, s; 10762e72742cSMatthew G. Knepley 10772e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Checking local point %D\n", rank, localPoint);CHKERRQ(ierr);} 10782e72742cSMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, localPoint, &supportSize);CHKERRQ(ierr); 10792e72742cSMatthew G. Knepley ierr = DMPlexGetSupport(dm, localPoint, &support);CHKERRQ(ierr); 10802e72742cSMatthew G. Knepley for (s = 0; s < supportSize; ++s) { 10812e72742cSMatthew G. Knepley const PetscInt face = support[s]; 10822e72742cSMatthew G. Knepley const PetscInt *cone; 10832e72742cSMatthew G. Knepley PetscInt coneSize, c, f, root; 10842e72742cSMatthew G. Knepley PetscBool isFace = PETSC_TRUE; 10852e72742cSMatthew G. Knepley 10862e72742cSMatthew G. Knepley /* Only add face once */ 10872e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Support point %D\n", rank, face);CHKERRQ(ierr);} 108866aa2a39SMatthew G. Knepley key.i = localPoint; 108966aa2a39SMatthew G. Knepley key.j = face; 109066aa2a39SMatthew G. Knepley ierr = PetscHMapIJGet(facehash, key, &f);CHKERRQ(ierr); 10912e72742cSMatthew G. Knepley if (f >= 0) continue; 10922e72742cSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, face, &coneSize);CHKERRQ(ierr); 10932e72742cSMatthew G. Knepley ierr = DMPlexGetCone(dm, face, &cone);CHKERRQ(ierr); 10942e72742cSMatthew G. Knepley /* If a cone point does not map to leaves on any proc, then do not put face in SF */ 10952e72742cSMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 10962e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Cone point %D\n", rank, cone[c]);CHKERRQ(ierr);} 10972e72742cSMatthew G. Knepley ierr = PetscHMapIGet(leafhash, cone[c], &root);CHKERRQ(ierr); 10982e72742cSMatthew G. Knepley if (!rootdegree[cone[c]] && (root < 0)) {isFace = PETSC_FALSE; break;} 10992e72742cSMatthew G. Knepley } 11002e72742cSMatthew G. Knepley if (isFace) { 110166aa2a39SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Found shared face %D\n", rank, face);CHKERRQ(ierr);} 110266aa2a39SMatthew G. Knepley ierr = PetscHMapIJSet(facehash, key, l);CHKERRQ(ierr); 11032e72742cSMatthew G. Knepley ierr = PetscSectionAddDof(candidateSection, localPoint, coneSize);CHKERRQ(ierr); 11047bffde78SMichael Lange } 11057bffde78SMichael Lange } 11062e72742cSMatthew G. Knepley } 11072e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);} 110866aa2a39SMatthew G. Knepley ierr = PetscHMapIJClear(facehash);CHKERRQ(ierr); 11097bffde78SMichael Lange ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr); 11107bffde78SMichael Lange ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr); 11117bffde78SMichael Lange ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr); 11122e72742cSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) { 11132e72742cSMatthew G. Knepley const PetscInt localPoint = localPoints[l]; 11142e72742cSMatthew G. Knepley const PetscInt *support; 11152e72742cSMatthew G. Knepley PetscInt supportSize, s, offset, idx = 0; 11162e72742cSMatthew G. Knepley 11172e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Checking local point %D\n", rank, localPoint);CHKERRQ(ierr);} 11182e72742cSMatthew G. Knepley ierr = PetscSectionGetOffset(candidateSection, localPoint, &offset);CHKERRQ(ierr); 11192e72742cSMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, localPoint, &supportSize);CHKERRQ(ierr); 11202e72742cSMatthew G. Knepley ierr = DMPlexGetSupport(dm, localPoint, &support);CHKERRQ(ierr); 11212e72742cSMatthew G. Knepley for (s = 0; s < supportSize; ++s) { 11222e72742cSMatthew G. Knepley const PetscInt face = support[s]; 11232e72742cSMatthew G. Knepley const PetscInt *cone; 11242e72742cSMatthew G. Knepley PetscInt coneSize, c, f, root; 11252e72742cSMatthew G. Knepley PetscBool isFace = PETSC_TRUE; 11262e72742cSMatthew G. Knepley 11272e72742cSMatthew G. Knepley /* Only add face once */ 11282e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Support point %D\n", rank, face);CHKERRQ(ierr);} 112966aa2a39SMatthew G. Knepley key.i = localPoint; 113066aa2a39SMatthew G. Knepley key.j = face; 113166aa2a39SMatthew G. Knepley ierr = PetscHMapIJGet(facehash, key, &f);CHKERRQ(ierr); 11322e72742cSMatthew G. Knepley if (f >= 0) continue; 11332e72742cSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, face, &coneSize);CHKERRQ(ierr); 11342e72742cSMatthew G. Knepley ierr = DMPlexGetCone(dm, face, &cone);CHKERRQ(ierr); 11352e72742cSMatthew G. Knepley /* If a cone point does not map to leaves on any proc, then do not put face in SF */ 11362e72742cSMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 11372e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Cone point %D\n", rank, cone[c]);CHKERRQ(ierr);} 11382e72742cSMatthew G. Knepley ierr = PetscHMapIGet(leafhash, cone[c], &root);CHKERRQ(ierr); 11392e72742cSMatthew G. Knepley if (!rootdegree[cone[c]] && (root < 0)) {isFace = PETSC_FALSE; break;} 11402e72742cSMatthew G. Knepley } 11412e72742cSMatthew G. Knepley if (isFace) { 114266aa2a39SMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Adding shared face %D at idx %D\n", rank, face, idx);CHKERRQ(ierr);} 114366aa2a39SMatthew G. Knepley ierr = PetscHMapIJSet(facehash, key, l);CHKERRQ(ierr); 11442e72742cSMatthew G. Knepley candidates[offset+idx].rank = -1; 11452e72742cSMatthew G. Knepley candidates[offset+idx++].index = coneSize-1; 11462e72742cSMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 11472e72742cSMatthew G. Knepley if (cone[c] == localPoint) continue; 11482e72742cSMatthew G. Knepley if (rootdegree[cone[c]]) { 11492e72742cSMatthew G. Knepley candidates[offset+idx].rank = rank; 11502e72742cSMatthew G. Knepley candidates[offset+idx++].index = cone[c]; 11512e72742cSMatthew G. Knepley } else { 11522e72742cSMatthew G. Knepley ierr = PetscHMapIGet(leafhash, cone[c], &root);CHKERRQ(ierr); 11532e72742cSMatthew G. Knepley if (root < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot locate local point %D in SF", cone[c]); 11542e72742cSMatthew G. Knepley candidates[offset+idx++] = remotePoints[root]; 11557bffde78SMichael Lange } 11567bffde78SMichael Lange } 11572e72742cSMatthew G. Knepley } 11582e72742cSMatthew G. Knepley } 11592e72742cSMatthew G. Knepley } 11602e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);} 116166aa2a39SMatthew G. Knepley ierr = PetscHMapIJDestroy(&facehash);CHKERRQ(ierr); 11622e72742cSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) candidateSection, NULL, "-petscsection_interp_candidate_view");CHKERRQ(ierr); 11632e72742cSMatthew G. Knepley ierr = SFNodeArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-petscsection_interp_candidate_view", "Candidates", NULL, candidatesSize, candidates);CHKERRQ(ierr); 11647bffde78SMichael Lange } 11657bffde78SMichael Lange /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */ 11662e72742cSMatthew G. Knepley /* Note that this section is indexed by offsets into leaves, not by point number */ 11677bffde78SMichael Lange { 11687bffde78SMichael Lange PetscSF sfMulti, sfInverse, sfCandidates; 11697bffde78SMichael Lange PetscInt *remoteOffsets; 11702e72742cSMatthew G. Knepley 11717bffde78SMichael Lange ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 11727bffde78SMichael Lange ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr); 11737bffde78SMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr); 11747bffde78SMichael Lange ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr); 11757bffde78SMichael Lange ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr); 11767bffde78SMichael Lange ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr); 11777bffde78SMichael Lange ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr); 11787bffde78SMichael Lange ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 11797bffde78SMichael Lange ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 11807bffde78SMichael Lange ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr); 11817bffde78SMichael Lange ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr); 11827bffde78SMichael Lange ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 11832e72742cSMatthew G. Knepley 11842e72742cSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) candidateSectionRemote, NULL, "-petscsection_interp_candidate_remote_view");CHKERRQ(ierr); 11852e72742cSMatthew G. Knepley ierr = SFNodeArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-petscsection_interp_candidate_remote_view", "Remote Candidates", NULL, candidatesRemoteSize, candidatesRemote);CHKERRQ(ierr); 11867bffde78SMichael Lange } 11872e72742cSMatthew G. Knepley /* */ 11887bffde78SMichael Lange { 11892e72742cSMatthew G. Knepley PetscInt idx; 11902e72742cSMatthew G. Knepley /* There is a section point for every leaf attached to a given root point */ 11912e72742cSMatthew G. Knepley for (r = 0, idx = 0; r < numRoots; ++r) { 11922e72742cSMatthew G. Knepley PetscInt deg; 11932e72742cSMatthew G. Knepley for (deg = 0; deg < rootdegree[r]; ++deg, ++idx) { 11942e72742cSMatthew G. Knepley PetscInt offset, dof, d; 11952e72742cSMatthew G. Knepley 11967bffde78SMichael Lange ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr); 11977bffde78SMichael Lange ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr); 11982e72742cSMatthew G. Knepley for (d = 0; d < dof; ++d) { 11992e72742cSMatthew G. Knepley const PetscInt sizeInd = offset+d; 12002e72742cSMatthew G. Knepley const PetscInt numPoints = candidatesRemote[sizeInd].index; 12012e72742cSMatthew G. Knepley const PetscInt *join = NULL; 12022e72742cSMatthew G. Knepley PetscInt points[1024], p, joinSize; 12032e72742cSMatthew G. Knepley 12042e72742cSMatthew G. Knepley points[0] = r; 12052e72742cSMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 12062e72742cSMatthew G. Knepley ierr = DMPlexMapToLocalPoint(roothash, localPoints, rank, candidatesRemote[offset+(++d)], &points[p+1]); 12072e72742cSMatthew G. Knepley if (ierr) {d += numPoints-1 - p; break;} /* We got a point not in our overlap */ 12082e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Checking local candidate %D\n", rank, points[p+1]);CHKERRQ(ierr);} 12097bffde78SMichael Lange } 12102e72742cSMatthew G. Knepley if (ierr) continue; 12112e72742cSMatthew G. Knepley ierr = DMPlexGetJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr); 12127bffde78SMichael Lange if (joinSize == 1) { 12132e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Adding face %D at idx %D\n", rank, join[0], sizeInd);CHKERRQ(ierr);} 12142e72742cSMatthew G. Knepley candidatesRemote[sizeInd].rank = rank; 12152e72742cSMatthew G. Knepley candidatesRemote[sizeInd].index = join[0]; 12167bffde78SMichael Lange } 12172e72742cSMatthew G. Knepley ierr = DMPlexRestoreJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr); 12187bffde78SMichael Lange } 12197bffde78SMichael Lange } 12207bffde78SMichael Lange } 12212e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);} 12227bffde78SMichael Lange } 12237bffde78SMichael Lange /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */ 12247bffde78SMichael Lange { 12257bffde78SMichael Lange PetscSF sfMulti, sfClaims, sfPointNew; 12267bffde78SMichael Lange PetscSFNode *remotePointsNew; 12272e72742cSMatthew G. Knepley PetscHMapI claimshash; 12282e72742cSMatthew G. Knepley PetscInt *remoteOffsets, *localPointsNew; 12292e72742cSMatthew G. Knepley PetscInt claimsSize, pStart, pEnd, root, numLocalNew, p, d; 12302e72742cSMatthew G. Knepley 12317bffde78SMichael Lange ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 12327bffde78SMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr); 12337bffde78SMichael Lange ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr); 12347bffde78SMichael Lange ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr); 12352e72742cSMatthew G. Knepley ierr = PetscSectionGetStorageSize(claimSection, &claimsSize);CHKERRQ(ierr); 12362e72742cSMatthew G. Knepley ierr = PetscMalloc1(claimsSize, &claims);CHKERRQ(ierr); 12377bffde78SMichael Lange ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 12387bffde78SMichael Lange ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 12397bffde78SMichael Lange ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr); 12407bffde78SMichael Lange ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 12412e72742cSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) claimSection, NULL, "-petscsection_interp_claim_view");CHKERRQ(ierr); 12422e72742cSMatthew G. Knepley ierr = SFNodeArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-petscsection_interp_claim_view", "Claims", NULL, claimsSize, claims);CHKERRQ(ierr); 12437bffde78SMichael Lange /* Walk the original section of local supports and add an SF entry for each updated item */ 1244e8f14785SLisandro Dalcin ierr = PetscHMapICreate(&claimshash);CHKERRQ(ierr); 12457bffde78SMichael Lange for (p = 0; p < numRoots; ++p) { 12462e72742cSMatthew G. Knepley PetscInt dof, offset; 12472e72742cSMatthew G. Knepley 12482e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Checking root for claims %D\n", rank, p);CHKERRQ(ierr);} 12497bffde78SMichael Lange ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr); 12507bffde78SMichael Lange ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr); 12512e72742cSMatthew G. Knepley for (d = 0; d < dof;) { 12522e72742cSMatthew G. Knepley if (claims[offset+d].rank >= 0) { 12532e72742cSMatthew G. Knepley const PetscInt faceInd = offset+d; 12542e72742cSMatthew G. Knepley const PetscInt numPoints = candidates[faceInd].index; 12552e72742cSMatthew G. Knepley const PetscInt *join = NULL; 12562e72742cSMatthew G. Knepley PetscInt joinSize, points[1024], c; 12572e72742cSMatthew G. Knepley 12582e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Found claim for remote point (%D, %D)\n", rank, claims[faceInd].rank, claims[faceInd].index);CHKERRQ(ierr);} 12592e72742cSMatthew G. Knepley points[0] = p; 12602e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] point %D\n", rank, points[0]);CHKERRQ(ierr);} 12612e72742cSMatthew G. Knepley for (c = 0, ++d; c < numPoints; ++c, ++d) { 1262e8f14785SLisandro Dalcin key.i = candidates[offset+d].index; 1263e8f14785SLisandro Dalcin key.j = candidates[offset+d].rank; 1264e8f14785SLisandro Dalcin ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr); 12652e72742cSMatthew G. Knepley points[c+1] = localPoints[root]; 12662e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] point %D\n", rank, points[c+1]);CHKERRQ(ierr);} 12672e72742cSMatthew G. Knepley } 12682e72742cSMatthew G. Knepley ierr = DMPlexGetJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr); 12692e72742cSMatthew G. Knepley if (joinSize == 1) { 12702e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Found local face %D\n", rank, join[0]);CHKERRQ(ierr);} 12712e72742cSMatthew G. Knepley ierr = PetscHMapISet(claimshash, join[0], faceInd);CHKERRQ(ierr); 12722e72742cSMatthew G. Knepley } 12732e72742cSMatthew G. Knepley ierr = DMPlexRestoreJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr); 12742e72742cSMatthew G. Knepley } else d += claims[offset+d].index+1; 12757bffde78SMichael Lange } 12767bffde78SMichael Lange } 12772e72742cSMatthew G. Knepley if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);} 12787bffde78SMichael Lange /* Create new pointSF from hashed claims */ 1279e8f14785SLisandro Dalcin ierr = PetscHMapIGetSize(claimshash, &numLocalNew);CHKERRQ(ierr); 12807bffde78SMichael Lange ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 12817bffde78SMichael Lange ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr); 12827bffde78SMichael Lange ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr); 12837bffde78SMichael Lange for (p = 0; p < numLeaves; ++p) { 12847bffde78SMichael Lange localPointsNew[p] = localPoints[p]; 12857bffde78SMichael Lange remotePointsNew[p].index = remotePoints[p].index; 12867bffde78SMichael Lange remotePointsNew[p].rank = remotePoints[p].rank; 12877bffde78SMichael Lange } 1288f3190fc2SToby Isaac p = numLeaves; 1289e8f14785SLisandro Dalcin ierr = PetscHMapIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr); 1290f3190fc2SToby Isaac ierr = PetscSortInt(numLocalNew, &localPointsNew[numLeaves]);CHKERRQ(ierr); 12917bffde78SMichael Lange for (p = numLeaves; p < numLeaves + numLocalNew; ++p) { 12922e72742cSMatthew G. Knepley PetscInt offset; 1293e8f14785SLisandro Dalcin ierr = PetscHMapIGet(claimshash, localPointsNew[p], &offset);CHKERRQ(ierr); 12947bffde78SMichael Lange remotePointsNew[p] = claims[offset]; 12957bffde78SMichael Lange } 12967bffde78SMichael Lange ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr); 12977bffde78SMichael Lange ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); 12987bffde78SMichael Lange ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr); 12997bffde78SMichael Lange ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr); 1300e8f14785SLisandro Dalcin ierr = PetscHMapIDestroy(&claimshash);CHKERRQ(ierr); 13017bffde78SMichael Lange } 1302e8f14785SLisandro Dalcin ierr = PetscHMapIDestroy(&leafhash);CHKERRQ(ierr); 1303e8f14785SLisandro Dalcin ierr = PetscHMapIJDestroy(&roothash);CHKERRQ(ierr); 13047bffde78SMichael Lange ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr); 13057bffde78SMichael Lange ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr); 13067bffde78SMichael Lange ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr); 13077bffde78SMichael Lange ierr = PetscFree(candidates);CHKERRQ(ierr); 13087bffde78SMichael Lange ierr = PetscFree(candidatesRemote);CHKERRQ(ierr); 13097bffde78SMichael Lange ierr = PetscFree(claims);CHKERRQ(ierr); 131025afeb17SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 13117bffde78SMichael Lange PetscFunctionReturn(0); 13127bffde78SMichael Lange } 13137bffde78SMichael Lange 13140c795ddcSMatthew G. Knepley /*@C 131580330477SMatthew G. Knepley DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc. 131680330477SMatthew G. Knepley 131780330477SMatthew G. Knepley Collective on DM 131880330477SMatthew G. Knepley 1319e56d480eSMatthew G. Knepley Input Parameters: 1320e56d480eSMatthew G. Knepley + dm - The DMPlex object with only cells and vertices 132110a67516SMatthew G. Knepley - dmInt - The interpolated DM 132280330477SMatthew G. Knepley 132380330477SMatthew G. Knepley Output Parameter: 13244e3744c5SMatthew G. Knepley . dmInt - The complete DMPlex object 132580330477SMatthew G. Knepley 132680330477SMatthew G. Knepley Level: intermediate 132780330477SMatthew G. Knepley 132895452b02SPatrick Sanan Notes: 132995452b02SPatrick Sanan It does not copy over the coordinates. 133043eeeb2dSStefano Zampini 133180330477SMatthew G. Knepley .keywords: mesh 133243eeeb2dSStefano Zampini .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates() 133380330477SMatthew G. Knepley @*/ 13349f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt) 13359f074e33SMatthew G Knepley { 13369a5b13a2SMatthew G Knepley DM idm, odm = dm; 13377bffde78SMichael Lange PetscSF sfPoint; 13382e1b13c2SMatthew G. Knepley PetscInt depth, dim, d; 133910a67516SMatthew G. Knepley const char *name; 1340b325530aSVaclav Hapla PetscBool flg=PETSC_TRUE; 13419f074e33SMatthew G Knepley PetscErrorCode ierr; 13429f074e33SMatthew G Knepley 13439f074e33SMatthew G Knepley PetscFunctionBegin; 134410a67516SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 134510a67516SMatthew G. Knepley PetscValidPointer(dmInt, 2); 1346a72f3261SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 13472e1b13c2SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1348c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1349b21b8912SMatthew G. Knepley if ((depth == dim) || (dim <= 1)) { 135076b791aaSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 135176b791aaSMatthew G. Knepley idm = dm; 1352b21b8912SMatthew G. Knepley } else { 13539a5b13a2SMatthew G Knepley for (d = 1; d < dim; ++d) { 13549a5b13a2SMatthew G Knepley /* Create interpolated mesh */ 135510a67516SMatthew G. Knepley ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr); 13569a5b13a2SMatthew G Knepley ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr); 1357c73cfb54SMatthew G. Knepley ierr = DMSetDimension(idm, dim);CHKERRQ(ierr); 13587bffde78SMichael Lange if (depth > 0) { 13597bffde78SMichael Lange ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr); 13607bffde78SMichael Lange ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr); 13612e72742cSMatthew G. Knepley ierr = DMPlexInterpolatePointSF(idm, sfPoint);CHKERRQ(ierr); 13627bffde78SMichael Lange } 13639a5b13a2SMatthew G Knepley if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 13649a5b13a2SMatthew G Knepley odm = idm; 13659f074e33SMatthew G Knepley } 136610a67516SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr); 136710a67516SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) idm, name);CHKERRQ(ierr); 136810a67516SMatthew G. Knepley ierr = DMPlexCopyCoordinates(dm, idm);CHKERRQ(ierr); 136910a67516SMatthew G. Knepley ierr = DMCopyLabels(dm, idm);CHKERRQ(ierr); 1370b325530aSVaclav Hapla ierr = PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_interpolate_orient_interfaces", &flg, NULL);CHKERRQ(ierr); 1371b325530aSVaclav Hapla if (flg) {ierr = DMPlexOrientInterface(idm);CHKERRQ(ierr);} 137284699f70SSatish Balay } 137343eeeb2dSStefano Zampini { 137443eeeb2dSStefano Zampini PetscBool isper; 137543eeeb2dSStefano Zampini const PetscReal *maxCell, *L; 137643eeeb2dSStefano Zampini const DMBoundaryType *bd; 137743eeeb2dSStefano Zampini 137843eeeb2dSStefano Zampini ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr); 137943eeeb2dSStefano Zampini ierr = DMSetPeriodicity(idm,isper,maxCell,L,bd);CHKERRQ(ierr); 138043eeeb2dSStefano Zampini } 13819a5b13a2SMatthew G Knepley *dmInt = idm; 1382a72f3261SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 13839f074e33SMatthew G Knepley PetscFunctionReturn(0); 13849f074e33SMatthew G Knepley } 138507b0a1fcSMatthew G Knepley 138680330477SMatthew G. Knepley /*@ 138780330477SMatthew G. Knepley DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices 138880330477SMatthew G. Knepley 138980330477SMatthew G. Knepley Collective on DM 139080330477SMatthew G. Knepley 139180330477SMatthew G. Knepley Input Parameter: 139280330477SMatthew G. Knepley . dmA - The DMPlex object with initial coordinates 139380330477SMatthew G. Knepley 139480330477SMatthew G. Knepley Output Parameter: 139580330477SMatthew G. Knepley . dmB - The DMPlex object with copied coordinates 139680330477SMatthew G. Knepley 139780330477SMatthew G. Knepley Level: intermediate 139880330477SMatthew G. Knepley 139980330477SMatthew G. Knepley Note: This is typically used when adding pieces other than vertices to a mesh 140080330477SMatthew G. Knepley 140180330477SMatthew G. Knepley .keywords: mesh 140265f90189SJed Brown .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection() 140380330477SMatthew G. Knepley @*/ 140407b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB) 140507b0a1fcSMatthew G Knepley { 140607b0a1fcSMatthew G Knepley Vec coordinatesA, coordinatesB; 140743eeeb2dSStefano Zampini VecType vtype; 140807b0a1fcSMatthew G Knepley PetscSection coordSectionA, coordSectionB; 140907b0a1fcSMatthew G Knepley PetscScalar *coordsA, *coordsB; 14100bedd6beSMatthew G. Knepley PetscInt spaceDim, Nf, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d; 141143eeeb2dSStefano Zampini PetscInt cStartA, cEndA, cStartB, cEndB, cS, cE; 141243eeeb2dSStefano Zampini PetscBool lc = PETSC_FALSE; 141307b0a1fcSMatthew G Knepley PetscErrorCode ierr; 141407b0a1fcSMatthew G Knepley 141507b0a1fcSMatthew G Knepley PetscFunctionBegin; 141643eeeb2dSStefano Zampini PetscValidHeaderSpecific(dmA, DM_CLASSID, 1); 141743eeeb2dSStefano Zampini PetscValidHeaderSpecific(dmB, DM_CLASSID, 2); 141876b791aaSMatthew G. Knepley if (dmA == dmB) PetscFunctionReturn(0); 141907b0a1fcSMatthew G Knepley ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr); 142007b0a1fcSMatthew G Knepley ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr); 142107b0a1fcSMatthew 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); 142243eeeb2dSStefano Zampini ierr = DMPlexGetHeightStratum(dmA, 0, &cStartA, &cEndA);CHKERRQ(ierr); 142343eeeb2dSStefano Zampini ierr = DMPlexGetHeightStratum(dmB, 0, &cStartB, &cEndB);CHKERRQ(ierr); 142469d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr); 142569d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr); 1426972bc18aSToby Isaac if (coordSectionA == coordSectionB) PetscFunctionReturn(0); 14270bedd6beSMatthew G. Knepley ierr = PetscSectionGetNumFields(coordSectionA, &Nf);CHKERRQ(ierr); 14280bedd6beSMatthew G. Knepley if (!Nf) PetscFunctionReturn(0); 14290bedd6beSMatthew G. Knepley if (Nf > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate fields must be 1, not %D", Nf); 1430df26b574SMatthew G. Knepley if (!coordSectionB) { 1431df26b574SMatthew G. Knepley PetscInt dim; 1432df26b574SMatthew G. Knepley 1433df26b574SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr); 1434df26b574SMatthew G. Knepley ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr); 1435df26b574SMatthew G. Knepley ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr); 1436df26b574SMatthew G. Knepley ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr); 1437df26b574SMatthew G. Knepley } 143807b0a1fcSMatthew G Knepley ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr); 143907b0a1fcSMatthew G Knepley ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr); 144007b0a1fcSMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr); 144143eeeb2dSStefano Zampini ierr = PetscSectionGetChart(coordSectionA, &cS, &cE);CHKERRQ(ierr); 144243eeeb2dSStefano Zampini if (cStartA <= cS && cS < cEndA) { /* localized coordinates */ 1443367003a6SStefano 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); 144443eeeb2dSStefano Zampini cS = cS - cStartA + cStartB; 144543eeeb2dSStefano Zampini cE = vEndB; 144643eeeb2dSStefano Zampini lc = PETSC_TRUE; 144743eeeb2dSStefano Zampini } else { 144843eeeb2dSStefano Zampini cS = vStartB; 144943eeeb2dSStefano Zampini cE = vEndB; 145043eeeb2dSStefano Zampini } 145143eeeb2dSStefano Zampini ierr = PetscSectionSetChart(coordSectionB, cS, cE);CHKERRQ(ierr); 145207b0a1fcSMatthew G Knepley for (v = vStartB; v < vEndB; ++v) { 145307b0a1fcSMatthew G Knepley ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr); 145407b0a1fcSMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr); 145507b0a1fcSMatthew G Knepley } 145643eeeb2dSStefano Zampini if (lc) { /* localized coordinates */ 145743eeeb2dSStefano Zampini PetscInt c; 145843eeeb2dSStefano Zampini 145943eeeb2dSStefano Zampini for (c = cS-cStartB; c < cEndB-cStartB; c++) { 146043eeeb2dSStefano Zampini PetscInt dof; 146143eeeb2dSStefano Zampini 146243eeeb2dSStefano Zampini ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr); 146343eeeb2dSStefano Zampini ierr = PetscSectionSetDof(coordSectionB, c + cStartB, dof);CHKERRQ(ierr); 146443eeeb2dSStefano Zampini ierr = PetscSectionSetFieldDof(coordSectionB, c + cStartB, 0, dof);CHKERRQ(ierr); 146543eeeb2dSStefano Zampini } 146643eeeb2dSStefano Zampini } 146707b0a1fcSMatthew G Knepley ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr); 146807b0a1fcSMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr); 146907b0a1fcSMatthew G Knepley ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr); 14708b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinatesB);CHKERRQ(ierr); 147107b0a1fcSMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr); 147207b0a1fcSMatthew G Knepley ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr); 147343eeeb2dSStefano Zampini ierr = VecGetBlockSize(coordinatesA, &d);CHKERRQ(ierr); 147443eeeb2dSStefano Zampini ierr = VecSetBlockSize(coordinatesB, d);CHKERRQ(ierr); 147543eeeb2dSStefano Zampini ierr = VecGetType(coordinatesA, &vtype);CHKERRQ(ierr); 147643eeeb2dSStefano Zampini ierr = VecSetType(coordinatesB, vtype);CHKERRQ(ierr); 147707b0a1fcSMatthew G Knepley ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr); 147807b0a1fcSMatthew G Knepley ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr); 147907b0a1fcSMatthew G Knepley for (v = 0; v < vEndB-vStartB; ++v) { 148043eeeb2dSStefano Zampini PetscInt offA, offB; 148143eeeb2dSStefano Zampini 148243eeeb2dSStefano Zampini ierr = PetscSectionGetOffset(coordSectionA, v + vStartA, &offA);CHKERRQ(ierr); 148343eeeb2dSStefano Zampini ierr = PetscSectionGetOffset(coordSectionB, v + vStartB, &offB);CHKERRQ(ierr); 148407b0a1fcSMatthew G Knepley for (d = 0; d < spaceDim; ++d) { 148543eeeb2dSStefano Zampini coordsB[offB+d] = coordsA[offA+d]; 148643eeeb2dSStefano Zampini } 148743eeeb2dSStefano Zampini } 148843eeeb2dSStefano Zampini if (lc) { /* localized coordinates */ 148943eeeb2dSStefano Zampini PetscInt c; 149043eeeb2dSStefano Zampini 149143eeeb2dSStefano Zampini for (c = cS-cStartB; c < cEndB-cStartB; c++) { 149243eeeb2dSStefano Zampini PetscInt dof, offA, offB; 149343eeeb2dSStefano Zampini 149443eeeb2dSStefano Zampini ierr = PetscSectionGetOffset(coordSectionA, c + cStartA, &offA);CHKERRQ(ierr); 149543eeeb2dSStefano Zampini ierr = PetscSectionGetOffset(coordSectionB, c + cStartB, &offB);CHKERRQ(ierr); 149643eeeb2dSStefano Zampini ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr); 149743eeeb2dSStefano Zampini ierr = PetscMemcpy(coordsB + offB,coordsA + offA,dof*sizeof(*coordsB));CHKERRQ(ierr); 149807b0a1fcSMatthew G Knepley } 149907b0a1fcSMatthew G Knepley } 150007b0a1fcSMatthew G Knepley ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr); 150107b0a1fcSMatthew G Knepley ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr); 150207b0a1fcSMatthew G Knepley ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr); 150307b0a1fcSMatthew G Knepley ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr); 150407b0a1fcSMatthew G Knepley PetscFunctionReturn(0); 150507b0a1fcSMatthew G Knepley } 15065c386225SMatthew G. Knepley 15074e3744c5SMatthew G. Knepley /*@ 15084e3744c5SMatthew G. Knepley DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh 15094e3744c5SMatthew G. Knepley 15104e3744c5SMatthew G. Knepley Collective on DM 15114e3744c5SMatthew G. Knepley 15124e3744c5SMatthew G. Knepley Input Parameter: 15134e3744c5SMatthew G. Knepley . dm - The complete DMPlex object 15144e3744c5SMatthew G. Knepley 15154e3744c5SMatthew G. Knepley Output Parameter: 15164e3744c5SMatthew G. Knepley . dmUnint - The DMPlex object with only cells and vertices 15174e3744c5SMatthew G. Knepley 15184e3744c5SMatthew G. Knepley Level: intermediate 15194e3744c5SMatthew G. Knepley 152095452b02SPatrick Sanan Notes: 152195452b02SPatrick Sanan It does not copy over the coordinates. 152243eeeb2dSStefano Zampini 15234e3744c5SMatthew G. Knepley .keywords: mesh 152443eeeb2dSStefano Zampini .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates() 15254e3744c5SMatthew G. Knepley @*/ 15264e3744c5SMatthew G. Knepley PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint) 15274e3744c5SMatthew G. Knepley { 15284e3744c5SMatthew G. Knepley DM udm; 1529c9f63434SStefano Zampini PetscInt dim, vStart, vEnd, cStart, cEnd, cMax, c, maxConeSize = 0, *cone; 15304e3744c5SMatthew G. Knepley PetscErrorCode ierr; 15314e3744c5SMatthew G. Knepley 15324e3744c5SMatthew G. Knepley PetscFunctionBegin; 153343eeeb2dSStefano Zampini PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 153443eeeb2dSStefano Zampini PetscValidPointer(dmUnint, 2); 1535c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 15364e3744c5SMatthew G. Knepley if (dim <= 1) { 15374e3744c5SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 1538595d4abbSMatthew G. Knepley *dmUnint = dm; 1539595d4abbSMatthew G. Knepley PetscFunctionReturn(0); 15404e3744c5SMatthew G. Knepley } 1541595d4abbSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1542595d4abbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1543c9f63434SStefano Zampini ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 15444e3744c5SMatthew G. Knepley ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr); 15454e3744c5SMatthew G. Knepley ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr); 1546c73cfb54SMatthew G. Knepley ierr = DMSetDimension(udm, dim);CHKERRQ(ierr); 15474e3744c5SMatthew G. Knepley ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr); 15484e3744c5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1549595d4abbSMatthew G. Knepley PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 15504e3744c5SMatthew G. Knepley 15514e3744c5SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 15524e3744c5SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 15534e3744c5SMatthew G. Knepley const PetscInt p = closure[cl]; 15544e3744c5SMatthew G. Knepley 15554e3744c5SMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) ++coneSize; 15564e3744c5SMatthew G. Knepley } 15574e3744c5SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 15584e3744c5SMatthew G. Knepley ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr); 1559595d4abbSMatthew G. Knepley maxConeSize = PetscMax(maxConeSize, coneSize); 15604e3744c5SMatthew G. Knepley } 15614e3744c5SMatthew G. Knepley ierr = DMSetUp(udm);CHKERRQ(ierr); 1562785e854fSJed Brown ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr); 15634e3744c5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1564595d4abbSMatthew G. Knepley PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 15654e3744c5SMatthew G. Knepley 15664e3744c5SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 15674e3744c5SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 15684e3744c5SMatthew G. Knepley const PetscInt p = closure[cl]; 15694e3744c5SMatthew G. Knepley 15704e3744c5SMatthew G. Knepley if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p; 15714e3744c5SMatthew G. Knepley } 15724e3744c5SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 15734e3744c5SMatthew G. Knepley ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr); 15744e3744c5SMatthew G. Knepley } 15754e3744c5SMatthew G. Knepley ierr = PetscFree(cone);CHKERRQ(ierr); 1576c9f63434SStefano Zampini ierr = DMPlexSetHybridBounds(udm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 15774e3744c5SMatthew G. Knepley ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr); 15784e3744c5SMatthew G. Knepley ierr = DMPlexStratify(udm);CHKERRQ(ierr); 15795c7de58aSMatthew G. Knepley /* Reduce SF */ 15805c7de58aSMatthew G. Knepley { 15815c7de58aSMatthew G. Knepley PetscSF sfPoint, sfPointUn; 15825c7de58aSMatthew G. Knepley const PetscSFNode *remotePoints; 15835c7de58aSMatthew G. Knepley const PetscInt *localPoints; 15845c7de58aSMatthew G. Knepley PetscSFNode *remotePointsUn; 15855c7de58aSMatthew G. Knepley PetscInt *localPointsUn; 15865c7de58aSMatthew G. Knepley PetscInt vEnd, numRoots, numLeaves, l; 15875c7de58aSMatthew G. Knepley PetscInt numLeavesUn = 0, n = 0; 15885c7de58aSMatthew G. Knepley PetscErrorCode ierr; 15895c7de58aSMatthew G. Knepley 15905c7de58aSMatthew G. Knepley /* Get original SF information */ 15915c7de58aSMatthew G. Knepley ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 15925c7de58aSMatthew G. Knepley ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr); 15935c7de58aSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr); 15945c7de58aSMatthew G. Knepley ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 15955c7de58aSMatthew G. Knepley /* Allocate space for cells and vertices */ 15965c7de58aSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++; 15975c7de58aSMatthew G. Knepley /* Fill in leaves */ 15985c7de58aSMatthew G. Knepley if (vEnd >= 0) { 15995c7de58aSMatthew G. Knepley ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr); 16005c7de58aSMatthew G. Knepley ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr); 16015c7de58aSMatthew G. Knepley for (l = 0; l < numLeaves; l++) { 16025c7de58aSMatthew G. Knepley if (localPoints[l] < vEnd) { 16035c7de58aSMatthew G. Knepley localPointsUn[n] = localPoints[l]; 16045c7de58aSMatthew G. Knepley remotePointsUn[n].rank = remotePoints[l].rank; 16055c7de58aSMatthew G. Knepley remotePointsUn[n].index = remotePoints[l].index; 16065c7de58aSMatthew G. Knepley ++n; 16075c7de58aSMatthew G. Knepley } 16085c7de58aSMatthew G. Knepley } 16095c7de58aSMatthew G. Knepley if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn); 16105c7de58aSMatthew G. Knepley ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr); 16115c7de58aSMatthew G. Knepley } 16125c7de58aSMatthew G. Knepley } 161343eeeb2dSStefano Zampini { 161443eeeb2dSStefano Zampini PetscBool isper; 161543eeeb2dSStefano Zampini const PetscReal *maxCell, *L; 161643eeeb2dSStefano Zampini const DMBoundaryType *bd; 161743eeeb2dSStefano Zampini 161843eeeb2dSStefano Zampini ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr); 161943eeeb2dSStefano Zampini ierr = DMSetPeriodicity(udm,isper,maxCell,L,bd);CHKERRQ(ierr); 162043eeeb2dSStefano Zampini } 162143eeeb2dSStefano Zampini 16224e3744c5SMatthew G. Knepley *dmUnint = udm; 16234e3744c5SMatthew G. Knepley PetscFunctionReturn(0); 16244e3744c5SMatthew G. Knepley } 1625