1552f7358SJed Brown #define PETSCDM_DLL 234541f0dSBarry Smith #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3552f7358SJed Brown #include <petscdmda.h> 40c312b8eSJed Brown #include <petscsf.h> 5552f7358SJed Brown 6552f7358SJed Brown #undef __FUNCT__ 7552f7358SJed Brown #define __FUNCT__ "DMSetFromOptions_Plex" 8552f7358SJed Brown PetscErrorCode DMSetFromOptions_Plex(DM dm) 9552f7358SJed Brown { 10552f7358SJed Brown DM_Plex *mesh = (DM_Plex*) dm->data; 11552f7358SJed Brown PetscErrorCode ierr; 12552f7358SJed Brown 13552f7358SJed Brown PetscFunctionBegin; 14552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 15552f7358SJed Brown ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr); 16552f7358SJed Brown /* Handle DMPlex refinement */ 17552f7358SJed Brown /* Handle associated vectors */ 18552f7358SJed Brown /* Handle viewing */ 190298fd71SBarry Smith ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr); 200298fd71SBarry Smith ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr); 21552f7358SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 22552f7358SJed Brown PetscFunctionReturn(0); 23552f7358SJed Brown } 24552f7358SJed Brown 25552f7358SJed Brown #undef __FUNCT__ 26552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareBoundary" 27552f7358SJed Brown /* 28552f7358SJed Brown Simple square boundary: 29552f7358SJed Brown 30552f7358SJed Brown 18--5-17--4--16 31552f7358SJed Brown | | | 32552f7358SJed Brown 6 10 3 33552f7358SJed Brown | | | 34552f7358SJed Brown 19-11-20--9--15 35552f7358SJed Brown | | | 36552f7358SJed Brown 7 8 2 37552f7358SJed Brown | | | 38552f7358SJed Brown 12--0-13--1--14 39552f7358SJed Brown */ 40552f7358SJed Brown PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 41552f7358SJed Brown { 42552f7358SJed Brown PetscInt numVertices = (edges[0]+1)*(edges[1]+1); 43552f7358SJed Brown PetscInt numEdges = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1]; 44552f7358SJed Brown PetscInt markerTop = 1; 45552f7358SJed Brown PetscInt markerBottom = 1; 46552f7358SJed Brown PetscInt markerRight = 1; 47552f7358SJed Brown PetscInt markerLeft = 1; 48552f7358SJed Brown PetscBool markerSeparate = PETSC_FALSE; 49552f7358SJed Brown Vec coordinates; 50552f7358SJed Brown PetscSection coordSection; 51552f7358SJed Brown PetscScalar *coords; 52552f7358SJed Brown PetscInt coordSize; 53552f7358SJed Brown PetscMPIInt rank; 54552f7358SJed Brown PetscInt v, vx, vy; 55552f7358SJed Brown PetscErrorCode ierr; 56552f7358SJed Brown 57552f7358SJed Brown PetscFunctionBegin; 580298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 59552f7358SJed Brown if (markerSeparate) { 60552f7358SJed Brown markerTop = 1; 61552f7358SJed Brown markerBottom = 0; 62552f7358SJed Brown markerRight = 0; 63552f7358SJed Brown markerLeft = 0; 64552f7358SJed Brown } 6582f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 66552f7358SJed Brown if (!rank) { 67552f7358SJed Brown PetscInt e, ex, ey; 68552f7358SJed Brown 69552f7358SJed Brown ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr); 70552f7358SJed Brown for (e = 0; e < numEdges; ++e) { 71552f7358SJed Brown ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 72552f7358SJed Brown } 73552f7358SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 74552f7358SJed Brown for (vx = 0; vx <= edges[0]; vx++) { 75552f7358SJed Brown for (ey = 0; ey < edges[1]; ey++) { 76552f7358SJed Brown PetscInt edge = vx*edges[1] + ey + edges[0]*(edges[1]+1); 77552f7358SJed Brown PetscInt vertex = ey*(edges[0]+1) + vx + numEdges; 78da80777bSKarl Rupp PetscInt cone[2]; 79552f7358SJed Brown 80da80777bSKarl Rupp cone[0] = vertex; cone[1] = vertex+edges[0]+1; 81552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 82552f7358SJed Brown if (vx == edges[0]) { 83552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 84552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 85552f7358SJed Brown if (ey == edges[1]-1) { 86552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 87552f7358SJed Brown } 88552f7358SJed Brown } else if (vx == 0) { 89552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 90552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 91552f7358SJed Brown if (ey == edges[1]-1) { 92552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 93552f7358SJed Brown } 94552f7358SJed Brown } 95552f7358SJed Brown } 96552f7358SJed Brown } 97552f7358SJed Brown for (vy = 0; vy <= edges[1]; vy++) { 98552f7358SJed Brown for (ex = 0; ex < edges[0]; ex++) { 99552f7358SJed Brown PetscInt edge = vy*edges[0] + ex; 100552f7358SJed Brown PetscInt vertex = vy*(edges[0]+1) + ex + numEdges; 101da80777bSKarl Rupp PetscInt cone[2]; 102552f7358SJed Brown 103da80777bSKarl Rupp cone[0] = vertex; cone[1] = vertex+1; 104552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 105552f7358SJed Brown if (vy == edges[1]) { 106552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 107552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 108552f7358SJed Brown if (ex == edges[0]-1) { 109552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 110552f7358SJed Brown } 111552f7358SJed Brown } else if (vy == 0) { 112552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 113552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 114552f7358SJed Brown if (ex == edges[0]-1) { 115552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 116552f7358SJed Brown } 117552f7358SJed Brown } 118552f7358SJed Brown } 119552f7358SJed Brown } 120552f7358SJed Brown } 121552f7358SJed Brown ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 122552f7358SJed Brown ierr = DMPlexStratify(dm);CHKERRQ(ierr); 123552f7358SJed Brown /* Build coordinates */ 124552f7358SJed Brown ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 125552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr); 126552f7358SJed Brown for (v = numEdges; v < numEdges+numVertices; ++v) { 127552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 128552f7358SJed Brown } 129552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 130552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 13182f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 132552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 133552f7358SJed Brown ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 134552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 135552f7358SJed Brown for (vy = 0; vy <= edges[1]; ++vy) { 136552f7358SJed Brown for (vx = 0; vx <= edges[0]; ++vx) { 137552f7358SJed Brown coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx; 138552f7358SJed Brown coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy; 139552f7358SJed Brown } 140552f7358SJed Brown } 141552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 142552f7358SJed Brown ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 143552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 144552f7358SJed Brown PetscFunctionReturn(0); 145552f7358SJed Brown } 146552f7358SJed Brown 147552f7358SJed Brown #undef __FUNCT__ 148552f7358SJed Brown #define __FUNCT__ "DMPlexCreateCubeBoundary" 149552f7358SJed Brown /* 150552f7358SJed Brown Simple cubic boundary: 151552f7358SJed Brown 152552f7358SJed Brown 2-------3 153552f7358SJed Brown /| /| 154552f7358SJed Brown 6-------7 | 155552f7358SJed Brown | | | | 156552f7358SJed Brown | 0-----|-1 157552f7358SJed Brown |/ |/ 158552f7358SJed Brown 4-------5 159552f7358SJed Brown */ 160552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[]) 161552f7358SJed Brown { 162552f7358SJed Brown PetscInt numVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1); 163552f7358SJed Brown PetscInt numFaces = 6; 164552f7358SJed Brown Vec coordinates; 165552f7358SJed Brown PetscSection coordSection; 166552f7358SJed Brown PetscScalar *coords; 167552f7358SJed Brown PetscInt coordSize; 168552f7358SJed Brown PetscMPIInt rank; 169552f7358SJed Brown PetscInt v, vx, vy, vz; 170552f7358SJed Brown PetscErrorCode ierr; 171552f7358SJed Brown 172552f7358SJed Brown PetscFunctionBegin; 17382f516ccSBarry Smith if ((faces[0] < 1) || (faces[1] < 1) || (faces[2] < 1)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Must have at least 1 face per side"); 17482f516ccSBarry Smith if ((faces[0] > 1) || (faces[1] > 1) || (faces[2] > 1)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Currently can't handle more than 1 face per side"); 175552f7358SJed Brown ierr = PetscMalloc(numVertices*2 * sizeof(PetscReal), &coords);CHKERRQ(ierr); 17682f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 177552f7358SJed Brown if (!rank) { 178552f7358SJed Brown PetscInt f; 179552f7358SJed Brown 180552f7358SJed Brown ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr); 181552f7358SJed Brown for (f = 0; f < numFaces; ++f) { 182552f7358SJed Brown ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 183552f7358SJed Brown } 184552f7358SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 185552f7358SJed Brown for (v = 0; v < numFaces+numVertices; ++v) { 186552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr); 187552f7358SJed Brown } 188552f7358SJed Brown { /* Side 0 (Front) */ 189da80777bSKarl Rupp PetscInt cone[4]; 190da80777bSKarl Rupp cone[0] = numFaces+4; cone[1] = numFaces+5; cone[2] = numFaces+7; cone[3] = numFaces+6; 191552f7358SJed Brown ierr = DMPlexSetCone(dm, 0, cone);CHKERRQ(ierr); 192552f7358SJed Brown } 193552f7358SJed Brown { /* Side 1 (Back) */ 194da80777bSKarl Rupp PetscInt cone[4]; 195da80777bSKarl Rupp cone[0] = numFaces+1; cone[1] = numFaces+0; cone[2] = numFaces+2; cone[3] = numFaces+3; 196552f7358SJed Brown ierr = DMPlexSetCone(dm, 1, cone);CHKERRQ(ierr); 197552f7358SJed Brown } 198552f7358SJed Brown { /* Side 2 (Bottom) */ 199da80777bSKarl Rupp PetscInt cone[4]; 200da80777bSKarl Rupp cone[0] = numFaces+0; cone[1] = numFaces+1; cone[2] = numFaces+5; cone[3] = numFaces+4; 201552f7358SJed Brown ierr = DMPlexSetCone(dm, 2, cone);CHKERRQ(ierr); 202552f7358SJed Brown } 203552f7358SJed Brown { /* Side 3 (Top) */ 204da80777bSKarl Rupp PetscInt cone[4]; 205da80777bSKarl Rupp cone[0] = numFaces+6; cone[1] = numFaces+7; cone[2] = numFaces+3; cone[3] = numFaces+2; 206552f7358SJed Brown ierr = DMPlexSetCone(dm, 3, cone);CHKERRQ(ierr); 207552f7358SJed Brown } 208552f7358SJed Brown { /* Side 4 (Left) */ 209da80777bSKarl Rupp PetscInt cone[4]; 210da80777bSKarl Rupp cone[0] = numFaces+0; cone[1] = numFaces+4; cone[2] = numFaces+6; cone[3] = numFaces+2; 211552f7358SJed Brown ierr = DMPlexSetCone(dm, 4, cone);CHKERRQ(ierr); 212552f7358SJed Brown } 213552f7358SJed Brown { /* Side 5 (Right) */ 214da80777bSKarl Rupp PetscInt cone[4]; 215da80777bSKarl Rupp cone[0] = numFaces+5; cone[1] = numFaces+1; cone[2] = numFaces+3; cone[3] = numFaces+7; 216552f7358SJed Brown ierr = DMPlexSetCone(dm, 5, cone);CHKERRQ(ierr); 217552f7358SJed Brown } 218552f7358SJed Brown } 219552f7358SJed Brown ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 220552f7358SJed Brown ierr = DMPlexStratify(dm);CHKERRQ(ierr); 221552f7358SJed Brown /* Build coordinates */ 222552f7358SJed Brown ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 223552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr); 224552f7358SJed Brown for (v = numFaces; v < numFaces+numVertices; ++v) { 225552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr); 226552f7358SJed Brown } 227552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 228552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 22982f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 230552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 231552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 232552f7358SJed Brown ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 233552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 234552f7358SJed Brown for (vz = 0; vz <= faces[2]; ++vz) { 235552f7358SJed Brown for (vy = 0; vy <= faces[1]; ++vy) { 236552f7358SJed Brown for (vx = 0; vx <= faces[0]; ++vx) { 237552f7358SJed Brown coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx; 238552f7358SJed Brown coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy; 239552f7358SJed Brown coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz; 240552f7358SJed Brown } 241552f7358SJed Brown } 242552f7358SJed Brown } 243552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 244552f7358SJed Brown ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 245552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 246552f7358SJed Brown PetscFunctionReturn(0); 247552f7358SJed Brown } 248552f7358SJed Brown 249552f7358SJed Brown #undef __FUNCT__ 250552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareMesh" 251552f7358SJed Brown /* 252552f7358SJed Brown Simple square mesh: 253552f7358SJed Brown 254552f7358SJed Brown 22--8-23--9--24 255552f7358SJed Brown | | | 256552f7358SJed Brown 13 2 14 3 15 257552f7358SJed Brown | | | 258552f7358SJed Brown 19--6-20--7--21 259552f7358SJed Brown | | | 260552f7358SJed Brown 10 0 11 1 12 261552f7358SJed Brown | | | 262552f7358SJed Brown 16--4-17--5--18 263552f7358SJed Brown */ 264552f7358SJed Brown PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 265552f7358SJed Brown { 266552f7358SJed Brown PetscInt markerTop = 1; 267552f7358SJed Brown PetscInt markerBottom = 1; 268552f7358SJed Brown PetscInt markerRight = 1; 269552f7358SJed Brown PetscInt markerLeft = 1; 270552f7358SJed Brown PetscBool markerSeparate = PETSC_FALSE; 271552f7358SJed Brown PetscMPIInt rank; 272552f7358SJed Brown PetscErrorCode ierr; 273552f7358SJed Brown 274552f7358SJed Brown PetscFunctionBegin; 27582f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 2760298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 277552f7358SJed Brown if (markerSeparate) { 278552f7358SJed Brown markerTop = 3; 279552f7358SJed Brown markerBottom = 1; 280552f7358SJed Brown markerRight = 2; 281552f7358SJed Brown markerLeft = 4; 282552f7358SJed Brown } 283552f7358SJed Brown { 284552f7358SJed Brown const PetscInt numXEdges = !rank ? edges[0] : 0; 285552f7358SJed Brown const PetscInt numYEdges = !rank ? edges[1] : 0; 286552f7358SJed Brown const PetscInt numXVertices = !rank ? edges[0]+1 : 0; 287552f7358SJed Brown const PetscInt numYVertices = !rank ? edges[1]+1 : 0; 288552f7358SJed Brown const PetscInt numTotXEdges = numXEdges*numYVertices; 289552f7358SJed Brown const PetscInt numTotYEdges = numYEdges*numXVertices; 290552f7358SJed Brown const PetscInt numVertices = numXVertices*numYVertices; 291552f7358SJed Brown const PetscInt numEdges = numTotXEdges + numTotYEdges; 292552f7358SJed Brown const PetscInt numFaces = numXEdges*numYEdges; 293552f7358SJed Brown const PetscInt firstVertex = numFaces; 294552f7358SJed Brown const PetscInt firstXEdge = numFaces + numVertices; 295552f7358SJed Brown const PetscInt firstYEdge = numFaces + numVertices + numTotXEdges; 296552f7358SJed Brown Vec coordinates; 297552f7358SJed Brown PetscSection coordSection; 298552f7358SJed Brown PetscScalar *coords; 299552f7358SJed Brown PetscInt coordSize; 300552f7358SJed Brown PetscInt v, vx, vy; 301552f7358SJed Brown PetscInt f, fx, fy, e, ex, ey; 302552f7358SJed Brown 303552f7358SJed Brown ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr); 304552f7358SJed Brown for (f = 0; f < numFaces; ++f) { 305552f7358SJed Brown ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 306552f7358SJed Brown } 307552f7358SJed Brown for (e = firstXEdge; e < firstXEdge+numEdges; ++e) { 308552f7358SJed Brown ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 309552f7358SJed Brown } 310552f7358SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 311552f7358SJed Brown /* Build faces */ 312552f7358SJed Brown for (fy = 0; fy < numYEdges; fy++) { 313552f7358SJed Brown for (fx = 0; fx < numXEdges; fx++) { 314552f7358SJed Brown const PetscInt face = fy*numXEdges + fx; 315552f7358SJed Brown const PetscInt edgeL = firstYEdge + fx*numYEdges + fy; 316552f7358SJed Brown const PetscInt edgeB = firstXEdge + fy*numXEdges + fx; 317552f7358SJed Brown const PetscInt ornt[4] = {0, 0, -2, -2}; 318da80777bSKarl Rupp PetscInt cone[4]; 319552f7358SJed Brown 320da80777bSKarl Rupp cone[0] = edgeB; cone[1] = edgeL+numYEdges; cone[2] = edgeB+numXEdges; cone[3] = edgeL; 321552f7358SJed Brown ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 322552f7358SJed Brown ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 323552f7358SJed Brown } 324552f7358SJed Brown } 325552f7358SJed Brown /* Build Y edges*/ 326552f7358SJed Brown for (vx = 0; vx < numXVertices; vx++) { 327552f7358SJed Brown for (ey = 0; ey < numYEdges; ey++) { 328552f7358SJed Brown const PetscInt edge = firstYEdge + vx*numYEdges + ey; 329552f7358SJed Brown const PetscInt vertex = firstVertex + ey*numXVertices + vx; 330da80777bSKarl Rupp PetscInt cone[2]; 331552f7358SJed Brown 332da80777bSKarl Rupp cone[0] = vertex; cone[1] = vertex+numXVertices; 333552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 334552f7358SJed Brown if (vx == numXVertices-1) { 335552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 336552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 337552f7358SJed Brown if (ey == numYEdges-1) { 338552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 339552f7358SJed Brown } 340552f7358SJed Brown } else if (vx == 0) { 341552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 342552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 343552f7358SJed Brown if (ey == numYEdges-1) { 344552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 345552f7358SJed Brown } 346552f7358SJed Brown } 347552f7358SJed Brown } 348552f7358SJed Brown } 349552f7358SJed Brown /* Build X edges*/ 350552f7358SJed Brown for (vy = 0; vy < numYVertices; vy++) { 351552f7358SJed Brown for (ex = 0; ex < numXEdges; ex++) { 352552f7358SJed Brown const PetscInt edge = firstXEdge + vy*numXEdges + ex; 353552f7358SJed Brown const PetscInt vertex = firstVertex + vy*numXVertices + ex; 354da80777bSKarl Rupp PetscInt cone[2]; 355552f7358SJed Brown 356da80777bSKarl Rupp cone[0] = vertex; cone[1] = vertex+1; 357552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 358552f7358SJed Brown if (vy == numYVertices-1) { 359552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 360552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 361552f7358SJed Brown if (ex == numXEdges-1) { 362552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 363552f7358SJed Brown } 364552f7358SJed Brown } else if (vy == 0) { 365552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 366552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 367552f7358SJed Brown if (ex == numXEdges-1) { 368552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 369552f7358SJed Brown } 370552f7358SJed Brown } 371552f7358SJed Brown } 372552f7358SJed Brown } 373552f7358SJed Brown ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 374552f7358SJed Brown ierr = DMPlexStratify(dm);CHKERRQ(ierr); 375552f7358SJed Brown /* Build coordinates */ 376552f7358SJed Brown ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 377552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr); 378552f7358SJed Brown for (v = firstVertex; v < firstVertex+numVertices; ++v) { 379552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 380552f7358SJed Brown } 381552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 382552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 38382f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 384552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 385552f7358SJed Brown ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 386552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 387552f7358SJed Brown for (vy = 0; vy < numYVertices; ++vy) { 388552f7358SJed Brown for (vx = 0; vx < numXVertices; ++vx) { 389552f7358SJed Brown coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx; 390552f7358SJed Brown coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy; 391552f7358SJed Brown } 392552f7358SJed Brown } 393552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 394552f7358SJed Brown ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 395552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 396552f7358SJed Brown } 397552f7358SJed Brown PetscFunctionReturn(0); 398552f7358SJed Brown } 399552f7358SJed Brown 400552f7358SJed Brown #undef __FUNCT__ 401552f7358SJed Brown #define __FUNCT__ "DMPlexCreateBoxMesh" 4026b75479cSMatthew G Knepley /*@ 4036b75479cSMatthew G Knepley DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box). 4046b75479cSMatthew G Knepley 4056b75479cSMatthew G Knepley Collective on MPI_Comm 4066b75479cSMatthew G Knepley 4076b75479cSMatthew G Knepley Input Parameters: 4086b75479cSMatthew G Knepley + comm - The communicator for the DM object 4096b75479cSMatthew G Knepley . dim - The spatial dimension 4106b75479cSMatthew G Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces) 4116b75479cSMatthew G Knepley 4126b75479cSMatthew G Knepley Output Parameter: 4136b75479cSMatthew G Knepley . dm - The DM object 4146b75479cSMatthew G Knepley 4156b75479cSMatthew G Knepley Level: beginner 4166b75479cSMatthew G Knepley 4176b75479cSMatthew G Knepley .keywords: DM, create 4186b75479cSMatthew G Knepley .seealso: DMSetType(), DMCreate() 4196b75479cSMatthew G Knepley @*/ 4200adebc6cSBarry Smith PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm) 4210adebc6cSBarry Smith { 422552f7358SJed Brown DM boundary; 423552f7358SJed Brown PetscErrorCode ierr; 424552f7358SJed Brown 425552f7358SJed Brown PetscFunctionBegin; 426552f7358SJed Brown PetscValidPointer(dm, 4); 427552f7358SJed Brown ierr = DMCreate(comm, &boundary);CHKERRQ(ierr); 428552f7358SJed Brown PetscValidLogicalCollectiveInt(boundary,dim,2); 429552f7358SJed Brown ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr); 430552f7358SJed Brown ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr); 431552f7358SJed Brown switch (dim) { 432552f7358SJed Brown case 2: 433552f7358SJed Brown { 434552f7358SJed Brown PetscReal lower[2] = {0.0, 0.0}; 435552f7358SJed Brown PetscReal upper[2] = {1.0, 1.0}; 436552f7358SJed Brown PetscInt edges[2] = {2, 2}; 437552f7358SJed Brown 438552f7358SJed Brown ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr); 439552f7358SJed Brown break; 440552f7358SJed Brown } 441552f7358SJed Brown case 3: 442552f7358SJed Brown { 443552f7358SJed Brown PetscReal lower[3] = {0.0, 0.0, 0.0}; 444552f7358SJed Brown PetscReal upper[3] = {1.0, 1.0, 1.0}; 445552f7358SJed Brown PetscInt faces[3] = {1, 1, 1}; 446552f7358SJed Brown 447552f7358SJed Brown ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr); 448552f7358SJed Brown break; 449552f7358SJed Brown } 450552f7358SJed Brown default: 451552f7358SJed Brown SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 452552f7358SJed Brown } 4530298fd71SBarry Smith ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr); 454552f7358SJed Brown ierr = DMDestroy(&boundary);CHKERRQ(ierr); 455552f7358SJed Brown PetscFunctionReturn(0); 456552f7358SJed Brown } 457552f7358SJed Brown 458552f7358SJed Brown #undef __FUNCT__ 459552f7358SJed Brown #define __FUNCT__ "DMPlexCreateHexBoxMesh" 460a6dfd86eSKarl Rupp PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DM *dm) 461a6dfd86eSKarl Rupp { 462552f7358SJed Brown PetscErrorCode ierr; 463552f7358SJed Brown 464552f7358SJed Brown PetscFunctionBegin; 465552f7358SJed Brown PetscValidPointer(dm, 4); 466552f7358SJed Brown ierr = DMCreate(comm, dm);CHKERRQ(ierr); 467552f7358SJed Brown PetscValidLogicalCollectiveInt(*dm,dim,2); 468552f7358SJed Brown ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 469552f7358SJed Brown ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 470552f7358SJed Brown switch (dim) { 471552f7358SJed Brown case 2: 472552f7358SJed Brown { 473552f7358SJed Brown PetscReal lower[2] = {0.0, 0.0}; 474552f7358SJed Brown PetscReal upper[2] = {1.0, 1.0}; 475552f7358SJed Brown 476552f7358SJed Brown ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells);CHKERRQ(ierr); 477552f7358SJed Brown break; 478552f7358SJed Brown } 479552f7358SJed Brown #if 0 480552f7358SJed Brown case 3: 481552f7358SJed Brown { 482552f7358SJed Brown PetscReal lower[3] = {0.0, 0.0, 0.0}; 483552f7358SJed Brown PetscReal upper[3] = {1.0, 1.0, 1.0}; 484552f7358SJed Brown 485552f7358SJed Brown ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells);CHKERRQ(ierr); 486552f7358SJed Brown break; 487552f7358SJed Brown } 488552f7358SJed Brown #endif 489552f7358SJed Brown default: 490552f7358SJed Brown SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 491552f7358SJed Brown } 492552f7358SJed Brown PetscFunctionReturn(0); 493552f7358SJed Brown } 494552f7358SJed Brown 495552f7358SJed Brown /* External function declarations here */ 496552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 497552f7358SJed Brown extern PetscErrorCode DMCreateMatrix_Plex(DM dm, MatType mtype, Mat *J); 498552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 499552f7358SJed Brown extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined); 500*9a922efaSMatthew G. Knepley extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 501552f7358SJed Brown extern PetscErrorCode DMSetUp_Plex(DM dm); 502552f7358SJed Brown extern PetscErrorCode DMDestroy_Plex(DM dm); 503552f7358SJed Brown extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 504552f7358SJed Brown extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 505552f7358SJed Brown extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS); 506552f7358SJed Brown 507552f7358SJed Brown #undef __FUNCT__ 508552f7358SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Plex" 509552f7358SJed Brown static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 510552f7358SJed Brown { 511552f7358SJed Brown PetscErrorCode ierr; 512552f7358SJed Brown 513552f7358SJed Brown PetscFunctionBegin; 514552f7358SJed Brown ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 515552f7358SJed Brown /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 516552f7358SJed Brown ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex);CHKERRQ(ierr); 517552f7358SJed Brown PetscFunctionReturn(0); 518552f7358SJed Brown } 519552f7358SJed Brown 520552f7358SJed Brown #undef __FUNCT__ 521552f7358SJed Brown #define __FUNCT__ "DMCreateLocalVector_Plex" 522552f7358SJed Brown static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 523552f7358SJed Brown { 524552f7358SJed Brown PetscErrorCode ierr; 525552f7358SJed Brown 526552f7358SJed Brown PetscFunctionBegin; 527552f7358SJed Brown ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 528552f7358SJed Brown ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 529552f7358SJed Brown PetscFunctionReturn(0); 530552f7358SJed Brown } 531552f7358SJed Brown 53238221697SMatthew G. Knepley #undef __FUNCT__ 533552f7358SJed Brown #define __FUNCT__ "DMInitialize_Plex" 534552f7358SJed Brown PetscErrorCode DMInitialize_Plex(DM dm) 535552f7358SJed Brown { 536552f7358SJed Brown PetscErrorCode ierr; 537552f7358SJed Brown 538552f7358SJed Brown PetscFunctionBegin; 539552f7358SJed Brown ierr = PetscStrallocpy(VECSTANDARD, (char**)&dm->vectype);CHKERRQ(ierr); 5408865f1eaSKarl Rupp 541552f7358SJed Brown dm->ops->view = DMView_Plex; 542552f7358SJed Brown dm->ops->setfromoptions = DMSetFromOptions_Plex; 54338221697SMatthew G. Knepley dm->ops->clone = DMClone_Plex; 544552f7358SJed Brown dm->ops->setup = DMSetUp_Plex; 545552f7358SJed Brown dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 546552f7358SJed Brown dm->ops->createlocalvector = DMCreateLocalVector_Plex; 547184d77edSJed Brown dm->ops->getlocaltoglobalmapping = NULL; 548184d77edSJed Brown dm->ops->getlocaltoglobalmappingblock = NULL; 5490298fd71SBarry Smith dm->ops->createfieldis = NULL; 550552f7358SJed Brown dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 551552f7358SJed Brown dm->ops->getcoloring = 0; 552552f7358SJed Brown dm->ops->creatematrix = DMCreateMatrix_Plex; 553552f7358SJed Brown dm->ops->createinterpolation = 0; 554552f7358SJed Brown dm->ops->getaggregates = 0; 555552f7358SJed Brown dm->ops->getinjection = 0; 556552f7358SJed Brown dm->ops->refine = DMRefine_Plex; 557552f7358SJed Brown dm->ops->coarsen = 0; 558552f7358SJed Brown dm->ops->refinehierarchy = 0; 559552f7358SJed Brown dm->ops->coarsenhierarchy = 0; 5600298fd71SBarry Smith dm->ops->globaltolocalbegin = NULL; 5610298fd71SBarry Smith dm->ops->globaltolocalend = NULL; 5620298fd71SBarry Smith dm->ops->localtoglobalbegin = NULL; 5630298fd71SBarry Smith dm->ops->localtoglobalend = NULL; 564552f7358SJed Brown dm->ops->destroy = DMDestroy_Plex; 565552f7358SJed Brown dm->ops->createsubdm = DMCreateSubDM_Plex; 566552f7358SJed Brown dm->ops->locatepoints = DMLocatePoints_Plex; 567552f7358SJed Brown PetscFunctionReturn(0); 568552f7358SJed Brown } 569552f7358SJed Brown 57063a16f15SMatthew G. Knepley #undef __FUNCT__ 57163a16f15SMatthew G. Knepley #define __FUNCT__ "DMClone_Plex" 57263a16f15SMatthew G. Knepley PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 57363a16f15SMatthew G. Knepley { 57463a16f15SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 57563a16f15SMatthew G. Knepley PetscErrorCode ierr; 57663a16f15SMatthew G. Knepley 57763a16f15SMatthew G. Knepley PetscFunctionBegin; 57863a16f15SMatthew G. Knepley mesh->refct++; 57963a16f15SMatthew G. Knepley (*newdm)->data = mesh; 58063a16f15SMatthew G. Knepley ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 58163a16f15SMatthew G. Knepley ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 58263a16f15SMatthew G. Knepley PetscFunctionReturn(0); 58363a16f15SMatthew G. Knepley } 58463a16f15SMatthew G. Knepley 5858818961aSMatthew G Knepley /*MC 5868818961aSMatthew G Knepley DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 5878818961aSMatthew G Knepley In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 5888818961aSMatthew G Knepley specified by a PetscSection object. Ownership in the global representation is determined by 5898818961aSMatthew G Knepley ownership of the underlying DMPlex points. This is specified by another PetscSection object. 5908818961aSMatthew G Knepley 5918818961aSMatthew G Knepley Level: intermediate 5928818961aSMatthew G Knepley 5938818961aSMatthew G Knepley .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 5948818961aSMatthew G Knepley M*/ 5958818961aSMatthew G Knepley 596552f7358SJed Brown #undef __FUNCT__ 597552f7358SJed Brown #define __FUNCT__ "DMCreate_Plex" 5988cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 599552f7358SJed Brown { 600552f7358SJed Brown DM_Plex *mesh; 601770b213bSMatthew G Knepley PetscInt unit, d; 602552f7358SJed Brown PetscErrorCode ierr; 603552f7358SJed Brown 604552f7358SJed Brown PetscFunctionBegin; 605552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 606552f7358SJed Brown ierr = PetscNewLog(dm, DM_Plex, &mesh);CHKERRQ(ierr); 607552f7358SJed Brown dm->data = mesh; 608552f7358SJed Brown 609552f7358SJed Brown mesh->refct = 1; 610552f7358SJed Brown mesh->dim = 0; 61182f516ccSBarry Smith ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 612552f7358SJed Brown mesh->maxConeSize = 0; 6130298fd71SBarry Smith mesh->cones = NULL; 6140298fd71SBarry Smith mesh->coneOrientations = NULL; 61582f516ccSBarry Smith ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 616552f7358SJed Brown mesh->maxSupportSize = 0; 6170298fd71SBarry Smith mesh->supports = NULL; 618552f7358SJed Brown mesh->refinementUniform = PETSC_TRUE; 619552f7358SJed Brown mesh->refinementLimit = -1.0; 620552f7358SJed Brown 6210298fd71SBarry Smith mesh->facesTmp = NULL; 622552f7358SJed Brown 6230298fd71SBarry Smith mesh->subpointMap = NULL; 624552f7358SJed Brown 6258865f1eaSKarl Rupp for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 626552f7358SJed Brown 6270298fd71SBarry Smith mesh->labels = NULL; 6280298fd71SBarry Smith mesh->globalVertexNumbers = NULL; 6290298fd71SBarry Smith mesh->globalCellNumbers = NULL; 6308865f1eaSKarl Rupp for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 631552f7358SJed Brown mesh->vtkCellHeight = 0; 632552f7358SJed Brown mesh->preallocCenterDim = -1; 633552f7358SJed Brown 6340298fd71SBarry Smith mesh->integrateResidualFEM = NULL; 6350298fd71SBarry Smith mesh->integrateJacobianActionFEM = NULL; 6360298fd71SBarry Smith mesh->integrateJacobianFEM = NULL; 637552f7358SJed Brown 638552f7358SJed Brown mesh->printSetValues = PETSC_FALSE; 639552f7358SJed Brown mesh->printFEM = 0; 640552f7358SJed Brown 641552f7358SJed Brown ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 642552f7358SJed Brown PetscFunctionReturn(0); 643552f7358SJed Brown } 644552f7358SJed Brown 645552f7358SJed Brown #undef __FUNCT__ 646552f7358SJed Brown #define __FUNCT__ "DMPlexCreate" 647552f7358SJed Brown /*@ 648552f7358SJed Brown DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 649552f7358SJed Brown 650552f7358SJed Brown Collective on MPI_Comm 651552f7358SJed Brown 652552f7358SJed Brown Input Parameter: 653552f7358SJed Brown . comm - The communicator for the DMPlex object 654552f7358SJed Brown 655552f7358SJed Brown Output Parameter: 656552f7358SJed Brown . mesh - The DMPlex object 657552f7358SJed Brown 658552f7358SJed Brown Level: beginner 659552f7358SJed Brown 660552f7358SJed Brown .keywords: DMPlex, create 661552f7358SJed Brown @*/ 662552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 663552f7358SJed Brown { 664552f7358SJed Brown PetscErrorCode ierr; 665552f7358SJed Brown 666552f7358SJed Brown PetscFunctionBegin; 667552f7358SJed Brown PetscValidPointer(mesh,2); 668552f7358SJed Brown ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 669552f7358SJed Brown ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 670552f7358SJed Brown PetscFunctionReturn(0); 671552f7358SJed Brown } 672552f7358SJed Brown 673552f7358SJed Brown #undef __FUNCT__ 6749298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildFromCellList_Private" 6759298eaa6SMatthew G Knepley /* 6769298eaa6SMatthew G Knepley This takes as input the common mesh generator output, a list of the vertices for each cell 6779298eaa6SMatthew G Knepley */ 6789298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 6799298eaa6SMatthew G Knepley { 6809298eaa6SMatthew G Knepley PetscInt *cone, c, p; 6819298eaa6SMatthew G Knepley PetscErrorCode ierr; 6829298eaa6SMatthew G Knepley 6839298eaa6SMatthew G Knepley PetscFunctionBegin; 6849298eaa6SMatthew G Knepley ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 6859298eaa6SMatthew G Knepley for (c = 0; c < numCells; ++c) { 6869298eaa6SMatthew G Knepley ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 6879298eaa6SMatthew G Knepley } 6889298eaa6SMatthew G Knepley ierr = DMSetUp(dm);CHKERRQ(ierr); 6899298eaa6SMatthew G Knepley ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 6909298eaa6SMatthew G Knepley for (c = 0; c < numCells; ++c) { 6919298eaa6SMatthew G Knepley for (p = 0; p < numCorners; ++p) { 6929298eaa6SMatthew G Knepley cone[p] = cells[c*numCorners+p]+numCells; 6939298eaa6SMatthew G Knepley } 6949298eaa6SMatthew G Knepley ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 6959298eaa6SMatthew G Knepley } 6969298eaa6SMatthew G Knepley ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 6979298eaa6SMatthew G Knepley ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 6989298eaa6SMatthew G Knepley ierr = DMPlexStratify(dm);CHKERRQ(ierr); 6999298eaa6SMatthew G Knepley PetscFunctionReturn(0); 7009298eaa6SMatthew G Knepley } 7019298eaa6SMatthew G Knepley 7029298eaa6SMatthew G Knepley #undef __FUNCT__ 7039298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildCoordinates_Private" 7049298eaa6SMatthew G Knepley /* 7059298eaa6SMatthew G Knepley This takes as input the coordinates for each vertex 7069298eaa6SMatthew G Knepley */ 7079298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 7089298eaa6SMatthew G Knepley { 7099298eaa6SMatthew G Knepley PetscSection coordSection; 7109298eaa6SMatthew G Knepley Vec coordinates; 7119298eaa6SMatthew G Knepley PetscScalar *coords; 7129298eaa6SMatthew G Knepley PetscInt coordSize, v, d; 7139298eaa6SMatthew G Knepley PetscErrorCode ierr; 7149298eaa6SMatthew G Knepley 7159298eaa6SMatthew G Knepley PetscFunctionBegin; 7169298eaa6SMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 7179298eaa6SMatthew G Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 7189298eaa6SMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 7199298eaa6SMatthew G Knepley ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 7209298eaa6SMatthew G Knepley for (v = numCells; v < numCells+numVertices; ++v) { 7219298eaa6SMatthew G Knepley ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 7229298eaa6SMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 7239298eaa6SMatthew G Knepley } 7249298eaa6SMatthew G Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 7259298eaa6SMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 7269298eaa6SMatthew G Knepley ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 7279298eaa6SMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 7289298eaa6SMatthew G Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 7299298eaa6SMatthew G Knepley ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 7309298eaa6SMatthew G Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 7319298eaa6SMatthew G Knepley for (v = 0; v < numVertices; ++v) { 7329298eaa6SMatthew G Knepley for (d = 0; d < spaceDim; ++d) { 7339298eaa6SMatthew G Knepley coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 7349298eaa6SMatthew G Knepley } 7359298eaa6SMatthew G Knepley } 7369298eaa6SMatthew G Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 7379298eaa6SMatthew G Knepley ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 7389298eaa6SMatthew G Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 7399298eaa6SMatthew G Knepley PetscFunctionReturn(0); 7409298eaa6SMatthew G Knepley } 7419298eaa6SMatthew G Knepley 7429298eaa6SMatthew G Knepley #undef __FUNCT__ 7439298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromCellList" 7449298eaa6SMatthew G Knepley /*@C 7459298eaa6SMatthew G Knepley DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 7469298eaa6SMatthew G Knepley 7479298eaa6SMatthew G Knepley Input Parameters: 7489298eaa6SMatthew G Knepley + comm - The communicator 7499298eaa6SMatthew G Knepley . dim - The topological dimension of the mesh 7509298eaa6SMatthew G Knepley . numCells - The number of cells 7519298eaa6SMatthew G Knepley . numVertices - The number of vertices 7529298eaa6SMatthew G Knepley . numCorners - The number of vertices for each cell 7539298eaa6SMatthew G Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 7549298eaa6SMatthew G Knepley . cells - An array of numCells*numCorners numbers, the vertices for each cell 7559298eaa6SMatthew G Knepley . spaceDim - The spatial dimension used for coordinates 7569298eaa6SMatthew G Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 7579298eaa6SMatthew G Knepley 7589298eaa6SMatthew G Knepley Output Parameter: 7599298eaa6SMatthew G Knepley . dm - The DM 7609298eaa6SMatthew G Knepley 7619298eaa6SMatthew G Knepley Note: Two triangles sharing a face 7629298eaa6SMatthew G Knepley $ 7639298eaa6SMatthew G Knepley $ 2 7649298eaa6SMatthew G Knepley $ / | \ 7659298eaa6SMatthew G Knepley $ / | \ 7669298eaa6SMatthew G Knepley $ / | \ 7679298eaa6SMatthew G Knepley $ 0 0 | 1 3 7689298eaa6SMatthew G Knepley $ \ | / 7699298eaa6SMatthew G Knepley $ \ | / 7709298eaa6SMatthew G Knepley $ \ | / 7719298eaa6SMatthew G Knepley $ 1 7729298eaa6SMatthew G Knepley would have input 7739298eaa6SMatthew G Knepley $ numCells = 2, numVertices = 4 7749298eaa6SMatthew G Knepley $ cells = [0 1 2 1 3 2] 7759298eaa6SMatthew G Knepley $ 7769298eaa6SMatthew G Knepley which would result in the DMPlex 7779298eaa6SMatthew G Knepley $ 7789298eaa6SMatthew G Knepley $ 4 7799298eaa6SMatthew G Knepley $ / | \ 7809298eaa6SMatthew G Knepley $ / | \ 7819298eaa6SMatthew G Knepley $ / | \ 7829298eaa6SMatthew G Knepley $ 2 0 | 1 5 7839298eaa6SMatthew G Knepley $ \ | / 7849298eaa6SMatthew G Knepley $ \ | / 7859298eaa6SMatthew G Knepley $ \ | / 7869298eaa6SMatthew G Knepley $ 3 7879298eaa6SMatthew G Knepley 7889298eaa6SMatthew G Knepley Level: beginner 7899298eaa6SMatthew G Knepley 7909298eaa6SMatthew G Knepley .seealso: DMPlexCreate() 7919298eaa6SMatthew G Knepley @*/ 7929298eaa6SMatthew G Knepley PetscErrorCode DMPlexCreateFromCellList(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const double vertexCoords[], DM *dm) 7939298eaa6SMatthew G Knepley { 7949298eaa6SMatthew G Knepley PetscErrorCode ierr; 7959298eaa6SMatthew G Knepley 7969298eaa6SMatthew G Knepley PetscFunctionBegin; 7979298eaa6SMatthew G Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 7989298eaa6SMatthew G Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 7999298eaa6SMatthew G Knepley ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 8009298eaa6SMatthew G Knepley ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 8019298eaa6SMatthew G Knepley if (interpolate) { 8029298eaa6SMatthew G Knepley DM idm; 8039298eaa6SMatthew G Knepley 8049298eaa6SMatthew G Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 8059298eaa6SMatthew G Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 8069298eaa6SMatthew G Knepley *dm = idm; 8079298eaa6SMatthew G Knepley } 8089298eaa6SMatthew G Knepley ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 8099298eaa6SMatthew G Knepley PetscFunctionReturn(0); 8109298eaa6SMatthew G Knepley } 8119298eaa6SMatthew G Knepley 8129298eaa6SMatthew G Knepley #undef __FUNCT__ 8139298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromDAG" 8149298eaa6SMatthew G Knepley /* 8159298eaa6SMatthew G Knepley This takes as input the raw Hasse Diagram data 8169298eaa6SMatthew G Knepley */ 8179298eaa6SMatthew G Knepley PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 8189298eaa6SMatthew G Knepley { 8199298eaa6SMatthew G Knepley Vec coordinates; 8209298eaa6SMatthew G Knepley PetscSection coordSection; 8219298eaa6SMatthew G Knepley PetscScalar *coords; 8229298eaa6SMatthew G Knepley PetscInt coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off; 8239298eaa6SMatthew G Knepley PetscErrorCode ierr; 8249298eaa6SMatthew G Knepley 8259298eaa6SMatthew G Knepley PetscFunctionBegin; 8269298eaa6SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 8279298eaa6SMatthew G Knepley for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 8289298eaa6SMatthew G Knepley ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 8299298eaa6SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 8309298eaa6SMatthew G Knepley ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 8319298eaa6SMatthew G Knepley } 8329298eaa6SMatthew G Knepley ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 8339298eaa6SMatthew G Knepley for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 8349298eaa6SMatthew G Knepley ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 8359298eaa6SMatthew G Knepley ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 8369298eaa6SMatthew G Knepley } 8379298eaa6SMatthew G Knepley ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 8389298eaa6SMatthew G Knepley ierr = DMPlexStratify(dm);CHKERRQ(ierr); 8399298eaa6SMatthew G Knepley /* Build coordinates */ 8409298eaa6SMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 8419298eaa6SMatthew G Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 8429298eaa6SMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 8439298eaa6SMatthew G Knepley ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 8449298eaa6SMatthew G Knepley for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 8459298eaa6SMatthew G Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 8466f8cbbeeSMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 8479298eaa6SMatthew G Knepley } 8489298eaa6SMatthew G Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 8499298eaa6SMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 8509298eaa6SMatthew G Knepley ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 8516f8cbbeeSMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 8529298eaa6SMatthew G Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 8539298eaa6SMatthew G Knepley ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 8549298eaa6SMatthew G Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 8559298eaa6SMatthew G Knepley for (v = 0; v < numPoints[0]; ++v) { 8569298eaa6SMatthew G Knepley PetscInt off; 8579298eaa6SMatthew G Knepley 8589298eaa6SMatthew G Knepley ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 8599298eaa6SMatthew G Knepley for (d = 0; d < dim; ++d) { 8609298eaa6SMatthew G Knepley coords[off+d] = vertexCoords[v*dim+d]; 8619298eaa6SMatthew G Knepley } 8629298eaa6SMatthew G Knepley } 8639298eaa6SMatthew G Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 8649298eaa6SMatthew G Knepley ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 8659298eaa6SMatthew G Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 8669298eaa6SMatthew G Knepley PetscFunctionReturn(0); 8679298eaa6SMatthew G Knepley } 868