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); 216113b454SMatthew 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__ 271df5d5c5SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateDoublet" 281df5d5c5SMatthew G. Knepley /*@ 291df5d5c5SMatthew G. Knepley DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement. 301df5d5c5SMatthew G. Knepley 311df5d5c5SMatthew G. Knepley Collective on MPI_Comm 321df5d5c5SMatthew G. Knepley 331df5d5c5SMatthew G. Knepley Input Parameters: 341df5d5c5SMatthew G. Knepley + comm - The communicator for the DM object 351df5d5c5SMatthew G. Knepley . dim - The spatial dimension 361df5d5c5SMatthew G. Knepley . simplex - Flag for simplicial cells, otherwise they are tensor product cells 371df5d5c5SMatthew G. Knepley . interpolate - Flag to create intermediate mesh pieces (edges, faces) 381df5d5c5SMatthew G. Knepley . refinementUniform - Flag for uniform parallel refinement 391df5d5c5SMatthew G. Knepley - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell 401df5d5c5SMatthew G. Knepley 411df5d5c5SMatthew G. Knepley Output Parameter: 421df5d5c5SMatthew G. Knepley . dm - The DM object 431df5d5c5SMatthew G. Knepley 441df5d5c5SMatthew G. Knepley Level: beginner 451df5d5c5SMatthew G. Knepley 461df5d5c5SMatthew G. Knepley .keywords: DM, create 471df5d5c5SMatthew G. Knepley .seealso: DMSetType(), DMCreate() 481df5d5c5SMatthew G. Knepley @*/ 491df5d5c5SMatthew G. Knepley PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm) 501df5d5c5SMatthew G. Knepley { 511df5d5c5SMatthew G. Knepley DM dm; 521df5d5c5SMatthew G. Knepley PetscInt p; 531df5d5c5SMatthew G. Knepley PetscMPIInt rank; 541df5d5c5SMatthew G. Knepley PetscErrorCode ierr; 551df5d5c5SMatthew G. Knepley 561df5d5c5SMatthew G. Knepley PetscFunctionBegin; 571df5d5c5SMatthew G. Knepley ierr = DMCreate(comm, &dm);CHKERRQ(ierr); 581df5d5c5SMatthew G. Knepley ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr); 591df5d5c5SMatthew G. Knepley ierr = DMPlexSetDimension(dm, dim);CHKERRQ(ierr); 601df5d5c5SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 61ce78fa2fSMatthew G. Knepley switch (dim) { 62ce78fa2fSMatthew G. Knepley case 2: 63ce78fa2fSMatthew G. Knepley if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);} 64ce78fa2fSMatthew G. Knepley else {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);} 65ce78fa2fSMatthew G. Knepley break; 66ce78fa2fSMatthew G. Knepley case 3: 67ce78fa2fSMatthew G. Knepley if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);} 68ce78fa2fSMatthew G. Knepley else {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);} 69ce78fa2fSMatthew G. Knepley break; 70ce78fa2fSMatthew G. Knepley default: 71ce78fa2fSMatthew G. Knepley SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 72ce78fa2fSMatthew G. Knepley } 731df5d5c5SMatthew G. Knepley if (rank) { 741df5d5c5SMatthew G. Knepley PetscInt numPoints[2] = {0, 0}; 751df5d5c5SMatthew G. Knepley ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 761df5d5c5SMatthew G. Knepley } else { 771df5d5c5SMatthew G. Knepley switch (dim) { 781df5d5c5SMatthew G. Knepley case 2: 791df5d5c5SMatthew G. Knepley if (simplex) { 801df5d5c5SMatthew G. Knepley PetscInt numPoints[2] = {4, 2}; 811df5d5c5SMatthew G. Knepley PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0}; 821df5d5c5SMatthew G. Knepley PetscInt cones[6] = {2, 3, 4, 5, 4, 3}; 831df5d5c5SMatthew G. Knepley PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 841df5d5c5SMatthew G. Knepley PetscScalar vertexCoords[8] = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5}; 851df5d5c5SMatthew G. Knepley PetscInt markerPoints[8] = {2, 1, 3, 1, 4, 1, 5, 1}; 861df5d5c5SMatthew G. Knepley 871df5d5c5SMatthew G. Knepley ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 881df5d5c5SMatthew G. Knepley for (p = 0; p < 4; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);} 891df5d5c5SMatthew G. Knepley } else { 901df5d5c5SMatthew G. Knepley PetscInt numPoints[2] = {6, 2}; 911df5d5c5SMatthew G. Knepley PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 921df5d5c5SMatthew G. Knepley PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4}; 931df5d5c5SMatthew G. Knepley PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 941df5d5c5SMatthew G. Knepley PetscScalar vertexCoords[12] = {-1.0, -0.5, 0.0, -0.5, 0.0, 0.5, -1.0, 0.5, 1.0, -0.5, 1.0, 0.5}; 951df5d5c5SMatthew G. Knepley 961df5d5c5SMatthew G. Knepley ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 971df5d5c5SMatthew G. Knepley } 981df5d5c5SMatthew G. Knepley break; 991df5d5c5SMatthew G. Knepley case 3: 1001df5d5c5SMatthew G. Knepley if (simplex) { 1011df5d5c5SMatthew G. Knepley PetscInt numPoints[2] = {5, 2}; 1021df5d5c5SMatthew G. Knepley PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 1031df5d5c5SMatthew G. Knepley PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 1041df5d5c5SMatthew G. Knepley PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 1051df5d5c5SMatthew G. Knepley PetscScalar vertexCoords[15] = {-1.0, 0.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0}; 1061df5d5c5SMatthew G. Knepley PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 1071df5d5c5SMatthew G. Knepley 1081df5d5c5SMatthew G. Knepley ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 1091df5d5c5SMatthew G. Knepley for (p = 0; p < 5; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);} 1101df5d5c5SMatthew G. Knepley } else { 1111df5d5c5SMatthew G. Knepley PetscInt numPoints[2] = {12, 2}; 1121df5d5c5SMatthew G. Knepley PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 1131df5d5c5SMatthew G. Knepley PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 1141df5d5c5SMatthew G. Knepley PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 1151df5d5c5SMatthew G. Knepley PetscScalar vertexCoords[36] = {-1.0, -0.5, -0.5, -1.0, 0.5, -0.5, 0.0, 0.5, -0.5, 0.0, -0.5, -0.5, 1161df5d5c5SMatthew G. Knepley -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 1171df5d5c5SMatthew G. Knepley 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5}; 1181df5d5c5SMatthew G. Knepley 1191df5d5c5SMatthew G. Knepley ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 1201df5d5c5SMatthew G. Knepley } 1211df5d5c5SMatthew G. Knepley break; 1221df5d5c5SMatthew G. Knepley default: 1231df5d5c5SMatthew G. Knepley SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 1241df5d5c5SMatthew G. Knepley } 1251df5d5c5SMatthew G. Knepley } 1261df5d5c5SMatthew G. Knepley *newdm = dm; 1271df5d5c5SMatthew G. Knepley if (refinementLimit > 0.0) { 1281df5d5c5SMatthew G. Knepley DM rdm; 1291df5d5c5SMatthew G. Knepley const char *name; 1301df5d5c5SMatthew G. Knepley 1311df5d5c5SMatthew G. Knepley ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr); 1321df5d5c5SMatthew G. Knepley ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr); 1331df5d5c5SMatthew G. Knepley ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr); 1341df5d5c5SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr); 1351df5d5c5SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) rdm, name);CHKERRQ(ierr); 1361df5d5c5SMatthew G. Knepley ierr = DMDestroy(newdm);CHKERRQ(ierr); 1371df5d5c5SMatthew G. Knepley *newdm = rdm; 1381df5d5c5SMatthew G. Knepley } 1391df5d5c5SMatthew G. Knepley if (interpolate) { 1401df5d5c5SMatthew G. Knepley DM idm; 1411df5d5c5SMatthew G. Knepley const char *name; 1421df5d5c5SMatthew G. Knepley 1431df5d5c5SMatthew G. Knepley ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr); 1441df5d5c5SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr); 1451df5d5c5SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) idm, name);CHKERRQ(ierr); 1461df5d5c5SMatthew G. Knepley ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr); 1471df5d5c5SMatthew G. Knepley ierr = DMPlexCopyLabels(*newdm, idm);CHKERRQ(ierr); 1481df5d5c5SMatthew G. Knepley ierr = DMDestroy(newdm);CHKERRQ(ierr); 1491df5d5c5SMatthew G. Knepley *newdm = idm; 1501df5d5c5SMatthew G. Knepley } 1511df5d5c5SMatthew G. Knepley { 1521df5d5c5SMatthew G. Knepley DM refinedMesh = NULL; 1531df5d5c5SMatthew G. Knepley DM distributedMesh = NULL; 1541df5d5c5SMatthew G. Knepley 1551df5d5c5SMatthew G. Knepley /* Distribute mesh over processes */ 15607104863SJed Brown ierr = DMPlexDistribute(*newdm, NULL, 0, NULL, &distributedMesh);CHKERRQ(ierr); 1571df5d5c5SMatthew G. Knepley if (distributedMesh) { 1581df5d5c5SMatthew G. Knepley ierr = DMDestroy(newdm);CHKERRQ(ierr); 1591df5d5c5SMatthew G. Knepley *newdm = distributedMesh; 1601df5d5c5SMatthew G. Knepley } 1611df5d5c5SMatthew G. Knepley if (refinementUniform) { 1621df5d5c5SMatthew G. Knepley ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr); 1631df5d5c5SMatthew G. Knepley ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr); 1641df5d5c5SMatthew G. Knepley if (refinedMesh) { 1651df5d5c5SMatthew G. Knepley ierr = DMDestroy(newdm);CHKERRQ(ierr); 1661df5d5c5SMatthew G. Knepley *newdm = refinedMesh; 1671df5d5c5SMatthew G. Knepley } 1681df5d5c5SMatthew G. Knepley } 1691df5d5c5SMatthew G. Knepley } 1701df5d5c5SMatthew G. Knepley PetscFunctionReturn(0); 1711df5d5c5SMatthew G. Knepley } 1721df5d5c5SMatthew G. Knepley 1731df5d5c5SMatthew G. Knepley #undef __FUNCT__ 174552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareBoundary" 175552f7358SJed Brown /* 176552f7358SJed Brown Simple square boundary: 177552f7358SJed Brown 178552f7358SJed Brown 18--5-17--4--16 179552f7358SJed Brown | | | 180552f7358SJed Brown 6 10 3 181552f7358SJed Brown | | | 182552f7358SJed Brown 19-11-20--9--15 183552f7358SJed Brown | | | 184552f7358SJed Brown 7 8 2 185552f7358SJed Brown | | | 186552f7358SJed Brown 12--0-13--1--14 187552f7358SJed Brown */ 188552f7358SJed Brown PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 189552f7358SJed Brown { 190552f7358SJed Brown PetscInt numVertices = (edges[0]+1)*(edges[1]+1); 191552f7358SJed Brown PetscInt numEdges = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1]; 192552f7358SJed Brown PetscInt markerTop = 1; 193552f7358SJed Brown PetscInt markerBottom = 1; 194552f7358SJed Brown PetscInt markerRight = 1; 195552f7358SJed Brown PetscInt markerLeft = 1; 196552f7358SJed Brown PetscBool markerSeparate = PETSC_FALSE; 197552f7358SJed Brown Vec coordinates; 198552f7358SJed Brown PetscSection coordSection; 199552f7358SJed Brown PetscScalar *coords; 200552f7358SJed Brown PetscInt coordSize; 201552f7358SJed Brown PetscMPIInt rank; 202552f7358SJed Brown PetscInt v, vx, vy; 203552f7358SJed Brown PetscErrorCode ierr; 204552f7358SJed Brown 205552f7358SJed Brown PetscFunctionBegin; 2060298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 207552f7358SJed Brown if (markerSeparate) { 208552f7358SJed Brown markerTop = 1; 209552f7358SJed Brown markerBottom = 0; 210552f7358SJed Brown markerRight = 0; 211552f7358SJed Brown markerLeft = 0; 212552f7358SJed Brown } 21382f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 214552f7358SJed Brown if (!rank) { 215552f7358SJed Brown PetscInt e, ex, ey; 216552f7358SJed Brown 217552f7358SJed Brown ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr); 218552f7358SJed Brown for (e = 0; e < numEdges; ++e) { 219552f7358SJed Brown ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 220552f7358SJed Brown } 221552f7358SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 222552f7358SJed Brown for (vx = 0; vx <= edges[0]; vx++) { 223552f7358SJed Brown for (ey = 0; ey < edges[1]; ey++) { 224552f7358SJed Brown PetscInt edge = vx*edges[1] + ey + edges[0]*(edges[1]+1); 225552f7358SJed Brown PetscInt vertex = ey*(edges[0]+1) + vx + numEdges; 226da80777bSKarl Rupp PetscInt cone[2]; 227552f7358SJed Brown 228da80777bSKarl Rupp cone[0] = vertex; cone[1] = vertex+edges[0]+1; 229552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 230552f7358SJed Brown if (vx == edges[0]) { 231552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 232552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 233552f7358SJed Brown if (ey == edges[1]-1) { 234552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 235552f7358SJed Brown } 236552f7358SJed Brown } else if (vx == 0) { 237552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 238552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 239552f7358SJed Brown if (ey == edges[1]-1) { 240552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 241552f7358SJed Brown } 242552f7358SJed Brown } 243552f7358SJed Brown } 244552f7358SJed Brown } 245552f7358SJed Brown for (vy = 0; vy <= edges[1]; vy++) { 246552f7358SJed Brown for (ex = 0; ex < edges[0]; ex++) { 247552f7358SJed Brown PetscInt edge = vy*edges[0] + ex; 248552f7358SJed Brown PetscInt vertex = vy*(edges[0]+1) + ex + numEdges; 249da80777bSKarl Rupp PetscInt cone[2]; 250552f7358SJed Brown 251da80777bSKarl Rupp cone[0] = vertex; cone[1] = vertex+1; 252552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 253552f7358SJed Brown if (vy == edges[1]) { 254552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 255552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 256552f7358SJed Brown if (ex == edges[0]-1) { 257552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 258552f7358SJed Brown } 259552f7358SJed Brown } else if (vy == 0) { 260552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 261552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 262552f7358SJed Brown if (ex == edges[0]-1) { 263552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 264552f7358SJed Brown } 265552f7358SJed Brown } 266552f7358SJed Brown } 267552f7358SJed Brown } 268552f7358SJed Brown } 269552f7358SJed Brown ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 270552f7358SJed Brown ierr = DMPlexStratify(dm);CHKERRQ(ierr); 271552f7358SJed Brown /* Build coordinates */ 272c2166f76SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 273552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr); 274552f7358SJed Brown for (v = numEdges; v < numEdges+numVertices; ++v) { 275552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 276552f7358SJed Brown } 277552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 278552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 27982f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 280552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2812eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 282552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 283552f7358SJed Brown for (vy = 0; vy <= edges[1]; ++vy) { 284552f7358SJed Brown for (vx = 0; vx <= edges[0]; ++vx) { 285552f7358SJed Brown coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx; 286552f7358SJed Brown coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy; 287552f7358SJed Brown } 288552f7358SJed Brown } 289552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 290552f7358SJed Brown ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 291552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 292552f7358SJed Brown PetscFunctionReturn(0); 293552f7358SJed Brown } 294552f7358SJed Brown 295552f7358SJed Brown #undef __FUNCT__ 296552f7358SJed Brown #define __FUNCT__ "DMPlexCreateCubeBoundary" 297552f7358SJed Brown /* 298552f7358SJed Brown Simple cubic boundary: 299552f7358SJed Brown 300552f7358SJed Brown 2-------3 301552f7358SJed Brown /| /| 302552f7358SJed Brown 6-------7 | 303552f7358SJed Brown | | | | 304552f7358SJed Brown | 0-----|-1 305552f7358SJed Brown |/ |/ 306552f7358SJed Brown 4-------5 307552f7358SJed Brown */ 308552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[]) 309552f7358SJed Brown { 310552f7358SJed Brown PetscInt numVertices = (faces[0]+1)*(faces[1]+1)*(faces[2]+1); 311552f7358SJed Brown PetscInt numFaces = 6; 312552f7358SJed Brown Vec coordinates; 313552f7358SJed Brown PetscSection coordSection; 314552f7358SJed Brown PetscScalar *coords; 315552f7358SJed Brown PetscInt coordSize; 316552f7358SJed Brown PetscMPIInt rank; 317552f7358SJed Brown PetscInt v, vx, vy, vz; 318552f7358SJed Brown PetscErrorCode ierr; 319552f7358SJed Brown 320552f7358SJed Brown PetscFunctionBegin; 32182f516ccSBarry 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"); 32282f516ccSBarry 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"); 32382f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 324552f7358SJed Brown if (!rank) { 325552f7358SJed Brown PetscInt f; 326552f7358SJed Brown 327552f7358SJed Brown ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr); 328552f7358SJed Brown for (f = 0; f < numFaces; ++f) { 329552f7358SJed Brown ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 330552f7358SJed Brown } 331552f7358SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 332552f7358SJed Brown for (v = 0; v < numFaces+numVertices; ++v) { 333552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr); 334552f7358SJed Brown } 335552f7358SJed Brown { /* Side 0 (Front) */ 336da80777bSKarl Rupp PetscInt cone[4]; 337da80777bSKarl Rupp cone[0] = numFaces+4; cone[1] = numFaces+5; cone[2] = numFaces+7; cone[3] = numFaces+6; 338552f7358SJed Brown ierr = DMPlexSetCone(dm, 0, cone);CHKERRQ(ierr); 339552f7358SJed Brown } 340552f7358SJed Brown { /* Side 1 (Back) */ 341da80777bSKarl Rupp PetscInt cone[4]; 342da80777bSKarl Rupp cone[0] = numFaces+1; cone[1] = numFaces+0; cone[2] = numFaces+2; cone[3] = numFaces+3; 343552f7358SJed Brown ierr = DMPlexSetCone(dm, 1, cone);CHKERRQ(ierr); 344552f7358SJed Brown } 345552f7358SJed Brown { /* Side 2 (Bottom) */ 346da80777bSKarl Rupp PetscInt cone[4]; 347da80777bSKarl Rupp cone[0] = numFaces+0; cone[1] = numFaces+1; cone[2] = numFaces+5; cone[3] = numFaces+4; 348552f7358SJed Brown ierr = DMPlexSetCone(dm, 2, cone);CHKERRQ(ierr); 349552f7358SJed Brown } 350552f7358SJed Brown { /* Side 3 (Top) */ 351da80777bSKarl Rupp PetscInt cone[4]; 352da80777bSKarl Rupp cone[0] = numFaces+6; cone[1] = numFaces+7; cone[2] = numFaces+3; cone[3] = numFaces+2; 353552f7358SJed Brown ierr = DMPlexSetCone(dm, 3, cone);CHKERRQ(ierr); 354552f7358SJed Brown } 355552f7358SJed Brown { /* Side 4 (Left) */ 356da80777bSKarl Rupp PetscInt cone[4]; 357da80777bSKarl Rupp cone[0] = numFaces+0; cone[1] = numFaces+4; cone[2] = numFaces+6; cone[3] = numFaces+2; 358552f7358SJed Brown ierr = DMPlexSetCone(dm, 4, cone);CHKERRQ(ierr); 359552f7358SJed Brown } 360552f7358SJed Brown { /* Side 5 (Right) */ 361da80777bSKarl Rupp PetscInt cone[4]; 362da80777bSKarl Rupp cone[0] = numFaces+5; cone[1] = numFaces+1; cone[2] = numFaces+3; cone[3] = numFaces+7; 363552f7358SJed Brown ierr = DMPlexSetCone(dm, 5, cone);CHKERRQ(ierr); 364552f7358SJed Brown } 365552f7358SJed Brown } 366552f7358SJed Brown ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 367552f7358SJed Brown ierr = DMPlexStratify(dm);CHKERRQ(ierr); 368552f7358SJed Brown /* Build coordinates */ 369c2166f76SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 370552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr); 371552f7358SJed Brown for (v = numFaces; v < numFaces+numVertices; ++v) { 372552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr); 373552f7358SJed Brown } 374552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 375552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 37682f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 377552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 378552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 3792eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 380552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 381552f7358SJed Brown for (vz = 0; vz <= faces[2]; ++vz) { 382552f7358SJed Brown for (vy = 0; vy <= faces[1]; ++vy) { 383552f7358SJed Brown for (vx = 0; vx <= faces[0]; ++vx) { 384552f7358SJed Brown coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx; 385552f7358SJed Brown coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy; 386552f7358SJed Brown coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz; 387552f7358SJed Brown } 388552f7358SJed Brown } 389552f7358SJed Brown } 390552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 391552f7358SJed Brown ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 392552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 393552f7358SJed Brown PetscFunctionReturn(0); 394552f7358SJed Brown } 395552f7358SJed Brown 396552f7358SJed Brown #undef __FUNCT__ 397552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareMesh" 398552f7358SJed Brown /* 399552f7358SJed Brown Simple square mesh: 400552f7358SJed Brown 401552f7358SJed Brown 22--8-23--9--24 402552f7358SJed Brown | | | 403552f7358SJed Brown 13 2 14 3 15 404552f7358SJed Brown | | | 405552f7358SJed Brown 19--6-20--7--21 406552f7358SJed Brown | | | 407552f7358SJed Brown 10 0 11 1 12 408552f7358SJed Brown | | | 409552f7358SJed Brown 16--4-17--5--18 410552f7358SJed Brown */ 411bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY) 412552f7358SJed Brown { 413552f7358SJed Brown PetscInt markerTop = 1; 414552f7358SJed Brown PetscInt markerBottom = 1; 415552f7358SJed Brown PetscInt markerRight = 1; 416552f7358SJed Brown PetscInt markerLeft = 1; 417552f7358SJed Brown PetscBool markerSeparate = PETSC_FALSE; 418552f7358SJed Brown PetscMPIInt rank; 419552f7358SJed Brown PetscErrorCode ierr; 420552f7358SJed Brown 421552f7358SJed Brown PetscFunctionBegin; 42282f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 4230298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 424552f7358SJed Brown if (markerSeparate) { 425552f7358SJed Brown markerTop = 3; 426552f7358SJed Brown markerBottom = 1; 427552f7358SJed Brown markerRight = 2; 428552f7358SJed Brown markerLeft = 4; 429552f7358SJed Brown } 430552f7358SJed Brown { 431552f7358SJed Brown const PetscInt numXEdges = !rank ? edges[0] : 0; 432552f7358SJed Brown const PetscInt numYEdges = !rank ? edges[1] : 0; 433*46675b8eSMatthew G. Knepley const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0; 434*46675b8eSMatthew G. Knepley const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0; 435552f7358SJed Brown const PetscInt numTotXEdges = numXEdges*numYVertices; 436552f7358SJed Brown const PetscInt numTotYEdges = numYEdges*numXVertices; 437552f7358SJed Brown const PetscInt numVertices = numXVertices*numYVertices; 438552f7358SJed Brown const PetscInt numEdges = numTotXEdges + numTotYEdges; 439552f7358SJed Brown const PetscInt numFaces = numXEdges*numYEdges; 440552f7358SJed Brown const PetscInt firstVertex = numFaces; 441552f7358SJed Brown const PetscInt firstXEdge = numFaces + numVertices; 442552f7358SJed Brown const PetscInt firstYEdge = numFaces + numVertices + numTotXEdges; 443552f7358SJed Brown Vec coordinates; 444552f7358SJed Brown PetscSection coordSection; 445552f7358SJed Brown PetscScalar *coords; 446552f7358SJed Brown PetscInt coordSize; 447552f7358SJed Brown PetscInt v, vx, vy; 448552f7358SJed Brown PetscInt f, fx, fy, e, ex, ey; 449552f7358SJed Brown 450552f7358SJed Brown ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr); 451552f7358SJed Brown for (f = 0; f < numFaces; ++f) { 452552f7358SJed Brown ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 453552f7358SJed Brown } 454552f7358SJed Brown for (e = firstXEdge; e < firstXEdge+numEdges; ++e) { 455552f7358SJed Brown ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 456552f7358SJed Brown } 457552f7358SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 458552f7358SJed Brown /* Build faces */ 459552f7358SJed Brown for (fy = 0; fy < numYEdges; fy++) { 460552f7358SJed Brown for (fx = 0; fx < numXEdges; fx++) { 461552f7358SJed Brown const PetscInt face = fy*numXEdges + fx; 462552f7358SJed Brown const PetscInt edgeL = firstYEdge + fx *numYEdges + fy; 463becd5721SMatthew G. Knepley const PetscInt edgeR = firstYEdge + ((fx+1)%numXVertices)*numYEdges + fy; 464552f7358SJed Brown const PetscInt edgeB = firstXEdge + fy *numXEdges + fx; 465becd5721SMatthew G. Knepley const PetscInt edgeT = firstXEdge + ((fy+1)%numYVertices)*numXEdges + fx; 466552f7358SJed Brown const PetscInt ornt[4] = {0, 0, -2, -2}; 467da80777bSKarl Rupp PetscInt cone[4]; 468552f7358SJed Brown 469becd5721SMatthew G. Knepley cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 470552f7358SJed Brown ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 471552f7358SJed Brown ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 472552f7358SJed Brown } 473552f7358SJed Brown } 474552f7358SJed Brown /* Build Y edges*/ 475552f7358SJed Brown for (vx = 0; vx < numXVertices; vx++) { 476552f7358SJed Brown for (ey = 0; ey < numYEdges; ey++) { 477*46675b8eSMatthew G. Knepley const PetscInt nextv = (bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : ((ey+1)%numYVertices)*numXVertices + vx; 478552f7358SJed Brown const PetscInt edge = firstYEdge + vx*numYEdges + ey; 479becd5721SMatthew G. Knepley const PetscInt vertexB = firstVertex + ey*numXVertices + vx; 480*46675b8eSMatthew G. Knepley const PetscInt vertexT = firstVertex + nextv; 481da80777bSKarl Rupp PetscInt cone[2]; 482552f7358SJed Brown 483becd5721SMatthew G. Knepley cone[0] = vertexB; cone[1] = vertexT; 484552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 485552f7358SJed Brown if (vx == numXVertices-1) { 486552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 487552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 488552f7358SJed Brown if (ey == numYEdges-1) { 489552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 490552f7358SJed Brown } 491552f7358SJed Brown } else if (vx == 0) { 492552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 493552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 494552f7358SJed Brown if (ey == numYEdges-1) { 495552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 496552f7358SJed Brown } 497552f7358SJed Brown } 498552f7358SJed Brown } 499552f7358SJed Brown } 500552f7358SJed Brown /* Build X edges*/ 501552f7358SJed Brown for (vy = 0; vy < numYVertices; vy++) { 502552f7358SJed Brown for (ex = 0; ex < numXEdges; ex++) { 503*46675b8eSMatthew G. Knepley const PetscInt nextv = (bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : vy*numXVertices + (ex+1)%numXVertices; 504552f7358SJed Brown const PetscInt edge = firstXEdge + vy*numXEdges + ex; 505becd5721SMatthew G. Knepley const PetscInt vertexL = firstVertex + vy*numXVertices + ex; 506*46675b8eSMatthew G. Knepley const PetscInt vertexR = firstVertex + nextv; 507da80777bSKarl Rupp PetscInt cone[2]; 508552f7358SJed Brown 509becd5721SMatthew G. Knepley cone[0] = vertexL; cone[1] = vertexR; 510552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 511552f7358SJed Brown if (vy == numYVertices-1) { 512552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 513552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 514552f7358SJed Brown if (ex == numXEdges-1) { 515552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 516552f7358SJed Brown } 517552f7358SJed Brown } else if (vy == 0) { 518552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 519552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 520552f7358SJed Brown if (ex == numXEdges-1) { 521552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 522552f7358SJed Brown } 523552f7358SJed Brown } 524552f7358SJed Brown } 525552f7358SJed Brown } 526552f7358SJed Brown ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 527552f7358SJed Brown ierr = DMPlexStratify(dm);CHKERRQ(ierr); 528552f7358SJed Brown /* Build coordinates */ 529c2166f76SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 530552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr); 531552f7358SJed Brown for (v = firstVertex; v < firstVertex+numVertices; ++v) { 532552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 533552f7358SJed Brown } 534552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 535552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 53682f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 537552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 5382eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 539552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 540552f7358SJed Brown for (vy = 0; vy < numYVertices; ++vy) { 541552f7358SJed Brown for (vx = 0; vx < numXVertices; ++vx) { 542552f7358SJed Brown coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx; 543552f7358SJed Brown coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy; 544552f7358SJed Brown } 545552f7358SJed Brown } 546552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 547552f7358SJed Brown ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 548552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 549552f7358SJed Brown } 550552f7358SJed Brown PetscFunctionReturn(0); 551552f7358SJed Brown } 552552f7358SJed Brown 553552f7358SJed Brown #undef __FUNCT__ 554552f7358SJed Brown #define __FUNCT__ "DMPlexCreateBoxMesh" 5556b75479cSMatthew G Knepley /*@ 5566b75479cSMatthew G Knepley DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box). 5576b75479cSMatthew G Knepley 5586b75479cSMatthew G Knepley Collective on MPI_Comm 5596b75479cSMatthew G Knepley 5606b75479cSMatthew G Knepley Input Parameters: 5616b75479cSMatthew G Knepley + comm - The communicator for the DM object 5626b75479cSMatthew G Knepley . dim - The spatial dimension 5636b75479cSMatthew G Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces) 5646b75479cSMatthew G Knepley 5656b75479cSMatthew G Knepley Output Parameter: 5666b75479cSMatthew G Knepley . dm - The DM object 5676b75479cSMatthew G Knepley 5686b75479cSMatthew G Knepley Level: beginner 5696b75479cSMatthew G Knepley 5706b75479cSMatthew G Knepley .keywords: DM, create 5716b75479cSMatthew G Knepley .seealso: DMSetType(), DMCreate() 5726b75479cSMatthew G Knepley @*/ 5730adebc6cSBarry Smith PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm) 5740adebc6cSBarry Smith { 575552f7358SJed Brown DM boundary; 576552f7358SJed Brown PetscErrorCode ierr; 577552f7358SJed Brown 578552f7358SJed Brown PetscFunctionBegin; 579552f7358SJed Brown PetscValidPointer(dm, 4); 580552f7358SJed Brown ierr = DMCreate(comm, &boundary);CHKERRQ(ierr); 581552f7358SJed Brown PetscValidLogicalCollectiveInt(boundary,dim,2); 582552f7358SJed Brown ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr); 583552f7358SJed Brown ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr); 584552f7358SJed Brown switch (dim) { 585552f7358SJed Brown case 2: 586552f7358SJed Brown { 587552f7358SJed Brown PetscReal lower[2] = {0.0, 0.0}; 588552f7358SJed Brown PetscReal upper[2] = {1.0, 1.0}; 589552f7358SJed Brown PetscInt edges[2] = {2, 2}; 590552f7358SJed Brown 591552f7358SJed Brown ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr); 592552f7358SJed Brown break; 593552f7358SJed Brown } 594552f7358SJed Brown case 3: 595552f7358SJed Brown { 596552f7358SJed Brown PetscReal lower[3] = {0.0, 0.0, 0.0}; 597552f7358SJed Brown PetscReal upper[3] = {1.0, 1.0, 1.0}; 598552f7358SJed Brown PetscInt faces[3] = {1, 1, 1}; 599552f7358SJed Brown 600552f7358SJed Brown ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr); 601552f7358SJed Brown break; 602552f7358SJed Brown } 603552f7358SJed Brown default: 604552f7358SJed Brown SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 605552f7358SJed Brown } 6060298fd71SBarry Smith ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr); 607552f7358SJed Brown ierr = DMDestroy(&boundary);CHKERRQ(ierr); 608552f7358SJed Brown PetscFunctionReturn(0); 609552f7358SJed Brown } 610552f7358SJed Brown 611552f7358SJed Brown #undef __FUNCT__ 612552f7358SJed Brown #define __FUNCT__ "DMPlexCreateHexBoxMesh" 613bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm) 614a6dfd86eSKarl Rupp { 615552f7358SJed Brown PetscErrorCode ierr; 616552f7358SJed Brown 617552f7358SJed Brown PetscFunctionBegin; 618552f7358SJed Brown PetscValidPointer(dm, 4); 619552f7358SJed Brown ierr = DMCreate(comm, dm);CHKERRQ(ierr); 620552f7358SJed Brown PetscValidLogicalCollectiveInt(*dm,dim,2); 621552f7358SJed Brown ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 622552f7358SJed Brown ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 623552f7358SJed Brown switch (dim) { 624552f7358SJed Brown case 2: 625552f7358SJed Brown { 626552f7358SJed Brown PetscReal lower[2] = {0.0, 0.0}; 627552f7358SJed Brown PetscReal upper[2] = {1.0, 1.0}; 628552f7358SJed Brown 629becd5721SMatthew G. Knepley ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr); 630552f7358SJed Brown break; 631552f7358SJed Brown } 632552f7358SJed Brown #if 0 633552f7358SJed Brown case 3: 634552f7358SJed Brown { 635552f7358SJed Brown PetscReal lower[3] = {0.0, 0.0, 0.0}; 636552f7358SJed Brown PetscReal upper[3] = {1.0, 1.0, 1.0}; 637552f7358SJed Brown 638becd5721SMatthew G. Knepley ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr); 639552f7358SJed Brown break; 640552f7358SJed Brown } 641552f7358SJed Brown #endif 642552f7358SJed Brown default: 643552f7358SJed Brown SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 644552f7358SJed Brown } 645552f7358SJed Brown PetscFunctionReturn(0); 646552f7358SJed Brown } 647552f7358SJed Brown 648552f7358SJed Brown /* External function declarations here */ 649552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 650b412c318SBarry Smith extern PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J); 651552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 652552f7358SJed Brown extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined); 6539a922efaSMatthew G. Knepley extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 654552f7358SJed Brown extern PetscErrorCode DMSetUp_Plex(DM dm); 655552f7358SJed Brown extern PetscErrorCode DMDestroy_Plex(DM dm); 656552f7358SJed Brown extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 657552f7358SJed Brown extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 658552f7358SJed Brown extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS); 659552f7358SJed Brown 660552f7358SJed Brown #undef __FUNCT__ 661552f7358SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Plex" 662552f7358SJed Brown static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 663552f7358SJed Brown { 664552f7358SJed Brown PetscErrorCode ierr; 665552f7358SJed Brown 666552f7358SJed Brown PetscFunctionBegin; 667552f7358SJed Brown ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 668552f7358SJed Brown /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 669552f7358SJed Brown ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex);CHKERRQ(ierr); 670552f7358SJed Brown PetscFunctionReturn(0); 671552f7358SJed Brown } 672552f7358SJed Brown 673552f7358SJed Brown #undef __FUNCT__ 674552f7358SJed Brown #define __FUNCT__ "DMCreateLocalVector_Plex" 675552f7358SJed Brown static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 676552f7358SJed Brown { 677552f7358SJed Brown PetscErrorCode ierr; 678552f7358SJed Brown 679552f7358SJed Brown PetscFunctionBegin; 680552f7358SJed Brown ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 681552f7358SJed Brown ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 682552f7358SJed Brown PetscFunctionReturn(0); 683552f7358SJed Brown } 684552f7358SJed Brown 68538221697SMatthew G. Knepley #undef __FUNCT__ 686552f7358SJed Brown #define __FUNCT__ "DMInitialize_Plex" 687552f7358SJed Brown PetscErrorCode DMInitialize_Plex(DM dm) 688552f7358SJed Brown { 689552f7358SJed Brown PetscFunctionBegin; 690552f7358SJed Brown dm->ops->view = DMView_Plex; 691552f7358SJed Brown dm->ops->setfromoptions = DMSetFromOptions_Plex; 69238221697SMatthew G. Knepley dm->ops->clone = DMClone_Plex; 693552f7358SJed Brown dm->ops->setup = DMSetUp_Plex; 694552f7358SJed Brown dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 695552f7358SJed Brown dm->ops->createlocalvector = DMCreateLocalVector_Plex; 696184d77edSJed Brown dm->ops->getlocaltoglobalmapping = NULL; 697184d77edSJed Brown dm->ops->getlocaltoglobalmappingblock = NULL; 6980298fd71SBarry Smith dm->ops->createfieldis = NULL; 699552f7358SJed Brown dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 700552f7358SJed Brown dm->ops->getcoloring = 0; 701552f7358SJed Brown dm->ops->creatematrix = DMCreateMatrix_Plex; 702552f7358SJed Brown dm->ops->createinterpolation = 0; 703552f7358SJed Brown dm->ops->getaggregates = 0; 704552f7358SJed Brown dm->ops->getinjection = 0; 705552f7358SJed Brown dm->ops->refine = DMRefine_Plex; 706552f7358SJed Brown dm->ops->coarsen = 0; 707552f7358SJed Brown dm->ops->refinehierarchy = 0; 708552f7358SJed Brown dm->ops->coarsenhierarchy = 0; 7090298fd71SBarry Smith dm->ops->globaltolocalbegin = NULL; 7100298fd71SBarry Smith dm->ops->globaltolocalend = NULL; 7110298fd71SBarry Smith dm->ops->localtoglobalbegin = NULL; 7120298fd71SBarry Smith dm->ops->localtoglobalend = NULL; 713552f7358SJed Brown dm->ops->destroy = DMDestroy_Plex; 714552f7358SJed Brown dm->ops->createsubdm = DMCreateSubDM_Plex; 715552f7358SJed Brown dm->ops->locatepoints = DMLocatePoints_Plex; 716552f7358SJed Brown PetscFunctionReturn(0); 717552f7358SJed Brown } 718552f7358SJed Brown 71963a16f15SMatthew G. Knepley #undef __FUNCT__ 72063a16f15SMatthew G. Knepley #define __FUNCT__ "DMClone_Plex" 72163a16f15SMatthew G. Knepley PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 72263a16f15SMatthew G. Knepley { 72363a16f15SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 72463a16f15SMatthew G. Knepley PetscErrorCode ierr; 72563a16f15SMatthew G. Knepley 72663a16f15SMatthew G. Knepley PetscFunctionBegin; 72763a16f15SMatthew G. Knepley mesh->refct++; 72863a16f15SMatthew G. Knepley (*newdm)->data = mesh; 72963a16f15SMatthew G. Knepley ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 73063a16f15SMatthew G. Knepley ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 73163a16f15SMatthew G. Knepley PetscFunctionReturn(0); 73263a16f15SMatthew G. Knepley } 73363a16f15SMatthew G. Knepley 7348818961aSMatthew G Knepley /*MC 7358818961aSMatthew G Knepley DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 7368818961aSMatthew G Knepley In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 7378818961aSMatthew G Knepley specified by a PetscSection object. Ownership in the global representation is determined by 7388818961aSMatthew G Knepley ownership of the underlying DMPlex points. This is specified by another PetscSection object. 7398818961aSMatthew G Knepley 7408818961aSMatthew G Knepley Level: intermediate 7418818961aSMatthew G Knepley 7428818961aSMatthew G Knepley .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 7438818961aSMatthew G Knepley M*/ 7448818961aSMatthew G Knepley 745552f7358SJed Brown #undef __FUNCT__ 746552f7358SJed Brown #define __FUNCT__ "DMCreate_Plex" 7478cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 748552f7358SJed Brown { 749552f7358SJed Brown DM_Plex *mesh; 750770b213bSMatthew G Knepley PetscInt unit, d; 751552f7358SJed Brown PetscErrorCode ierr; 752552f7358SJed Brown 753552f7358SJed Brown PetscFunctionBegin; 754552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 755b00a9115SJed Brown ierr = PetscNewLog(dm,&mesh);CHKERRQ(ierr); 756552f7358SJed Brown dm->data = mesh; 757552f7358SJed Brown 758552f7358SJed Brown mesh->refct = 1; 759552f7358SJed Brown mesh->dim = 0; 76082f516ccSBarry Smith ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 761552f7358SJed Brown mesh->maxConeSize = 0; 7620298fd71SBarry Smith mesh->cones = NULL; 7630298fd71SBarry Smith mesh->coneOrientations = NULL; 76482f516ccSBarry Smith ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 765552f7358SJed Brown mesh->maxSupportSize = 0; 7660298fd71SBarry Smith mesh->supports = NULL; 767552f7358SJed Brown mesh->refinementUniform = PETSC_TRUE; 768552f7358SJed Brown mesh->refinementLimit = -1.0; 769552f7358SJed Brown 7700298fd71SBarry Smith mesh->facesTmp = NULL; 771552f7358SJed Brown 7720298fd71SBarry Smith mesh->subpointMap = NULL; 773552f7358SJed Brown 7748865f1eaSKarl Rupp for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 775552f7358SJed Brown 7760298fd71SBarry Smith mesh->labels = NULL; 7778af19771SMatthew G. Knepley mesh->depthLabel = NULL; 7780298fd71SBarry Smith mesh->globalVertexNumbers = NULL; 7790298fd71SBarry Smith mesh->globalCellNumbers = NULL; 7808865f1eaSKarl Rupp for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 781552f7358SJed Brown mesh->vtkCellHeight = 0; 782552f7358SJed Brown mesh->preallocCenterDim = -1; 783552f7358SJed Brown 784552f7358SJed Brown mesh->printSetValues = PETSC_FALSE; 785552f7358SJed Brown mesh->printFEM = 0; 7866113b454SMatthew G. Knepley mesh->printTol = 1.0e-10; 787552f7358SJed Brown 788552f7358SJed Brown ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 789552f7358SJed Brown PetscFunctionReturn(0); 790552f7358SJed Brown } 791552f7358SJed Brown 792552f7358SJed Brown #undef __FUNCT__ 793552f7358SJed Brown #define __FUNCT__ "DMPlexCreate" 794552f7358SJed Brown /*@ 795552f7358SJed Brown DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 796552f7358SJed Brown 797552f7358SJed Brown Collective on MPI_Comm 798552f7358SJed Brown 799552f7358SJed Brown Input Parameter: 800552f7358SJed Brown . comm - The communicator for the DMPlex object 801552f7358SJed Brown 802552f7358SJed Brown Output Parameter: 803552f7358SJed Brown . mesh - The DMPlex object 804552f7358SJed Brown 805552f7358SJed Brown Level: beginner 806552f7358SJed Brown 807552f7358SJed Brown .keywords: DMPlex, create 808552f7358SJed Brown @*/ 809552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 810552f7358SJed Brown { 811552f7358SJed Brown PetscErrorCode ierr; 812552f7358SJed Brown 813552f7358SJed Brown PetscFunctionBegin; 814552f7358SJed Brown PetscValidPointer(mesh,2); 815552f7358SJed Brown ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 816552f7358SJed Brown ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 817552f7358SJed Brown PetscFunctionReturn(0); 818552f7358SJed Brown } 819552f7358SJed Brown 820552f7358SJed Brown #undef __FUNCT__ 8219298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildFromCellList_Private" 8229298eaa6SMatthew G Knepley /* 8239298eaa6SMatthew G Knepley This takes as input the common mesh generator output, a list of the vertices for each cell 8249298eaa6SMatthew G Knepley */ 8259298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 8269298eaa6SMatthew G Knepley { 8279298eaa6SMatthew G Knepley PetscInt *cone, c, p; 8289298eaa6SMatthew G Knepley PetscErrorCode ierr; 8299298eaa6SMatthew G Knepley 8309298eaa6SMatthew G Knepley PetscFunctionBegin; 8319298eaa6SMatthew G Knepley ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 8329298eaa6SMatthew G Knepley for (c = 0; c < numCells; ++c) { 8339298eaa6SMatthew G Knepley ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 8349298eaa6SMatthew G Knepley } 8359298eaa6SMatthew G Knepley ierr = DMSetUp(dm);CHKERRQ(ierr); 8369298eaa6SMatthew G Knepley ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 8379298eaa6SMatthew G Knepley for (c = 0; c < numCells; ++c) { 8389298eaa6SMatthew G Knepley for (p = 0; p < numCorners; ++p) { 8399298eaa6SMatthew G Knepley cone[p] = cells[c*numCorners+p]+numCells; 8409298eaa6SMatthew G Knepley } 8419298eaa6SMatthew G Knepley ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 8429298eaa6SMatthew G Knepley } 8439298eaa6SMatthew G Knepley ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 8449298eaa6SMatthew G Knepley ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 8459298eaa6SMatthew G Knepley ierr = DMPlexStratify(dm);CHKERRQ(ierr); 8469298eaa6SMatthew G Knepley PetscFunctionReturn(0); 8479298eaa6SMatthew G Knepley } 8489298eaa6SMatthew G Knepley 8499298eaa6SMatthew G Knepley #undef __FUNCT__ 8509298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildCoordinates_Private" 8519298eaa6SMatthew G Knepley /* 8529298eaa6SMatthew G Knepley This takes as input the coordinates for each vertex 8539298eaa6SMatthew G Knepley */ 8549298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 8559298eaa6SMatthew G Knepley { 8569298eaa6SMatthew G Knepley PetscSection coordSection; 8579298eaa6SMatthew G Knepley Vec coordinates; 8589298eaa6SMatthew G Knepley PetscScalar *coords; 8599298eaa6SMatthew G Knepley PetscInt coordSize, v, d; 8609298eaa6SMatthew G Knepley PetscErrorCode ierr; 8619298eaa6SMatthew G Knepley 8629298eaa6SMatthew G Knepley PetscFunctionBegin; 863c2166f76SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 8649298eaa6SMatthew G Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 8659298eaa6SMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 8669298eaa6SMatthew G Knepley ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 8679298eaa6SMatthew G Knepley for (v = numCells; v < numCells+numVertices; ++v) { 8689298eaa6SMatthew G Knepley ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 8699298eaa6SMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 8709298eaa6SMatthew G Knepley } 8719298eaa6SMatthew G Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 8729298eaa6SMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 8739298eaa6SMatthew G Knepley ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 8749298eaa6SMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 8759298eaa6SMatthew G Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 8762eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 8779298eaa6SMatthew G Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 8789298eaa6SMatthew G Knepley for (v = 0; v < numVertices; ++v) { 8799298eaa6SMatthew G Knepley for (d = 0; d < spaceDim; ++d) { 8809298eaa6SMatthew G Knepley coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 8819298eaa6SMatthew G Knepley } 8829298eaa6SMatthew G Knepley } 8839298eaa6SMatthew G Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 8849298eaa6SMatthew G Knepley ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 8859298eaa6SMatthew G Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 8869298eaa6SMatthew G Knepley PetscFunctionReturn(0); 8879298eaa6SMatthew G Knepley } 8889298eaa6SMatthew G Knepley 8899298eaa6SMatthew G Knepley #undef __FUNCT__ 8909298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromCellList" 8919298eaa6SMatthew G Knepley /*@C 8929298eaa6SMatthew G Knepley DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 8939298eaa6SMatthew G Knepley 8949298eaa6SMatthew G Knepley Input Parameters: 8959298eaa6SMatthew G Knepley + comm - The communicator 8969298eaa6SMatthew G Knepley . dim - The topological dimension of the mesh 8979298eaa6SMatthew G Knepley . numCells - The number of cells 8989298eaa6SMatthew G Knepley . numVertices - The number of vertices 8999298eaa6SMatthew G Knepley . numCorners - The number of vertices for each cell 9009298eaa6SMatthew G Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 9019298eaa6SMatthew G Knepley . cells - An array of numCells*numCorners numbers, the vertices for each cell 9029298eaa6SMatthew G Knepley . spaceDim - The spatial dimension used for coordinates 9039298eaa6SMatthew G Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 9049298eaa6SMatthew G Knepley 9059298eaa6SMatthew G Knepley Output Parameter: 9069298eaa6SMatthew G Knepley . dm - The DM 9079298eaa6SMatthew G Knepley 9089298eaa6SMatthew G Knepley Note: Two triangles sharing a face 9099298eaa6SMatthew G Knepley $ 9109298eaa6SMatthew G Knepley $ 2 9119298eaa6SMatthew G Knepley $ / | \ 9129298eaa6SMatthew G Knepley $ / | \ 9139298eaa6SMatthew G Knepley $ / | \ 9149298eaa6SMatthew G Knepley $ 0 0 | 1 3 9159298eaa6SMatthew G Knepley $ \ | / 9169298eaa6SMatthew G Knepley $ \ | / 9179298eaa6SMatthew G Knepley $ \ | / 9189298eaa6SMatthew G Knepley $ 1 9199298eaa6SMatthew G Knepley would have input 9209298eaa6SMatthew G Knepley $ numCells = 2, numVertices = 4 9219298eaa6SMatthew G Knepley $ cells = [0 1 2 1 3 2] 9229298eaa6SMatthew G Knepley $ 9239298eaa6SMatthew G Knepley which would result in the DMPlex 9249298eaa6SMatthew G Knepley $ 9259298eaa6SMatthew G Knepley $ 4 9269298eaa6SMatthew G Knepley $ / | \ 9279298eaa6SMatthew G Knepley $ / | \ 9289298eaa6SMatthew G Knepley $ / | \ 9299298eaa6SMatthew G Knepley $ 2 0 | 1 5 9309298eaa6SMatthew G Knepley $ \ | / 9319298eaa6SMatthew G Knepley $ \ | / 9329298eaa6SMatthew G Knepley $ \ | / 9339298eaa6SMatthew G Knepley $ 3 9349298eaa6SMatthew G Knepley 9359298eaa6SMatthew G Knepley Level: beginner 9369298eaa6SMatthew G Knepley 937939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromDAG(), DMPlexCreate() 9389298eaa6SMatthew G Knepley @*/ 9399298eaa6SMatthew 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) 9409298eaa6SMatthew G Knepley { 9419298eaa6SMatthew G Knepley PetscErrorCode ierr; 9429298eaa6SMatthew G Knepley 9439298eaa6SMatthew G Knepley PetscFunctionBegin; 9449298eaa6SMatthew G Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 9459298eaa6SMatthew G Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 9469298eaa6SMatthew G Knepley ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 9479298eaa6SMatthew G Knepley ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 9489298eaa6SMatthew G Knepley if (interpolate) { 9499298eaa6SMatthew G Knepley DM idm; 9509298eaa6SMatthew G Knepley 9519298eaa6SMatthew G Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 9529298eaa6SMatthew G Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 9539298eaa6SMatthew G Knepley *dm = idm; 9549298eaa6SMatthew G Knepley } 9559298eaa6SMatthew G Knepley ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 9569298eaa6SMatthew G Knepley PetscFunctionReturn(0); 9579298eaa6SMatthew G Knepley } 9589298eaa6SMatthew G Knepley 9599298eaa6SMatthew G Knepley #undef __FUNCT__ 9609298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromDAG" 961939f6067SMatthew G. Knepley /*@ 962939f6067SMatthew G. Knepley DMPlexCreateFromDAG - This takes as input the adjacency-list representation of the Directed Acyclic Graph (Hasse Diagram) encoding a mesh, and produces a DM 963939f6067SMatthew G. Knepley 964939f6067SMatthew G. Knepley Input Parameters: 965939f6067SMatthew G. Knepley + dm - The empty DM object, usually from DMCreate() and DMPlexSetDimension() 966939f6067SMatthew G. Knepley . depth - The depth of the DAG 967939f6067SMatthew G. Knepley . numPoints - The number of points at each depth 968939f6067SMatthew G. Knepley . coneSize - The cone size of each point 969939f6067SMatthew G. Knepley . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 970939f6067SMatthew G. Knepley . coneOrientations - The orientation of each cone point 971939f6067SMatthew G. Knepley - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex 972939f6067SMatthew G. Knepley 973939f6067SMatthew G. Knepley Output Parameter: 974939f6067SMatthew G. Knepley . dm - The DM 975939f6067SMatthew G. Knepley 976939f6067SMatthew G. Knepley Note: Two triangles sharing a face would have input 977939f6067SMatthew G. Knepley $ depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 978939f6067SMatthew G. Knepley $ cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 979939f6067SMatthew G. Knepley $ vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 980939f6067SMatthew G. Knepley $ 981939f6067SMatthew G. Knepley which would result in the DMPlex 982939f6067SMatthew G. Knepley $ 983939f6067SMatthew G. Knepley $ 4 984939f6067SMatthew G. Knepley $ / | \ 985939f6067SMatthew G. Knepley $ / | \ 986939f6067SMatthew G. Knepley $ / | \ 987939f6067SMatthew G. Knepley $ 2 0 | 1 5 988939f6067SMatthew G. Knepley $ \ | / 989939f6067SMatthew G. Knepley $ \ | / 990939f6067SMatthew G. Knepley $ \ | / 991939f6067SMatthew G. Knepley $ 3 992939f6067SMatthew G. Knepley $ 993939f6067SMatthew G. Knepley $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList() 994939f6067SMatthew G. Knepley 995939f6067SMatthew G. Knepley Level: advanced 996939f6067SMatthew G. Knepley 997939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromCellList(), DMPlexCreate() 998939f6067SMatthew G. Knepley @*/ 9999298eaa6SMatthew G Knepley PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 10009298eaa6SMatthew G Knepley { 10019298eaa6SMatthew G Knepley Vec coordinates; 10029298eaa6SMatthew G Knepley PetscSection coordSection; 10039298eaa6SMatthew G Knepley PetscScalar *coords; 10049298eaa6SMatthew G Knepley PetscInt coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off; 10059298eaa6SMatthew G Knepley PetscErrorCode ierr; 10069298eaa6SMatthew G Knepley 10079298eaa6SMatthew G Knepley PetscFunctionBegin; 10089298eaa6SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 10099298eaa6SMatthew G Knepley for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 10109298eaa6SMatthew G Knepley ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 10119298eaa6SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 10129298eaa6SMatthew G Knepley ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 10139298eaa6SMatthew G Knepley } 10149298eaa6SMatthew G Knepley ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 10159298eaa6SMatthew G Knepley for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 10169298eaa6SMatthew G Knepley ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 10179298eaa6SMatthew G Knepley ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 10189298eaa6SMatthew G Knepley } 10199298eaa6SMatthew G Knepley ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 10209298eaa6SMatthew G Knepley ierr = DMPlexStratify(dm);CHKERRQ(ierr); 10219298eaa6SMatthew G Knepley /* Build coordinates */ 1022c2166f76SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 10239298eaa6SMatthew G Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 10249298eaa6SMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 10259298eaa6SMatthew G Knepley ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 10269298eaa6SMatthew G Knepley for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 10279298eaa6SMatthew G Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 10286f8cbbeeSMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 10299298eaa6SMatthew G Knepley } 10309298eaa6SMatthew G Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 10319298eaa6SMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 10329298eaa6SMatthew G Knepley ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 10336f8cbbeeSMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 10349298eaa6SMatthew G Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 10352eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 10369298eaa6SMatthew G Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 10379298eaa6SMatthew G Knepley for (v = 0; v < numPoints[0]; ++v) { 10389298eaa6SMatthew G Knepley PetscInt off; 10399298eaa6SMatthew G Knepley 10409298eaa6SMatthew G Knepley ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 10419298eaa6SMatthew G Knepley for (d = 0; d < dim; ++d) { 10429298eaa6SMatthew G Knepley coords[off+d] = vertexCoords[v*dim+d]; 10439298eaa6SMatthew G Knepley } 10449298eaa6SMatthew G Knepley } 10459298eaa6SMatthew G Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 10469298eaa6SMatthew G Knepley ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 10479298eaa6SMatthew G Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 10489298eaa6SMatthew G Knepley PetscFunctionReturn(0); 10499298eaa6SMatthew G Knepley } 1050