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); 21*6113b454SMatthew G. Knepley ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", PETSC_FALSE, &mesh->printTol, NULL);CHKERRQ(ierr); 22552f7358SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 23552f7358SJed Brown PetscFunctionReturn(0); 24552f7358SJed Brown } 25552f7358SJed Brown 26552f7358SJed Brown #undef __FUNCT__ 27552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareBoundary" 28552f7358SJed Brown /* 29552f7358SJed Brown Simple square boundary: 30552f7358SJed Brown 31552f7358SJed Brown 18--5-17--4--16 32552f7358SJed Brown | | | 33552f7358SJed Brown 6 10 3 34552f7358SJed Brown | | | 35552f7358SJed Brown 19-11-20--9--15 36552f7358SJed Brown | | | 37552f7358SJed Brown 7 8 2 38552f7358SJed Brown | | | 39552f7358SJed Brown 12--0-13--1--14 40552f7358SJed Brown */ 41552f7358SJed Brown PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 42552f7358SJed Brown { 43552f7358SJed Brown PetscInt numVertices = (edges[0]+1)*(edges[1]+1); 44552f7358SJed Brown PetscInt numEdges = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1]; 45552f7358SJed Brown PetscInt markerTop = 1; 46552f7358SJed Brown PetscInt markerBottom = 1; 47552f7358SJed Brown PetscInt markerRight = 1; 48552f7358SJed Brown PetscInt markerLeft = 1; 49552f7358SJed Brown PetscBool markerSeparate = PETSC_FALSE; 50552f7358SJed Brown Vec coordinates; 51552f7358SJed Brown PetscSection coordSection; 52552f7358SJed Brown PetscScalar *coords; 53552f7358SJed Brown PetscInt coordSize; 54552f7358SJed Brown PetscMPIInt rank; 55552f7358SJed Brown PetscInt v, vx, vy; 56552f7358SJed Brown PetscErrorCode ierr; 57552f7358SJed Brown 58552f7358SJed Brown PetscFunctionBegin; 590298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 60552f7358SJed Brown if (markerSeparate) { 61552f7358SJed Brown markerTop = 1; 62552f7358SJed Brown markerBottom = 0; 63552f7358SJed Brown markerRight = 0; 64552f7358SJed Brown markerLeft = 0; 65552f7358SJed Brown } 6682f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 67552f7358SJed Brown if (!rank) { 68552f7358SJed Brown PetscInt e, ex, ey; 69552f7358SJed Brown 70552f7358SJed Brown ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr); 71552f7358SJed Brown for (e = 0; e < numEdges; ++e) { 72552f7358SJed Brown ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 73552f7358SJed Brown } 74552f7358SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 75552f7358SJed Brown for (vx = 0; vx <= edges[0]; vx++) { 76552f7358SJed Brown for (ey = 0; ey < edges[1]; ey++) { 77552f7358SJed Brown PetscInt edge = vx*edges[1] + ey + edges[0]*(edges[1]+1); 78552f7358SJed Brown PetscInt vertex = ey*(edges[0]+1) + vx + numEdges; 79da80777bSKarl Rupp PetscInt cone[2]; 80552f7358SJed Brown 81da80777bSKarl Rupp cone[0] = vertex; cone[1] = vertex+edges[0]+1; 82552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 83552f7358SJed Brown if (vx == edges[0]) { 84552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 85552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 86552f7358SJed Brown if (ey == edges[1]-1) { 87552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 88552f7358SJed Brown } 89552f7358SJed Brown } else if (vx == 0) { 90552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 91552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 92552f7358SJed Brown if (ey == edges[1]-1) { 93552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 94552f7358SJed Brown } 95552f7358SJed Brown } 96552f7358SJed Brown } 97552f7358SJed Brown } 98552f7358SJed Brown for (vy = 0; vy <= edges[1]; vy++) { 99552f7358SJed Brown for (ex = 0; ex < edges[0]; ex++) { 100552f7358SJed Brown PetscInt edge = vy*edges[0] + ex; 101552f7358SJed Brown PetscInt vertex = vy*(edges[0]+1) + ex + numEdges; 102da80777bSKarl Rupp PetscInt cone[2]; 103552f7358SJed Brown 104da80777bSKarl Rupp cone[0] = vertex; cone[1] = vertex+1; 105552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 106552f7358SJed Brown if (vy == edges[1]) { 107552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 108552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 109552f7358SJed Brown if (ex == edges[0]-1) { 110552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 111552f7358SJed Brown } 112552f7358SJed Brown } else if (vy == 0) { 113552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 114552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 115552f7358SJed Brown if (ex == edges[0]-1) { 116552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 117552f7358SJed Brown } 118552f7358SJed Brown } 119552f7358SJed Brown } 120552f7358SJed Brown } 121552f7358SJed Brown } 122552f7358SJed Brown ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 123552f7358SJed Brown ierr = DMPlexStratify(dm);CHKERRQ(ierr); 124552f7358SJed Brown /* Build coordinates */ 125552f7358SJed Brown ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 126552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr); 127552f7358SJed Brown for (v = numEdges; v < numEdges+numVertices; ++v) { 128552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 129552f7358SJed Brown } 130552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 131552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 13282f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 133552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 134552f7358SJed Brown ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 135552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 136552f7358SJed Brown for (vy = 0; vy <= edges[1]; ++vy) { 137552f7358SJed Brown for (vx = 0; vx <= edges[0]; ++vx) { 138552f7358SJed Brown coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx; 139552f7358SJed Brown coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy; 140552f7358SJed Brown } 141552f7358SJed Brown } 142552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 143552f7358SJed Brown ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 144552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 145552f7358SJed Brown PetscFunctionReturn(0); 146552f7358SJed Brown } 147552f7358SJed Brown 148552f7358SJed Brown #undef __FUNCT__ 149552f7358SJed Brown #define __FUNCT__ "DMPlexCreateCubeBoundary" 150552f7358SJed Brown /* 151552f7358SJed Brown Simple cubic boundary: 152552f7358SJed Brown 153552f7358SJed Brown 2-------3 154552f7358SJed Brown /| /| 155552f7358SJed Brown 6-------7 | 156552f7358SJed Brown | | | | 157552f7358SJed Brown | 0-----|-1 158552f7358SJed Brown |/ |/ 159552f7358SJed Brown 4-------5 160552f7358SJed Brown */ 161552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[]) 162552f7358SJed Brown { 163552f7358SJed Brown PetscInt numVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1); 164552f7358SJed Brown PetscInt numFaces = 6; 165552f7358SJed Brown Vec coordinates; 166552f7358SJed Brown PetscSection coordSection; 167552f7358SJed Brown PetscScalar *coords; 168552f7358SJed Brown PetscInt coordSize; 169552f7358SJed Brown PetscMPIInt rank; 170552f7358SJed Brown PetscInt v, vx, vy, vz; 171552f7358SJed Brown PetscErrorCode ierr; 172552f7358SJed Brown 173552f7358SJed Brown PetscFunctionBegin; 17482f516ccSBarry 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"); 17582f516ccSBarry 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"); 176552f7358SJed Brown ierr = PetscMalloc(numVertices*2 * sizeof(PetscReal), &coords);CHKERRQ(ierr); 17782f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 178552f7358SJed Brown if (!rank) { 179552f7358SJed Brown PetscInt f; 180552f7358SJed Brown 181552f7358SJed Brown ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr); 182552f7358SJed Brown for (f = 0; f < numFaces; ++f) { 183552f7358SJed Brown ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 184552f7358SJed Brown } 185552f7358SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 186552f7358SJed Brown for (v = 0; v < numFaces+numVertices; ++v) { 187552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr); 188552f7358SJed Brown } 189552f7358SJed Brown { /* Side 0 (Front) */ 190da80777bSKarl Rupp PetscInt cone[4]; 191da80777bSKarl Rupp cone[0] = numFaces+4; cone[1] = numFaces+5; cone[2] = numFaces+7; cone[3] = numFaces+6; 192552f7358SJed Brown ierr = DMPlexSetCone(dm, 0, cone);CHKERRQ(ierr); 193552f7358SJed Brown } 194552f7358SJed Brown { /* Side 1 (Back) */ 195da80777bSKarl Rupp PetscInt cone[4]; 196da80777bSKarl Rupp cone[0] = numFaces+1; cone[1] = numFaces+0; cone[2] = numFaces+2; cone[3] = numFaces+3; 197552f7358SJed Brown ierr = DMPlexSetCone(dm, 1, cone);CHKERRQ(ierr); 198552f7358SJed Brown } 199552f7358SJed Brown { /* Side 2 (Bottom) */ 200da80777bSKarl Rupp PetscInt cone[4]; 201da80777bSKarl Rupp cone[0] = numFaces+0; cone[1] = numFaces+1; cone[2] = numFaces+5; cone[3] = numFaces+4; 202552f7358SJed Brown ierr = DMPlexSetCone(dm, 2, cone);CHKERRQ(ierr); 203552f7358SJed Brown } 204552f7358SJed Brown { /* Side 3 (Top) */ 205da80777bSKarl Rupp PetscInt cone[4]; 206da80777bSKarl Rupp cone[0] = numFaces+6; cone[1] = numFaces+7; cone[2] = numFaces+3; cone[3] = numFaces+2; 207552f7358SJed Brown ierr = DMPlexSetCone(dm, 3, cone);CHKERRQ(ierr); 208552f7358SJed Brown } 209552f7358SJed Brown { /* Side 4 (Left) */ 210da80777bSKarl Rupp PetscInt cone[4]; 211da80777bSKarl Rupp cone[0] = numFaces+0; cone[1] = numFaces+4; cone[2] = numFaces+6; cone[3] = numFaces+2; 212552f7358SJed Brown ierr = DMPlexSetCone(dm, 4, cone);CHKERRQ(ierr); 213552f7358SJed Brown } 214552f7358SJed Brown { /* Side 5 (Right) */ 215da80777bSKarl Rupp PetscInt cone[4]; 216da80777bSKarl Rupp cone[0] = numFaces+5; cone[1] = numFaces+1; cone[2] = numFaces+3; cone[3] = numFaces+7; 217552f7358SJed Brown ierr = DMPlexSetCone(dm, 5, cone);CHKERRQ(ierr); 218552f7358SJed Brown } 219552f7358SJed Brown } 220552f7358SJed Brown ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 221552f7358SJed Brown ierr = DMPlexStratify(dm);CHKERRQ(ierr); 222552f7358SJed Brown /* Build coordinates */ 223552f7358SJed Brown ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 224552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr); 225552f7358SJed Brown for (v = numFaces; v < numFaces+numVertices; ++v) { 226552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr); 227552f7358SJed Brown } 228552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 229552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 23082f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 231552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 232552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 233552f7358SJed Brown ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 234552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 235552f7358SJed Brown for (vz = 0; vz <= faces[2]; ++vz) { 236552f7358SJed Brown for (vy = 0; vy <= faces[1]; ++vy) { 237552f7358SJed Brown for (vx = 0; vx <= faces[0]; ++vx) { 238552f7358SJed Brown coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx; 239552f7358SJed Brown coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy; 240552f7358SJed Brown coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz; 241552f7358SJed Brown } 242552f7358SJed Brown } 243552f7358SJed Brown } 244552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 245552f7358SJed Brown ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 246552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 247552f7358SJed Brown PetscFunctionReturn(0); 248552f7358SJed Brown } 249552f7358SJed Brown 250552f7358SJed Brown #undef __FUNCT__ 251552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareMesh" 252552f7358SJed Brown /* 253552f7358SJed Brown Simple square mesh: 254552f7358SJed Brown 255552f7358SJed Brown 22--8-23--9--24 256552f7358SJed Brown | | | 257552f7358SJed Brown 13 2 14 3 15 258552f7358SJed Brown | | | 259552f7358SJed Brown 19--6-20--7--21 260552f7358SJed Brown | | | 261552f7358SJed Brown 10 0 11 1 12 262552f7358SJed Brown | | | 263552f7358SJed Brown 16--4-17--5--18 264552f7358SJed Brown */ 265552f7358SJed Brown PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 266552f7358SJed Brown { 267552f7358SJed Brown PetscInt markerTop = 1; 268552f7358SJed Brown PetscInt markerBottom = 1; 269552f7358SJed Brown PetscInt markerRight = 1; 270552f7358SJed Brown PetscInt markerLeft = 1; 271552f7358SJed Brown PetscBool markerSeparate = PETSC_FALSE; 272552f7358SJed Brown PetscMPIInt rank; 273552f7358SJed Brown PetscErrorCode ierr; 274552f7358SJed Brown 275552f7358SJed Brown PetscFunctionBegin; 27682f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 2770298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 278552f7358SJed Brown if (markerSeparate) { 279552f7358SJed Brown markerTop = 3; 280552f7358SJed Brown markerBottom = 1; 281552f7358SJed Brown markerRight = 2; 282552f7358SJed Brown markerLeft = 4; 283552f7358SJed Brown } 284552f7358SJed Brown { 285552f7358SJed Brown const PetscInt numXEdges = !rank ? edges[0] : 0; 286552f7358SJed Brown const PetscInt numYEdges = !rank ? edges[1] : 0; 287552f7358SJed Brown const PetscInt numXVertices = !rank ? edges[0]+1 : 0; 288552f7358SJed Brown const PetscInt numYVertices = !rank ? edges[1]+1 : 0; 289552f7358SJed Brown const PetscInt numTotXEdges = numXEdges*numYVertices; 290552f7358SJed Brown const PetscInt numTotYEdges = numYEdges*numXVertices; 291552f7358SJed Brown const PetscInt numVertices = numXVertices*numYVertices; 292552f7358SJed Brown const PetscInt numEdges = numTotXEdges + numTotYEdges; 293552f7358SJed Brown const PetscInt numFaces = numXEdges*numYEdges; 294552f7358SJed Brown const PetscInt firstVertex = numFaces; 295552f7358SJed Brown const PetscInt firstXEdge = numFaces + numVertices; 296552f7358SJed Brown const PetscInt firstYEdge = numFaces + numVertices + numTotXEdges; 297552f7358SJed Brown Vec coordinates; 298552f7358SJed Brown PetscSection coordSection; 299552f7358SJed Brown PetscScalar *coords; 300552f7358SJed Brown PetscInt coordSize; 301552f7358SJed Brown PetscInt v, vx, vy; 302552f7358SJed Brown PetscInt f, fx, fy, e, ex, ey; 303552f7358SJed Brown 304552f7358SJed Brown ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr); 305552f7358SJed Brown for (f = 0; f < numFaces; ++f) { 306552f7358SJed Brown ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 307552f7358SJed Brown } 308552f7358SJed Brown for (e = firstXEdge; e < firstXEdge+numEdges; ++e) { 309552f7358SJed Brown ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 310552f7358SJed Brown } 311552f7358SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 312552f7358SJed Brown /* Build faces */ 313552f7358SJed Brown for (fy = 0; fy < numYEdges; fy++) { 314552f7358SJed Brown for (fx = 0; fx < numXEdges; fx++) { 315552f7358SJed Brown const PetscInt face = fy*numXEdges + fx; 316552f7358SJed Brown const PetscInt edgeL = firstYEdge + fx*numYEdges + fy; 317552f7358SJed Brown const PetscInt edgeB = firstXEdge + fy*numXEdges + fx; 318552f7358SJed Brown const PetscInt ornt[4] = {0, 0, -2, -2}; 319da80777bSKarl Rupp PetscInt cone[4]; 320552f7358SJed Brown 321da80777bSKarl Rupp cone[0] = edgeB; cone[1] = edgeL+numYEdges; cone[2] = edgeB+numXEdges; cone[3] = edgeL; 322552f7358SJed Brown ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 323552f7358SJed Brown ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 324552f7358SJed Brown } 325552f7358SJed Brown } 326552f7358SJed Brown /* Build Y edges*/ 327552f7358SJed Brown for (vx = 0; vx < numXVertices; vx++) { 328552f7358SJed Brown for (ey = 0; ey < numYEdges; ey++) { 329552f7358SJed Brown const PetscInt edge = firstYEdge + vx*numYEdges + ey; 330552f7358SJed Brown const PetscInt vertex = firstVertex + ey*numXVertices + vx; 331da80777bSKarl Rupp PetscInt cone[2]; 332552f7358SJed Brown 333da80777bSKarl Rupp cone[0] = vertex; cone[1] = vertex+numXVertices; 334552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 335552f7358SJed Brown if (vx == numXVertices-1) { 336552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 337552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 338552f7358SJed Brown if (ey == numYEdges-1) { 339552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 340552f7358SJed Brown } 341552f7358SJed Brown } else if (vx == 0) { 342552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 343552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 344552f7358SJed Brown if (ey == numYEdges-1) { 345552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 346552f7358SJed Brown } 347552f7358SJed Brown } 348552f7358SJed Brown } 349552f7358SJed Brown } 350552f7358SJed Brown /* Build X edges*/ 351552f7358SJed Brown for (vy = 0; vy < numYVertices; vy++) { 352552f7358SJed Brown for (ex = 0; ex < numXEdges; ex++) { 353552f7358SJed Brown const PetscInt edge = firstXEdge + vy*numXEdges + ex; 354552f7358SJed Brown const PetscInt vertex = firstVertex + vy*numXVertices + ex; 355da80777bSKarl Rupp PetscInt cone[2]; 356552f7358SJed Brown 357da80777bSKarl Rupp cone[0] = vertex; cone[1] = vertex+1; 358552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 359552f7358SJed Brown if (vy == numYVertices-1) { 360552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 361552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 362552f7358SJed Brown if (ex == numXEdges-1) { 363552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 364552f7358SJed Brown } 365552f7358SJed Brown } else if (vy == 0) { 366552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 367552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 368552f7358SJed Brown if (ex == numXEdges-1) { 369552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 370552f7358SJed Brown } 371552f7358SJed Brown } 372552f7358SJed Brown } 373552f7358SJed Brown } 374552f7358SJed Brown ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 375552f7358SJed Brown ierr = DMPlexStratify(dm);CHKERRQ(ierr); 376552f7358SJed Brown /* Build coordinates */ 377552f7358SJed Brown ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 378552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr); 379552f7358SJed Brown for (v = firstVertex; v < firstVertex+numVertices; ++v) { 380552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 381552f7358SJed Brown } 382552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 383552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 38482f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 385552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 386552f7358SJed Brown ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 387552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 388552f7358SJed Brown for (vy = 0; vy < numYVertices; ++vy) { 389552f7358SJed Brown for (vx = 0; vx < numXVertices; ++vx) { 390552f7358SJed Brown coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx; 391552f7358SJed Brown coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy; 392552f7358SJed Brown } 393552f7358SJed Brown } 394552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 395552f7358SJed Brown ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 396552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 397552f7358SJed Brown } 398552f7358SJed Brown PetscFunctionReturn(0); 399552f7358SJed Brown } 400552f7358SJed Brown 401552f7358SJed Brown #undef __FUNCT__ 402552f7358SJed Brown #define __FUNCT__ "DMPlexCreateBoxMesh" 4036b75479cSMatthew G Knepley /*@ 4046b75479cSMatthew G Knepley DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box). 4056b75479cSMatthew G Knepley 4066b75479cSMatthew G Knepley Collective on MPI_Comm 4076b75479cSMatthew G Knepley 4086b75479cSMatthew G Knepley Input Parameters: 4096b75479cSMatthew G Knepley + comm - The communicator for the DM object 4106b75479cSMatthew G Knepley . dim - The spatial dimension 4116b75479cSMatthew G Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces) 4126b75479cSMatthew G Knepley 4136b75479cSMatthew G Knepley Output Parameter: 4146b75479cSMatthew G Knepley . dm - The DM object 4156b75479cSMatthew G Knepley 4166b75479cSMatthew G Knepley Level: beginner 4176b75479cSMatthew G Knepley 4186b75479cSMatthew G Knepley .keywords: DM, create 4196b75479cSMatthew G Knepley .seealso: DMSetType(), DMCreate() 4206b75479cSMatthew G Knepley @*/ 4210adebc6cSBarry Smith PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm) 4220adebc6cSBarry Smith { 423552f7358SJed Brown DM boundary; 424552f7358SJed Brown PetscErrorCode ierr; 425552f7358SJed Brown 426552f7358SJed Brown PetscFunctionBegin; 427552f7358SJed Brown PetscValidPointer(dm, 4); 428552f7358SJed Brown ierr = DMCreate(comm, &boundary);CHKERRQ(ierr); 429552f7358SJed Brown PetscValidLogicalCollectiveInt(boundary,dim,2); 430552f7358SJed Brown ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr); 431552f7358SJed Brown ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr); 432552f7358SJed Brown switch (dim) { 433552f7358SJed Brown case 2: 434552f7358SJed Brown { 435552f7358SJed Brown PetscReal lower[2] = {0.0, 0.0}; 436552f7358SJed Brown PetscReal upper[2] = {1.0, 1.0}; 437552f7358SJed Brown PetscInt edges[2] = {2, 2}; 438552f7358SJed Brown 439552f7358SJed Brown ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr); 440552f7358SJed Brown break; 441552f7358SJed Brown } 442552f7358SJed Brown case 3: 443552f7358SJed Brown { 444552f7358SJed Brown PetscReal lower[3] = {0.0, 0.0, 0.0}; 445552f7358SJed Brown PetscReal upper[3] = {1.0, 1.0, 1.0}; 446552f7358SJed Brown PetscInt faces[3] = {1, 1, 1}; 447552f7358SJed Brown 448552f7358SJed Brown ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr); 449552f7358SJed Brown break; 450552f7358SJed Brown } 451552f7358SJed Brown default: 452552f7358SJed Brown SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 453552f7358SJed Brown } 4540298fd71SBarry Smith ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr); 455552f7358SJed Brown ierr = DMDestroy(&boundary);CHKERRQ(ierr); 456552f7358SJed Brown PetscFunctionReturn(0); 457552f7358SJed Brown } 458552f7358SJed Brown 459552f7358SJed Brown #undef __FUNCT__ 460552f7358SJed Brown #define __FUNCT__ "DMPlexCreateHexBoxMesh" 461a6dfd86eSKarl Rupp PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DM *dm) 462a6dfd86eSKarl Rupp { 463552f7358SJed Brown PetscErrorCode ierr; 464552f7358SJed Brown 465552f7358SJed Brown PetscFunctionBegin; 466552f7358SJed Brown PetscValidPointer(dm, 4); 467552f7358SJed Brown ierr = DMCreate(comm, dm);CHKERRQ(ierr); 468552f7358SJed Brown PetscValidLogicalCollectiveInt(*dm,dim,2); 469552f7358SJed Brown ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 470552f7358SJed Brown ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 471552f7358SJed Brown switch (dim) { 472552f7358SJed Brown case 2: 473552f7358SJed Brown { 474552f7358SJed Brown PetscReal lower[2] = {0.0, 0.0}; 475552f7358SJed Brown PetscReal upper[2] = {1.0, 1.0}; 476552f7358SJed Brown 477552f7358SJed Brown ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells);CHKERRQ(ierr); 478552f7358SJed Brown break; 479552f7358SJed Brown } 480552f7358SJed Brown #if 0 481552f7358SJed Brown case 3: 482552f7358SJed Brown { 483552f7358SJed Brown PetscReal lower[3] = {0.0, 0.0, 0.0}; 484552f7358SJed Brown PetscReal upper[3] = {1.0, 1.0, 1.0}; 485552f7358SJed Brown 486552f7358SJed Brown ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells);CHKERRQ(ierr); 487552f7358SJed Brown break; 488552f7358SJed Brown } 489552f7358SJed Brown #endif 490552f7358SJed Brown default: 491552f7358SJed Brown SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 492552f7358SJed Brown } 493552f7358SJed Brown PetscFunctionReturn(0); 494552f7358SJed Brown } 495552f7358SJed Brown 496552f7358SJed Brown /* External function declarations here */ 497552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 498552f7358SJed Brown extern PetscErrorCode DMCreateMatrix_Plex(DM dm, MatType mtype, Mat *J); 499552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 500552f7358SJed Brown extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined); 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 532552f7358SJed Brown 533552f7358SJed Brown #undef __FUNCT__ 534552f7358SJed Brown #define __FUNCT__ "DMInitialize_Plex" 535552f7358SJed Brown PetscErrorCode DMInitialize_Plex(DM dm) 536552f7358SJed Brown { 537552f7358SJed Brown PetscErrorCode ierr; 538552f7358SJed Brown 539552f7358SJed Brown PetscFunctionBegin; 540552f7358SJed Brown ierr = PetscStrallocpy(VECSTANDARD, (char**)&dm->vectype);CHKERRQ(ierr); 5418865f1eaSKarl Rupp 542552f7358SJed Brown dm->ops->view = DMView_Plex; 543552f7358SJed Brown dm->ops->setfromoptions = DMSetFromOptions_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 5708818961aSMatthew G Knepley /*MC 5718818961aSMatthew G Knepley DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 5728818961aSMatthew G Knepley In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 5738818961aSMatthew G Knepley specified by a PetscSection object. Ownership in the global representation is determined by 5748818961aSMatthew G Knepley ownership of the underlying DMPlex points. This is specified by another PetscSection object. 5758818961aSMatthew G Knepley 5768818961aSMatthew G Knepley Level: intermediate 5778818961aSMatthew G Knepley 5788818961aSMatthew G Knepley .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 5798818961aSMatthew G Knepley M*/ 5808818961aSMatthew G Knepley 581552f7358SJed Brown #undef __FUNCT__ 582552f7358SJed Brown #define __FUNCT__ "DMCreate_Plex" 5838cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 584552f7358SJed Brown { 585552f7358SJed Brown DM_Plex *mesh; 586770b213bSMatthew G Knepley PetscInt unit, d; 587552f7358SJed Brown PetscErrorCode ierr; 588552f7358SJed Brown 589552f7358SJed Brown PetscFunctionBegin; 590552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 591552f7358SJed Brown ierr = PetscNewLog(dm, DM_Plex, &mesh);CHKERRQ(ierr); 592552f7358SJed Brown dm->data = mesh; 593552f7358SJed Brown 594552f7358SJed Brown mesh->refct = 1; 595552f7358SJed Brown mesh->dim = 0; 59682f516ccSBarry Smith ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 597552f7358SJed Brown mesh->maxConeSize = 0; 5980298fd71SBarry Smith mesh->cones = NULL; 5990298fd71SBarry Smith mesh->coneOrientations = NULL; 60082f516ccSBarry Smith ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 601552f7358SJed Brown mesh->maxSupportSize = 0; 6020298fd71SBarry Smith mesh->supports = NULL; 603552f7358SJed Brown mesh->refinementUniform = PETSC_TRUE; 604552f7358SJed Brown mesh->refinementLimit = -1.0; 605552f7358SJed Brown 6060298fd71SBarry Smith mesh->facesTmp = NULL; 607552f7358SJed Brown 6080298fd71SBarry Smith mesh->subpointMap = NULL; 609552f7358SJed Brown 6108865f1eaSKarl Rupp for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 611552f7358SJed Brown 6120298fd71SBarry Smith mesh->labels = NULL; 6130298fd71SBarry Smith mesh->globalVertexNumbers = NULL; 6140298fd71SBarry Smith mesh->globalCellNumbers = NULL; 6158865f1eaSKarl Rupp for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 616552f7358SJed Brown mesh->vtkCellHeight = 0; 617552f7358SJed Brown mesh->preallocCenterDim = -1; 618552f7358SJed Brown 619552f7358SJed Brown mesh->printSetValues = PETSC_FALSE; 620552f7358SJed Brown mesh->printFEM = 0; 621*6113b454SMatthew G. Knepley mesh->printTol = 1.0e-10; 622552f7358SJed Brown 623552f7358SJed Brown ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 624552f7358SJed Brown PetscFunctionReturn(0); 625552f7358SJed Brown } 626552f7358SJed Brown 627552f7358SJed Brown #undef __FUNCT__ 628552f7358SJed Brown #define __FUNCT__ "DMPlexCreate" 629552f7358SJed Brown /*@ 630552f7358SJed Brown DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 631552f7358SJed Brown 632552f7358SJed Brown Collective on MPI_Comm 633552f7358SJed Brown 634552f7358SJed Brown Input Parameter: 635552f7358SJed Brown . comm - The communicator for the DMPlex object 636552f7358SJed Brown 637552f7358SJed Brown Output Parameter: 638552f7358SJed Brown . mesh - The DMPlex object 639552f7358SJed Brown 640552f7358SJed Brown Level: beginner 641552f7358SJed Brown 642552f7358SJed Brown .keywords: DMPlex, create 643552f7358SJed Brown @*/ 644552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 645552f7358SJed Brown { 646552f7358SJed Brown PetscErrorCode ierr; 647552f7358SJed Brown 648552f7358SJed Brown PetscFunctionBegin; 649552f7358SJed Brown PetscValidPointer(mesh,2); 650552f7358SJed Brown ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 651552f7358SJed Brown ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 652552f7358SJed Brown PetscFunctionReturn(0); 653552f7358SJed Brown } 654552f7358SJed Brown 655552f7358SJed Brown #undef __FUNCT__ 656552f7358SJed Brown #define __FUNCT__ "DMPlexClone" 657552f7358SJed Brown /*@ 658552f7358SJed Brown DMPlexClone - Creates a DMPlex object with the same mesh as the original. 659552f7358SJed Brown 660552f7358SJed Brown Collective on MPI_Comm 661552f7358SJed Brown 662552f7358SJed Brown Input Parameter: 663552f7358SJed Brown . dm - The original DMPlex object 664552f7358SJed Brown 665552f7358SJed Brown Output Parameter: 666552f7358SJed Brown . newdm - The new DMPlex object 667552f7358SJed Brown 668552f7358SJed Brown Level: beginner 669552f7358SJed Brown 670552f7358SJed Brown .keywords: DMPlex, create 671552f7358SJed Brown @*/ 672552f7358SJed Brown PetscErrorCode DMPlexClone(DM dm, DM *newdm) 673552f7358SJed Brown { 674552f7358SJed Brown DM_Plex *mesh; 6751c161d0cSMatthew G. Knepley Vec coords; 676552f7358SJed Brown void *ctx; 677552f7358SJed Brown PetscErrorCode ierr; 678552f7358SJed Brown 679552f7358SJed Brown PetscFunctionBegin; 680552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 681552f7358SJed Brown PetscValidPointer(newdm,2); 68282f516ccSBarry Smith ierr = DMCreate(PetscObjectComm((PetscObject)dm), newdm);CHKERRQ(ierr); 683552f7358SJed Brown ierr = PetscSFDestroy(&(*newdm)->sf);CHKERRQ(ierr); 684552f7358SJed Brown ierr = PetscObjectReference((PetscObject) dm->sf);CHKERRQ(ierr); 685552f7358SJed Brown (*newdm)->sf = dm->sf; 686552f7358SJed Brown mesh = (DM_Plex*) dm->data; 687552f7358SJed Brown mesh->refct++; 688552f7358SJed Brown (*newdm)->data = mesh; 689552f7358SJed Brown ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 690552f7358SJed Brown ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 691552f7358SJed Brown ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); 692552f7358SJed Brown ierr = DMSetApplicationContext(*newdm, ctx);CHKERRQ(ierr); 6931c161d0cSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coords);CHKERRQ(ierr); 6941c161d0cSMatthew G. Knepley if (coords) { 695f5e6d292SJed Brown ierr = DMSetCoordinatesLocal(*newdm, coords);CHKERRQ(ierr); 696f5e6d292SJed Brown } else { 6971c161d0cSMatthew G. Knepley ierr = DMGetCoordinates(dm, &coords);CHKERRQ(ierr); 698fb07c940SMatthew G. Knepley if (coords) {ierr = DMSetCoordinates(*newdm, coords);CHKERRQ(ierr);} 6991c161d0cSMatthew G. Knepley } 700552f7358SJed Brown PetscFunctionReturn(0); 701552f7358SJed Brown } 7029298eaa6SMatthew G Knepley 7039298eaa6SMatthew G Knepley #undef __FUNCT__ 7049298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildFromCellList_Private" 7059298eaa6SMatthew G Knepley /* 7069298eaa6SMatthew G Knepley This takes as input the common mesh generator output, a list of the vertices for each cell 7079298eaa6SMatthew G Knepley */ 7089298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 7099298eaa6SMatthew G Knepley { 7109298eaa6SMatthew G Knepley PetscInt *cone, c, p; 7119298eaa6SMatthew G Knepley PetscErrorCode ierr; 7129298eaa6SMatthew G Knepley 7139298eaa6SMatthew G Knepley PetscFunctionBegin; 7149298eaa6SMatthew G Knepley ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 7159298eaa6SMatthew G Knepley for (c = 0; c < numCells; ++c) { 7169298eaa6SMatthew G Knepley ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 7179298eaa6SMatthew G Knepley } 7189298eaa6SMatthew G Knepley ierr = DMSetUp(dm);CHKERRQ(ierr); 7199298eaa6SMatthew G Knepley ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 7209298eaa6SMatthew G Knepley for (c = 0; c < numCells; ++c) { 7219298eaa6SMatthew G Knepley for (p = 0; p < numCorners; ++p) { 7229298eaa6SMatthew G Knepley cone[p] = cells[c*numCorners+p]+numCells; 7239298eaa6SMatthew G Knepley } 7249298eaa6SMatthew G Knepley ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 7259298eaa6SMatthew G Knepley } 7269298eaa6SMatthew G Knepley ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 7279298eaa6SMatthew G Knepley ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 7289298eaa6SMatthew G Knepley ierr = DMPlexStratify(dm);CHKERRQ(ierr); 7299298eaa6SMatthew G Knepley PetscFunctionReturn(0); 7309298eaa6SMatthew G Knepley } 7319298eaa6SMatthew G Knepley 7329298eaa6SMatthew G Knepley #undef __FUNCT__ 7339298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildCoordinates_Private" 7349298eaa6SMatthew G Knepley /* 7359298eaa6SMatthew G Knepley This takes as input the coordinates for each vertex 7369298eaa6SMatthew G Knepley */ 7379298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 7389298eaa6SMatthew G Knepley { 7399298eaa6SMatthew G Knepley PetscSection coordSection; 7409298eaa6SMatthew G Knepley Vec coordinates; 7419298eaa6SMatthew G Knepley PetscScalar *coords; 7429298eaa6SMatthew G Knepley PetscInt coordSize, v, d; 7439298eaa6SMatthew G Knepley PetscErrorCode ierr; 7449298eaa6SMatthew G Knepley 7459298eaa6SMatthew G Knepley PetscFunctionBegin; 7469298eaa6SMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 7479298eaa6SMatthew G Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 7489298eaa6SMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 7499298eaa6SMatthew G Knepley ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 7509298eaa6SMatthew G Knepley for (v = numCells; v < numCells+numVertices; ++v) { 7519298eaa6SMatthew G Knepley ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 7529298eaa6SMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 7539298eaa6SMatthew G Knepley } 7549298eaa6SMatthew G Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 7559298eaa6SMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 7569298eaa6SMatthew G Knepley ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 7579298eaa6SMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 7589298eaa6SMatthew G Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 7599298eaa6SMatthew G Knepley ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 7609298eaa6SMatthew G Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 7619298eaa6SMatthew G Knepley for (v = 0; v < numVertices; ++v) { 7629298eaa6SMatthew G Knepley for (d = 0; d < spaceDim; ++d) { 7639298eaa6SMatthew G Knepley coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 7649298eaa6SMatthew G Knepley } 7659298eaa6SMatthew G Knepley } 7669298eaa6SMatthew G Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 7679298eaa6SMatthew G Knepley ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 7689298eaa6SMatthew G Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 7699298eaa6SMatthew G Knepley PetscFunctionReturn(0); 7709298eaa6SMatthew G Knepley } 7719298eaa6SMatthew G Knepley 7729298eaa6SMatthew G Knepley #undef __FUNCT__ 7739298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromCellList" 7749298eaa6SMatthew G Knepley /*@C 7759298eaa6SMatthew G Knepley DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 7769298eaa6SMatthew G Knepley 7779298eaa6SMatthew G Knepley Input Parameters: 7789298eaa6SMatthew G Knepley + comm - The communicator 7799298eaa6SMatthew G Knepley . dim - The topological dimension of the mesh 7809298eaa6SMatthew G Knepley . numCells - The number of cells 7819298eaa6SMatthew G Knepley . numVertices - The number of vertices 7829298eaa6SMatthew G Knepley . numCorners - The number of vertices for each cell 7839298eaa6SMatthew G Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 7849298eaa6SMatthew G Knepley . cells - An array of numCells*numCorners numbers, the vertices for each cell 7859298eaa6SMatthew G Knepley . spaceDim - The spatial dimension used for coordinates 7869298eaa6SMatthew G Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 7879298eaa6SMatthew G Knepley 7889298eaa6SMatthew G Knepley Output Parameter: 7899298eaa6SMatthew G Knepley . dm - The DM 7909298eaa6SMatthew G Knepley 7919298eaa6SMatthew G Knepley Note: Two triangles sharing a face 7929298eaa6SMatthew G Knepley $ 7939298eaa6SMatthew G Knepley $ 2 7949298eaa6SMatthew G Knepley $ / | \ 7959298eaa6SMatthew G Knepley $ / | \ 7969298eaa6SMatthew G Knepley $ / | \ 7979298eaa6SMatthew G Knepley $ 0 0 | 1 3 7989298eaa6SMatthew G Knepley $ \ | / 7999298eaa6SMatthew G Knepley $ \ | / 8009298eaa6SMatthew G Knepley $ \ | / 8019298eaa6SMatthew G Knepley $ 1 8029298eaa6SMatthew G Knepley would have input 8039298eaa6SMatthew G Knepley $ numCells = 2, numVertices = 4 8049298eaa6SMatthew G Knepley $ cells = [0 1 2 1 3 2] 8059298eaa6SMatthew G Knepley $ 8069298eaa6SMatthew G Knepley which would result in the DMPlex 8079298eaa6SMatthew G Knepley $ 8089298eaa6SMatthew G Knepley $ 4 8099298eaa6SMatthew G Knepley $ / | \ 8109298eaa6SMatthew G Knepley $ / | \ 8119298eaa6SMatthew G Knepley $ / | \ 8129298eaa6SMatthew G Knepley $ 2 0 | 1 5 8139298eaa6SMatthew G Knepley $ \ | / 8149298eaa6SMatthew G Knepley $ \ | / 8159298eaa6SMatthew G Knepley $ \ | / 8169298eaa6SMatthew G Knepley $ 3 8179298eaa6SMatthew G Knepley 8189298eaa6SMatthew G Knepley Level: beginner 8199298eaa6SMatthew G Knepley 8209298eaa6SMatthew G Knepley .seealso: DMPlexCreate() 8219298eaa6SMatthew G Knepley @*/ 8229298eaa6SMatthew 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) 8239298eaa6SMatthew G Knepley { 8249298eaa6SMatthew G Knepley PetscErrorCode ierr; 8259298eaa6SMatthew G Knepley 8269298eaa6SMatthew G Knepley PetscFunctionBegin; 8279298eaa6SMatthew G Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 8289298eaa6SMatthew G Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 8299298eaa6SMatthew G Knepley ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 8309298eaa6SMatthew G Knepley ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 8319298eaa6SMatthew G Knepley if (interpolate) { 8329298eaa6SMatthew G Knepley DM idm; 8339298eaa6SMatthew G Knepley 8349298eaa6SMatthew G Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 8359298eaa6SMatthew G Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 8369298eaa6SMatthew G Knepley *dm = idm; 8379298eaa6SMatthew G Knepley } 8389298eaa6SMatthew G Knepley ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 8399298eaa6SMatthew G Knepley PetscFunctionReturn(0); 8409298eaa6SMatthew G Knepley } 8419298eaa6SMatthew G Knepley 8429298eaa6SMatthew G Knepley #undef __FUNCT__ 8439298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromDAG" 8449298eaa6SMatthew G Knepley /* 8459298eaa6SMatthew G Knepley This takes as input the raw Hasse Diagram data 8469298eaa6SMatthew G Knepley */ 8479298eaa6SMatthew G Knepley PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 8489298eaa6SMatthew G Knepley { 8499298eaa6SMatthew G Knepley Vec coordinates; 8509298eaa6SMatthew G Knepley PetscSection coordSection; 8519298eaa6SMatthew G Knepley PetscScalar *coords; 8529298eaa6SMatthew G Knepley PetscInt coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off; 8539298eaa6SMatthew G Knepley PetscErrorCode ierr; 8549298eaa6SMatthew G Knepley 8559298eaa6SMatthew G Knepley PetscFunctionBegin; 8569298eaa6SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 8579298eaa6SMatthew G Knepley for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 8589298eaa6SMatthew G Knepley ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 8599298eaa6SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 8609298eaa6SMatthew G Knepley ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 8619298eaa6SMatthew G Knepley } 8629298eaa6SMatthew G Knepley ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 8639298eaa6SMatthew G Knepley for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 8649298eaa6SMatthew G Knepley ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 8659298eaa6SMatthew G Knepley ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 8669298eaa6SMatthew G Knepley } 8679298eaa6SMatthew G Knepley ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 8689298eaa6SMatthew G Knepley ierr = DMPlexStratify(dm);CHKERRQ(ierr); 8699298eaa6SMatthew G Knepley /* Build coordinates */ 8709298eaa6SMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 8719298eaa6SMatthew G Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 8729298eaa6SMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 8739298eaa6SMatthew G Knepley ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 8749298eaa6SMatthew G Knepley for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 8759298eaa6SMatthew G Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 8766f8cbbeeSMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 8779298eaa6SMatthew G Knepley } 8789298eaa6SMatthew G Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 8799298eaa6SMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 8809298eaa6SMatthew G Knepley ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 8816f8cbbeeSMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 8829298eaa6SMatthew G Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 8839298eaa6SMatthew G Knepley ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 8849298eaa6SMatthew G Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 8859298eaa6SMatthew G Knepley for (v = 0; v < numPoints[0]; ++v) { 8869298eaa6SMatthew G Knepley PetscInt off; 8879298eaa6SMatthew G Knepley 8889298eaa6SMatthew G Knepley ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 8899298eaa6SMatthew G Knepley for (d = 0; d < dim; ++d) { 8909298eaa6SMatthew G Knepley coords[off+d] = vertexCoords[v*dim+d]; 8919298eaa6SMatthew G Knepley } 8929298eaa6SMatthew G Knepley } 8939298eaa6SMatthew G Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 8949298eaa6SMatthew G Knepley ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 8959298eaa6SMatthew G Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 8969298eaa6SMatthew G Knepley PetscFunctionReturn(0); 8979298eaa6SMatthew G Knepley } 898