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__ 71df5d5c5SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateDoublet" 81df5d5c5SMatthew G. Knepley /*@ 91df5d5c5SMatthew G. Knepley DMPlexCreateDoublet - Creates a mesh of two cells of the specified type, optionally with later refinement. 101df5d5c5SMatthew G. Knepley 111df5d5c5SMatthew G. Knepley Collective on MPI_Comm 121df5d5c5SMatthew G. Knepley 131df5d5c5SMatthew G. Knepley Input Parameters: 141df5d5c5SMatthew G. Knepley + comm - The communicator for the DM object 151df5d5c5SMatthew G. Knepley . dim - The spatial dimension 161df5d5c5SMatthew G. Knepley . simplex - Flag for simplicial cells, otherwise they are tensor product cells 171df5d5c5SMatthew G. Knepley . interpolate - Flag to create intermediate mesh pieces (edges, faces) 181df5d5c5SMatthew G. Knepley . refinementUniform - Flag for uniform parallel refinement 191df5d5c5SMatthew G. Knepley - refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell 201df5d5c5SMatthew G. Knepley 211df5d5c5SMatthew G. Knepley Output Parameter: 221df5d5c5SMatthew G. Knepley . dm - The DM object 231df5d5c5SMatthew G. Knepley 241df5d5c5SMatthew G. Knepley Level: beginner 251df5d5c5SMatthew G. Knepley 261df5d5c5SMatthew G. Knepley .keywords: DM, create 271df5d5c5SMatthew G. Knepley .seealso: DMSetType(), DMCreate() 281df5d5c5SMatthew G. Knepley @*/ 291df5d5c5SMatthew G. Knepley PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscBool refinementUniform, PetscReal refinementLimit, DM *newdm) 301df5d5c5SMatthew G. Knepley { 311df5d5c5SMatthew G. Knepley DM dm; 321df5d5c5SMatthew G. Knepley PetscInt p; 331df5d5c5SMatthew G. Knepley PetscMPIInt rank; 341df5d5c5SMatthew G. Knepley PetscErrorCode ierr; 351df5d5c5SMatthew G. Knepley 361df5d5c5SMatthew G. Knepley PetscFunctionBegin; 371df5d5c5SMatthew G. Knepley ierr = DMCreate(comm, &dm);CHKERRQ(ierr); 381df5d5c5SMatthew G. Knepley ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr); 391df5d5c5SMatthew G. Knepley ierr = DMPlexSetDimension(dm, dim);CHKERRQ(ierr); 401df5d5c5SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 41ce78fa2fSMatthew G. Knepley switch (dim) { 42ce78fa2fSMatthew G. Knepley case 2: 43ce78fa2fSMatthew G. Knepley if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);} 44ce78fa2fSMatthew G. Knepley else {ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);} 45ce78fa2fSMatthew G. Knepley break; 46ce78fa2fSMatthew G. Knepley case 3: 47ce78fa2fSMatthew G. Knepley if (simplex) {ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);} 48ce78fa2fSMatthew G. Knepley else {ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);} 49ce78fa2fSMatthew G. Knepley break; 50ce78fa2fSMatthew G. Knepley default: 51ce78fa2fSMatthew G. Knepley SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 52ce78fa2fSMatthew G. Knepley } 531df5d5c5SMatthew G. Knepley if (rank) { 541df5d5c5SMatthew G. Knepley PetscInt numPoints[2] = {0, 0}; 551df5d5c5SMatthew G. Knepley ierr = DMPlexCreateFromDAG(dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 561df5d5c5SMatthew G. Knepley } else { 571df5d5c5SMatthew G. Knepley switch (dim) { 581df5d5c5SMatthew G. Knepley case 2: 591df5d5c5SMatthew G. Knepley if (simplex) { 601df5d5c5SMatthew G. Knepley PetscInt numPoints[2] = {4, 2}; 611df5d5c5SMatthew G. Knepley PetscInt coneSize[6] = {3, 3, 0, 0, 0, 0}; 621df5d5c5SMatthew G. Knepley PetscInt cones[6] = {2, 3, 4, 5, 4, 3}; 631df5d5c5SMatthew G. Knepley PetscInt coneOrientations[6] = {0, 0, 0, 0, 0, 0}; 641df5d5c5SMatthew G. Knepley PetscScalar vertexCoords[8] = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5}; 651df5d5c5SMatthew G. Knepley PetscInt markerPoints[8] = {2, 1, 3, 1, 4, 1, 5, 1}; 661df5d5c5SMatthew G. Knepley 671df5d5c5SMatthew G. Knepley ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 681df5d5c5SMatthew G. Knepley for (p = 0; p < 4; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);} 691df5d5c5SMatthew G. Knepley } else { 701df5d5c5SMatthew G. Knepley PetscInt numPoints[2] = {6, 2}; 711df5d5c5SMatthew G. Knepley PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 721df5d5c5SMatthew G. Knepley PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4}; 731df5d5c5SMatthew G. Knepley PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 741df5d5c5SMatthew 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}; 751df5d5c5SMatthew G. Knepley 761df5d5c5SMatthew G. Knepley ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 771df5d5c5SMatthew G. Knepley } 781df5d5c5SMatthew G. Knepley break; 791df5d5c5SMatthew G. Knepley case 3: 801df5d5c5SMatthew G. Knepley if (simplex) { 811df5d5c5SMatthew G. Knepley PetscInt numPoints[2] = {5, 2}; 821df5d5c5SMatthew G. Knepley PetscInt coneSize[7] = {4, 4, 0, 0, 0, 0, 0}; 831df5d5c5SMatthew G. Knepley PetscInt cones[8] = {4, 3, 5, 2, 5, 3, 4, 6}; 841df5d5c5SMatthew G. Knepley PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 851df5d5c5SMatthew 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}; 861df5d5c5SMatthew G. Knepley PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 871df5d5c5SMatthew G. Knepley 881df5d5c5SMatthew G. Knepley ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 891df5d5c5SMatthew G. Knepley for (p = 0; p < 5; ++p) {ierr = DMPlexSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);} 901df5d5c5SMatthew G. Knepley } else { 911df5d5c5SMatthew G. Knepley PetscInt numPoints[2] = {12, 2}; 921df5d5c5SMatthew G. Knepley PetscInt coneSize[14] = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 931df5d5c5SMatthew G. Knepley PetscInt cones[16] = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8}; 941df5d5c5SMatthew G. Knepley PetscInt coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; 951df5d5c5SMatthew 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, 961df5d5c5SMatthew 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, 971df5d5c5SMatthew 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}; 981df5d5c5SMatthew G. Knepley 991df5d5c5SMatthew G. Knepley ierr = DMPlexCreateFromDAG(dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 1001df5d5c5SMatthew G. Knepley } 1011df5d5c5SMatthew G. Knepley break; 1021df5d5c5SMatthew G. Knepley default: 1031df5d5c5SMatthew G. Knepley SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 1041df5d5c5SMatthew G. Knepley } 1051df5d5c5SMatthew G. Knepley } 1061df5d5c5SMatthew G. Knepley *newdm = dm; 1071df5d5c5SMatthew G. Knepley if (refinementLimit > 0.0) { 1081df5d5c5SMatthew G. Knepley DM rdm; 1091df5d5c5SMatthew G. Knepley const char *name; 1101df5d5c5SMatthew G. Knepley 1111df5d5c5SMatthew G. Knepley ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr); 1121df5d5c5SMatthew G. Knepley ierr = DMPlexSetRefinementLimit(*newdm, refinementLimit);CHKERRQ(ierr); 1131df5d5c5SMatthew G. Knepley ierr = DMRefine(*newdm, comm, &rdm);CHKERRQ(ierr); 1141df5d5c5SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr); 1151df5d5c5SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) rdm, name);CHKERRQ(ierr); 1161df5d5c5SMatthew G. Knepley ierr = DMDestroy(newdm);CHKERRQ(ierr); 1171df5d5c5SMatthew G. Knepley *newdm = rdm; 1181df5d5c5SMatthew G. Knepley } 1191df5d5c5SMatthew G. Knepley if (interpolate) { 1201df5d5c5SMatthew G. Knepley DM idm; 1211df5d5c5SMatthew G. Knepley const char *name; 1221df5d5c5SMatthew G. Knepley 1231df5d5c5SMatthew G. Knepley ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr); 1241df5d5c5SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr); 1251df5d5c5SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) idm, name);CHKERRQ(ierr); 1261df5d5c5SMatthew G. Knepley ierr = DMPlexCopyCoordinates(*newdm, idm);CHKERRQ(ierr); 1271df5d5c5SMatthew G. Knepley ierr = DMPlexCopyLabels(*newdm, idm);CHKERRQ(ierr); 1281df5d5c5SMatthew G. Knepley ierr = DMDestroy(newdm);CHKERRQ(ierr); 1291df5d5c5SMatthew G. Knepley *newdm = idm; 1301df5d5c5SMatthew G. Knepley } 1311df5d5c5SMatthew G. Knepley { 1321df5d5c5SMatthew G. Knepley DM refinedMesh = NULL; 1331df5d5c5SMatthew G. Knepley DM distributedMesh = NULL; 1341df5d5c5SMatthew G. Knepley 1351df5d5c5SMatthew G. Knepley /* Distribute mesh over processes */ 13607104863SJed Brown ierr = DMPlexDistribute(*newdm, NULL, 0, NULL, &distributedMesh);CHKERRQ(ierr); 1371df5d5c5SMatthew G. Knepley if (distributedMesh) { 1381df5d5c5SMatthew G. Knepley ierr = DMDestroy(newdm);CHKERRQ(ierr); 1391df5d5c5SMatthew G. Knepley *newdm = distributedMesh; 1401df5d5c5SMatthew G. Knepley } 1411df5d5c5SMatthew G. Knepley if (refinementUniform) { 1421df5d5c5SMatthew G. Knepley ierr = DMPlexSetRefinementUniform(*newdm, refinementUniform);CHKERRQ(ierr); 1431df5d5c5SMatthew G. Knepley ierr = DMRefine(*newdm, comm, &refinedMesh);CHKERRQ(ierr); 1441df5d5c5SMatthew G. Knepley if (refinedMesh) { 1451df5d5c5SMatthew G. Knepley ierr = DMDestroy(newdm);CHKERRQ(ierr); 1461df5d5c5SMatthew G. Knepley *newdm = refinedMesh; 1471df5d5c5SMatthew G. Knepley } 1481df5d5c5SMatthew G. Knepley } 1491df5d5c5SMatthew G. Knepley } 1501df5d5c5SMatthew G. Knepley PetscFunctionReturn(0); 1511df5d5c5SMatthew G. Knepley } 1521df5d5c5SMatthew G. Knepley 1531df5d5c5SMatthew G. Knepley #undef __FUNCT__ 154552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareBoundary" 15526492d91SMatthew G. Knepley /*@ 15626492d91SMatthew G. Knepley DMPlexCreateSquareBoundary - Creates a 1D mesh the is the boundary of a square lattice. 157552f7358SJed Brown 15826492d91SMatthew G. Knepley Collective on MPI_Comm 15926492d91SMatthew G. Knepley 16026492d91SMatthew G. Knepley Input Parameters: 16126492d91SMatthew G. Knepley + comm - The communicator for the DM object 16226492d91SMatthew G. Knepley . lower - The lower left corner coordinates 16326492d91SMatthew G. Knepley . upper - The upper right corner coordinates 16426492d91SMatthew G. Knepley - edges - The number of cells in each direction 16526492d91SMatthew G. Knepley 16626492d91SMatthew G. Knepley Output Parameter: 16726492d91SMatthew G. Knepley . dm - The DM object 16826492d91SMatthew G. Knepley 16926492d91SMatthew G. Knepley Note: Here is the numbering returned for 2 cells in each direction: 17026492d91SMatthew G. Knepley $ 18--5-17--4--16 17126492d91SMatthew G. Knepley $ | | | 17226492d91SMatthew G. Knepley $ 6 10 3 17326492d91SMatthew G. Knepley $ | | | 17426492d91SMatthew G. Knepley $ 19-11-20--9--15 17526492d91SMatthew G. Knepley $ | | | 17626492d91SMatthew G. Knepley $ 7 8 2 17726492d91SMatthew G. Knepley $ | | | 17826492d91SMatthew G. Knepley $ 12--0-13--1--14 17926492d91SMatthew G. Knepley 18026492d91SMatthew G. Knepley Level: beginner 18126492d91SMatthew G. Knepley 18226492d91SMatthew G. Knepley .keywords: DM, create 18326492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate() 18426492d91SMatthew G. Knepley @*/ 185552f7358SJed Brown PetscErrorCode DMPlexCreateSquareBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[]) 186552f7358SJed Brown { 187552f7358SJed Brown PetscInt numVertices = (edges[0]+1)*(edges[1]+1); 188552f7358SJed Brown PetscInt numEdges = edges[0]*(edges[1]+1) + (edges[0]+1)*edges[1]; 189552f7358SJed Brown PetscInt markerTop = 1; 190552f7358SJed Brown PetscInt markerBottom = 1; 191552f7358SJed Brown PetscInt markerRight = 1; 192552f7358SJed Brown PetscInt markerLeft = 1; 193552f7358SJed Brown PetscBool markerSeparate = PETSC_FALSE; 194552f7358SJed Brown Vec coordinates; 195552f7358SJed Brown PetscSection coordSection; 196552f7358SJed Brown PetscScalar *coords; 197552f7358SJed Brown PetscInt coordSize; 198552f7358SJed Brown PetscMPIInt rank; 199552f7358SJed Brown PetscInt v, vx, vy; 200552f7358SJed Brown PetscErrorCode ierr; 201552f7358SJed Brown 202552f7358SJed Brown PetscFunctionBegin; 2030298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 204552f7358SJed Brown if (markerSeparate) { 205552f7358SJed Brown markerTop = 1; 206552f7358SJed Brown markerBottom = 0; 207552f7358SJed Brown markerRight = 0; 208552f7358SJed Brown markerLeft = 0; 209552f7358SJed Brown } 21082f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 211552f7358SJed Brown if (!rank) { 212552f7358SJed Brown PetscInt e, ex, ey; 213552f7358SJed Brown 214552f7358SJed Brown ierr = DMPlexSetChart(dm, 0, numEdges+numVertices);CHKERRQ(ierr); 215552f7358SJed Brown for (e = 0; e < numEdges; ++e) { 216552f7358SJed Brown ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 217552f7358SJed Brown } 218552f7358SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 219552f7358SJed Brown for (vx = 0; vx <= edges[0]; vx++) { 220552f7358SJed Brown for (ey = 0; ey < edges[1]; ey++) { 221552f7358SJed Brown PetscInt edge = vx*edges[1] + ey + edges[0]*(edges[1]+1); 222552f7358SJed Brown PetscInt vertex = ey*(edges[0]+1) + vx + numEdges; 223da80777bSKarl Rupp PetscInt cone[2]; 224552f7358SJed Brown 225da80777bSKarl Rupp cone[0] = vertex; cone[1] = vertex+edges[0]+1; 226552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 227552f7358SJed Brown if (vx == edges[0]) { 228552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 229552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 230552f7358SJed Brown if (ey == edges[1]-1) { 231552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 232552f7358SJed Brown } 233552f7358SJed Brown } else if (vx == 0) { 234552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 235552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 236552f7358SJed Brown if (ey == edges[1]-1) { 237552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 238552f7358SJed Brown } 239552f7358SJed Brown } 240552f7358SJed Brown } 241552f7358SJed Brown } 242552f7358SJed Brown for (vy = 0; vy <= edges[1]; vy++) { 243552f7358SJed Brown for (ex = 0; ex < edges[0]; ex++) { 244552f7358SJed Brown PetscInt edge = vy*edges[0] + ex; 245552f7358SJed Brown PetscInt vertex = vy*(edges[0]+1) + ex + numEdges; 246da80777bSKarl Rupp PetscInt cone[2]; 247552f7358SJed Brown 248da80777bSKarl Rupp cone[0] = vertex; cone[1] = vertex+1; 249552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 250552f7358SJed Brown if (vy == edges[1]) { 251552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 252552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 253552f7358SJed Brown if (ex == edges[0]-1) { 254552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 255552f7358SJed Brown } 256552f7358SJed Brown } else if (vy == 0) { 257552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 258552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 259552f7358SJed Brown if (ex == edges[0]-1) { 260552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 261552f7358SJed Brown } 262552f7358SJed Brown } 263552f7358SJed Brown } 264552f7358SJed Brown } 265552f7358SJed Brown } 266552f7358SJed Brown ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 267552f7358SJed Brown ierr = DMPlexStratify(dm);CHKERRQ(ierr); 268552f7358SJed Brown /* Build coordinates */ 269c2166f76SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 270552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numEdges, numEdges + numVertices);CHKERRQ(ierr); 271552f7358SJed Brown for (v = numEdges; v < numEdges+numVertices; ++v) { 272552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 273552f7358SJed Brown } 274552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 275552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 27682f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 277552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 2782eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 279552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 280552f7358SJed Brown for (vy = 0; vy <= edges[1]; ++vy) { 281552f7358SJed Brown for (vx = 0; vx <= edges[0]; ++vx) { 282552f7358SJed Brown coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/edges[0])*vx; 283552f7358SJed Brown coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/edges[1])*vy; 284552f7358SJed Brown } 285552f7358SJed Brown } 286552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 287552f7358SJed Brown ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 288552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 289552f7358SJed Brown PetscFunctionReturn(0); 290552f7358SJed Brown } 291552f7358SJed Brown 292552f7358SJed Brown #undef __FUNCT__ 293552f7358SJed Brown #define __FUNCT__ "DMPlexCreateCubeBoundary" 29426492d91SMatthew G. Knepley /*@ 29526492d91SMatthew G. Knepley DMPlexCreateCubeBoundary - Creates a 2D mesh the is the boundary of a cubic lattice. 296552f7358SJed Brown 29726492d91SMatthew G. Knepley Collective on MPI_Comm 29826492d91SMatthew G. Knepley 29926492d91SMatthew G. Knepley Input Parameters: 30026492d91SMatthew G. Knepley + comm - The communicator for the DM object 30126492d91SMatthew G. Knepley . lower - The lower left front corner coordinates 30226492d91SMatthew G. Knepley . upper - The upper right back corner coordinates 30326492d91SMatthew G. Knepley - edges - The number of cells in each direction 30426492d91SMatthew G. Knepley 30526492d91SMatthew G. Knepley Output Parameter: 30626492d91SMatthew G. Knepley . dm - The DM object 30726492d91SMatthew G. Knepley 30826492d91SMatthew G. Knepley Level: beginner 30926492d91SMatthew G. Knepley 31026492d91SMatthew G. Knepley .keywords: DM, create 31126492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMSetType(), DMCreate() 31226492d91SMatthew G. Knepley @*/ 313552f7358SJed Brown PetscErrorCode DMPlexCreateCubeBoundary(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt faces[]) 314552f7358SJed Brown { 315*9e8abbc3SMichael Lange PetscInt vertices[3], numVertices; 3167b59f5a9SMichael Lange PetscInt numFaces = 2*faces[0]*faces[1] + 2*faces[1]*faces[2] + 2*faces[0]*faces[2]; 317552f7358SJed Brown Vec coordinates; 318552f7358SJed Brown PetscSection coordSection; 319552f7358SJed Brown PetscScalar *coords; 320552f7358SJed Brown PetscInt coordSize; 321552f7358SJed Brown PetscMPIInt rank; 322552f7358SJed Brown PetscInt v, vx, vy, vz; 3237b59f5a9SMichael Lange PetscInt voffset, iface=0, cone[4]; 324552f7358SJed Brown PetscErrorCode ierr; 325552f7358SJed Brown 326552f7358SJed Brown PetscFunctionBegin; 32782f516ccSBarry 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"); 32882f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 329*9e8abbc3SMichael Lange vertices[0] = faces[0]+1; vertices[1] = faces[1]+1; vertices[2] = faces[2]+1; 330*9e8abbc3SMichael Lange numVertices = vertices[0]*vertices[1]*vertices[2]; 331552f7358SJed Brown if (!rank) { 332552f7358SJed Brown PetscInt f; 333552f7358SJed Brown 334552f7358SJed Brown ierr = DMPlexSetChart(dm, 0, numFaces+numVertices);CHKERRQ(ierr); 335552f7358SJed Brown for (f = 0; f < numFaces; ++f) { 336552f7358SJed Brown ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 337552f7358SJed Brown } 338552f7358SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 339552f7358SJed Brown for (v = 0; v < numFaces+numVertices; ++v) { 340552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", v, 1);CHKERRQ(ierr); 341552f7358SJed Brown } 3427b59f5a9SMichael Lange 3437b59f5a9SMichael Lange /* Side 0 (Top) */ 3447b59f5a9SMichael Lange for (vy = 0; vy < faces[1]; vy++) { 3457b59f5a9SMichael Lange for (vx = 0; vx < faces[0]; vx++) { 3467b59f5a9SMichael Lange voffset = numFaces + vertices[0]*vertices[1]*(vertices[2]-1) + vy*vertices[0] + vx; 3477b59f5a9SMichael Lange cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]+1; cone[3] = voffset+vertices[0]; 3487b59f5a9SMichael Lange ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 3497b59f5a9SMichael Lange iface++; 350552f7358SJed Brown } 351552f7358SJed Brown } 3527b59f5a9SMichael Lange 3537b59f5a9SMichael Lange /* Side 1 (Bottom) */ 3547b59f5a9SMichael Lange for (vy = 0; vy < faces[1]; vy++) { 3557b59f5a9SMichael Lange for (vx = 0; vx < faces[0]; vx++) { 3567b59f5a9SMichael Lange voffset = numFaces + vy*(faces[0]+1) + vx; 3577b59f5a9SMichael Lange cone[0] = voffset+1; cone[1] = voffset; cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]+1; 3587b59f5a9SMichael Lange ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 3597b59f5a9SMichael Lange iface++; 360552f7358SJed Brown } 361552f7358SJed Brown } 3627b59f5a9SMichael Lange 3637b59f5a9SMichael Lange /* Side 2 (Front) */ 3647b59f5a9SMichael Lange for (vz = 0; vz < faces[2]; vz++) { 3657b59f5a9SMichael Lange for (vx = 0; vx < faces[0]; vx++) { 3667b59f5a9SMichael Lange voffset = numFaces + vz*vertices[0]*vertices[1] + vx; 3677b59f5a9SMichael Lange cone[0] = voffset; cone[1] = voffset+1; cone[2] = voffset+vertices[0]*vertices[1]+1; cone[3] = voffset+vertices[0]*vertices[1]; 3687b59f5a9SMichael Lange ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 3697b59f5a9SMichael Lange iface++; 370552f7358SJed Brown } 3717b59f5a9SMichael Lange } 3727b59f5a9SMichael Lange 3737b59f5a9SMichael Lange /* Side 3 (Back) */ 3747b59f5a9SMichael Lange for (vz = 0; vz < faces[2]; vz++) { 3757b59f5a9SMichael Lange for (vx = 0; vx < faces[0]; vx++) { 3767b59f5a9SMichael Lange voffset = numFaces + vz*vertices[0]*vertices[1] + vertices[0]*(vertices[1]-1) + vx; 3777b59f5a9SMichael Lange cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset+vertices[0]*vertices[1]+1; 3787b59f5a9SMichael Lange cone[2] = voffset+1; cone[3] = voffset; 3797b59f5a9SMichael Lange ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 3807b59f5a9SMichael Lange iface++; 3817b59f5a9SMichael Lange } 3827b59f5a9SMichael Lange } 3837b59f5a9SMichael Lange 3847b59f5a9SMichael Lange /* Side 4 (Left) */ 3857b59f5a9SMichael Lange for (vz = 0; vz < faces[2]; vz++) { 3867b59f5a9SMichael Lange for (vy = 0; vy < faces[1]; vy++) { 3877b59f5a9SMichael Lange voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0]; 3887b59f5a9SMichael Lange cone[0] = voffset; cone[1] = voffset+vertices[0]*vertices[1]; 3897b59f5a9SMichael Lange cone[2] = voffset+vertices[0]*vertices[1]+vertices[0]; cone[3] = voffset+vertices[0]; 3907b59f5a9SMichael Lange ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 3917b59f5a9SMichael Lange iface++; 3927b59f5a9SMichael Lange } 3937b59f5a9SMichael Lange } 3947b59f5a9SMichael Lange 3957b59f5a9SMichael Lange /* Side 5 (Right) */ 3967b59f5a9SMichael Lange for (vz = 0; vz < faces[2]; vz++) { 3977b59f5a9SMichael Lange for (vy = 0; vy < faces[1]; vy++) { 3987b59f5a9SMichael Lange voffset = numFaces + vz*vertices[0]*vertices[1] + vy*vertices[0] + vx; 3997b59f5a9SMichael Lange cone[0] = voffset+vertices[0]*vertices[1]; cone[1] = voffset; 4007b59f5a9SMichael Lange cone[2] = voffset+vertices[0]; cone[3] = voffset+vertices[0]*vertices[1]+vertices[0]; 4017b59f5a9SMichael Lange ierr = DMPlexSetCone(dm, iface, cone);CHKERRQ(ierr); 4027b59f5a9SMichael Lange iface++; 4037b59f5a9SMichael Lange } 404552f7358SJed Brown } 405552f7358SJed Brown } 406552f7358SJed Brown ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 407552f7358SJed Brown ierr = DMPlexStratify(dm);CHKERRQ(ierr); 408552f7358SJed Brown /* Build coordinates */ 409c2166f76SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 410552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numFaces, numFaces + numVertices);CHKERRQ(ierr); 411552f7358SJed Brown for (v = numFaces; v < numFaces+numVertices; ++v) { 412552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, 3);CHKERRQ(ierr); 413552f7358SJed Brown } 414552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 415552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 41682f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 417552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 418552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 4192eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 420552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 421552f7358SJed Brown for (vz = 0; vz <= faces[2]; ++vz) { 422552f7358SJed Brown for (vy = 0; vy <= faces[1]; ++vy) { 423552f7358SJed Brown for (vx = 0; vx <= faces[0]; ++vx) { 424552f7358SJed Brown coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+0] = lower[0] + ((upper[0] - lower[0])/faces[0])*vx; 425552f7358SJed Brown coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+1] = lower[1] + ((upper[1] - lower[1])/faces[1])*vy; 426552f7358SJed Brown coords[((vz*(faces[1]+1)+vy)*(faces[0]+1)+vx)*3+2] = lower[2] + ((upper[2] - lower[2])/faces[2])*vz; 427552f7358SJed Brown } 428552f7358SJed Brown } 429552f7358SJed Brown } 430552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 431552f7358SJed Brown ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 432552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 433552f7358SJed Brown PetscFunctionReturn(0); 434552f7358SJed Brown } 435552f7358SJed Brown 436552f7358SJed Brown #undef __FUNCT__ 437552f7358SJed Brown #define __FUNCT__ "DMPlexCreateSquareMesh" 43826492d91SMatthew G. Knepley /*@ 43926492d91SMatthew G. Knepley DMPlexCreateSquareMesh - Creates a 2D mesh for a square lattice. 440552f7358SJed Brown 44126492d91SMatthew G. Knepley Collective on MPI_Comm 44226492d91SMatthew G. Knepley 44326492d91SMatthew G. Knepley Input Parameters: 44426492d91SMatthew G. Knepley + comm - The communicator for the DM object 44526492d91SMatthew G. Knepley . lower - The lower left corner coordinates 44626492d91SMatthew G. Knepley . upper - The upper right corner coordinates 447fbf5b483SMatthew G. Knepley . edges - The number of cells in each direction 448fbf5b483SMatthew G. Knepley . bdX - The boundary type for the X direction 449fbf5b483SMatthew G. Knepley - bdY - The boundary type for the Y direction 45026492d91SMatthew G. Knepley 45126492d91SMatthew G. Knepley Output Parameter: 45226492d91SMatthew G. Knepley . dm - The DM object 45326492d91SMatthew G. Knepley 45426492d91SMatthew G. Knepley Note: Here is the numbering returned for 2 cells in each direction: 45526492d91SMatthew G. Knepley $ 22--8-23--9--24 45626492d91SMatthew G. Knepley $ | | | 45726492d91SMatthew G. Knepley $ 13 2 14 3 15 45826492d91SMatthew G. Knepley $ | | | 45926492d91SMatthew G. Knepley $ 19--6-20--7--21 46026492d91SMatthew G. Knepley $ | | | 46126492d91SMatthew G. Knepley $ 10 0 11 1 12 46226492d91SMatthew G. Knepley $ | | | 46326492d91SMatthew G. Knepley $ 16--4-17--5--18 46426492d91SMatthew G. Knepley 46526492d91SMatthew G. Knepley Level: beginner 46626492d91SMatthew G. Knepley 46726492d91SMatthew G. Knepley .keywords: DM, create 46826492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMPlexCreateSquareBoundary(), DMPlexCreateCubeBoundary(), DMSetType(), DMCreate() 46926492d91SMatthew G. Knepley @*/ 470bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateSquareMesh(DM dm, const PetscReal lower[], const PetscReal upper[], const PetscInt edges[], DMBoundaryType bdX, DMBoundaryType bdY) 471552f7358SJed Brown { 472552f7358SJed Brown PetscInt markerTop = 1; 473552f7358SJed Brown PetscInt markerBottom = 1; 474552f7358SJed Brown PetscInt markerRight = 1; 475552f7358SJed Brown PetscInt markerLeft = 1; 476552f7358SJed Brown PetscBool markerSeparate = PETSC_FALSE; 477552f7358SJed Brown PetscMPIInt rank; 478552f7358SJed Brown PetscErrorCode ierr; 479552f7358SJed Brown 480552f7358SJed Brown PetscFunctionBegin; 48182f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 4820298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_plex_separate_marker", &markerSeparate, NULL);CHKERRQ(ierr); 483552f7358SJed Brown if (markerSeparate) { 484552f7358SJed Brown markerTop = 3; 485552f7358SJed Brown markerBottom = 1; 486552f7358SJed Brown markerRight = 2; 487552f7358SJed Brown markerLeft = 4; 488552f7358SJed Brown } 489552f7358SJed Brown { 490552f7358SJed Brown const PetscInt numXEdges = !rank ? edges[0] : 0; 491552f7358SJed Brown const PetscInt numYEdges = !rank ? edges[1] : 0; 49246675b8eSMatthew G. Knepley const PetscInt numXVertices = !rank ? (bdX == DM_BOUNDARY_PERIODIC || bdX == DM_BOUNDARY_TWIST ? edges[0] : edges[0]+1) : 0; 49346675b8eSMatthew G. Knepley const PetscInt numYVertices = !rank ? (bdY == DM_BOUNDARY_PERIODIC || bdY == DM_BOUNDARY_TWIST ? edges[1] : edges[1]+1) : 0; 494552f7358SJed Brown const PetscInt numTotXEdges = numXEdges*numYVertices; 495552f7358SJed Brown const PetscInt numTotYEdges = numYEdges*numXVertices; 496552f7358SJed Brown const PetscInt numVertices = numXVertices*numYVertices; 497552f7358SJed Brown const PetscInt numEdges = numTotXEdges + numTotYEdges; 498552f7358SJed Brown const PetscInt numFaces = numXEdges*numYEdges; 499552f7358SJed Brown const PetscInt firstVertex = numFaces; 500552f7358SJed Brown const PetscInt firstXEdge = numFaces + numVertices; 501552f7358SJed Brown const PetscInt firstYEdge = numFaces + numVertices + numTotXEdges; 502552f7358SJed Brown Vec coordinates; 503552f7358SJed Brown PetscSection coordSection; 504552f7358SJed Brown PetscScalar *coords; 505552f7358SJed Brown PetscInt coordSize; 506552f7358SJed Brown PetscInt v, vx, vy; 507552f7358SJed Brown PetscInt f, fx, fy, e, ex, ey; 508552f7358SJed Brown 509552f7358SJed Brown ierr = DMPlexSetChart(dm, 0, numFaces+numEdges+numVertices);CHKERRQ(ierr); 510552f7358SJed Brown for (f = 0; f < numFaces; ++f) { 511552f7358SJed Brown ierr = DMPlexSetConeSize(dm, f, 4);CHKERRQ(ierr); 512552f7358SJed Brown } 513552f7358SJed Brown for (e = firstXEdge; e < firstXEdge+numEdges; ++e) { 514552f7358SJed Brown ierr = DMPlexSetConeSize(dm, e, 2);CHKERRQ(ierr); 515552f7358SJed Brown } 516552f7358SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 517552f7358SJed Brown /* Build faces */ 518552f7358SJed Brown for (fy = 0; fy < numYEdges; fy++) { 519552f7358SJed Brown for (fx = 0; fx < numXEdges; fx++) { 520552f7358SJed Brown const PetscInt face = fy*numXEdges + fx; 521552f7358SJed Brown const PetscInt edgeL = firstYEdge + fx *numYEdges + fy; 522becd5721SMatthew G. Knepley const PetscInt edgeR = firstYEdge + ((fx+1)%numXVertices)*numYEdges + fy; 523552f7358SJed Brown const PetscInt edgeB = firstXEdge + fy *numXEdges + fx; 524becd5721SMatthew G. Knepley const PetscInt edgeT = firstXEdge + ((fy+1)%numYVertices)*numXEdges + fx; 525552f7358SJed Brown const PetscInt ornt[4] = {0, 0, -2, -2}; 526da80777bSKarl Rupp PetscInt cone[4]; 527552f7358SJed Brown 528becd5721SMatthew G. Knepley cone[0] = edgeB; cone[1] = edgeR; cone[2] = edgeT; cone[3] = edgeL; 529552f7358SJed Brown ierr = DMPlexSetCone(dm, face, cone);CHKERRQ(ierr); 530552f7358SJed Brown ierr = DMPlexSetConeOrientation(dm, face, ornt);CHKERRQ(ierr); 531552f7358SJed Brown } 532552f7358SJed Brown } 533552f7358SJed Brown /* Build Y edges*/ 534552f7358SJed Brown for (vx = 0; vx < numXVertices; vx++) { 535552f7358SJed Brown for (ey = 0; ey < numYEdges; ey++) { 53646675b8eSMatthew G. Knepley const PetscInt nextv = (bdY == DM_BOUNDARY_TWIST && ey == numYEdges-1) ? (numXVertices-vx-1) : ((ey+1)%numYVertices)*numXVertices + vx; 537552f7358SJed Brown const PetscInt edge = firstYEdge + vx*numYEdges + ey; 538becd5721SMatthew G. Knepley const PetscInt vertexB = firstVertex + ey*numXVertices + vx; 53946675b8eSMatthew G. Knepley const PetscInt vertexT = firstVertex + nextv; 540da80777bSKarl Rupp PetscInt cone[2]; 541552f7358SJed Brown 542becd5721SMatthew G. Knepley cone[0] = vertexB; cone[1] = vertexT; 543552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 544c4b39bc1SMatthew G. Knepley if ((bdX != DM_BOUNDARY_PERIODIC) && (bdX != DM_BOUNDARY_TWIST)) { 545552f7358SJed Brown if (vx == numXVertices-1) { 546552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerRight);CHKERRQ(ierr); 547552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerRight);CHKERRQ(ierr); 548552f7358SJed Brown if (ey == numYEdges-1) { 549552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerRight);CHKERRQ(ierr); 550552f7358SJed Brown } 551552f7358SJed Brown } else if (vx == 0) { 552552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerLeft);CHKERRQ(ierr); 553552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerLeft);CHKERRQ(ierr); 554552f7358SJed Brown if (ey == numYEdges-1) { 555552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerLeft);CHKERRQ(ierr); 556552f7358SJed Brown } 557552f7358SJed Brown } 558552f7358SJed Brown } 559552f7358SJed Brown } 560c4b39bc1SMatthew G. Knepley } 561552f7358SJed Brown /* Build X edges*/ 562552f7358SJed Brown for (vy = 0; vy < numYVertices; vy++) { 563552f7358SJed Brown for (ex = 0; ex < numXEdges; ex++) { 56446675b8eSMatthew G. Knepley const PetscInt nextv = (bdX == DM_BOUNDARY_TWIST && ex == numXEdges-1) ? (numYVertices-vy-1)*numXVertices : vy*numXVertices + (ex+1)%numXVertices; 565552f7358SJed Brown const PetscInt edge = firstXEdge + vy*numXEdges + ex; 566becd5721SMatthew G. Knepley const PetscInt vertexL = firstVertex + vy*numXVertices + ex; 56746675b8eSMatthew G. Knepley const PetscInt vertexR = firstVertex + nextv; 568da80777bSKarl Rupp PetscInt cone[2]; 569552f7358SJed Brown 570becd5721SMatthew G. Knepley cone[0] = vertexL; cone[1] = vertexR; 571552f7358SJed Brown ierr = DMPlexSetCone(dm, edge, cone);CHKERRQ(ierr); 572c4b39bc1SMatthew G. Knepley if ((bdY != DM_BOUNDARY_PERIODIC) && (bdY != DM_BOUNDARY_TWIST)) { 573552f7358SJed Brown if (vy == numYVertices-1) { 574552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerTop);CHKERRQ(ierr); 575552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerTop);CHKERRQ(ierr); 576552f7358SJed Brown if (ex == numXEdges-1) { 577552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerTop);CHKERRQ(ierr); 578552f7358SJed Brown } 579552f7358SJed Brown } else if (vy == 0) { 580552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", edge, markerBottom);CHKERRQ(ierr); 581552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[0], markerBottom);CHKERRQ(ierr); 582552f7358SJed Brown if (ex == numXEdges-1) { 583552f7358SJed Brown ierr = DMPlexSetLabelValue(dm, "marker", cone[1], markerBottom);CHKERRQ(ierr); 584552f7358SJed Brown } 585552f7358SJed Brown } 586552f7358SJed Brown } 587552f7358SJed Brown } 588c4b39bc1SMatthew G. Knepley } 589552f7358SJed Brown ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 590552f7358SJed Brown ierr = DMPlexStratify(dm);CHKERRQ(ierr); 591552f7358SJed Brown /* Build coordinates */ 592c2166f76SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 593552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numVertices);CHKERRQ(ierr); 594552f7358SJed Brown for (v = firstVertex; v < firstVertex+numVertices; ++v) { 595552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr); 596552f7358SJed Brown } 597552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 598552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 59982f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 600552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 6012eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 602552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 603552f7358SJed Brown for (vy = 0; vy < numYVertices; ++vy) { 604552f7358SJed Brown for (vx = 0; vx < numXVertices; ++vx) { 605552f7358SJed Brown coords[(vy*numXVertices+vx)*2+0] = lower[0] + ((upper[0] - lower[0])/numXEdges)*vx; 606552f7358SJed Brown coords[(vy*numXVertices+vx)*2+1] = lower[1] + ((upper[1] - lower[1])/numYEdges)*vy; 607552f7358SJed Brown } 608552f7358SJed Brown } 609552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 610552f7358SJed Brown ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 611552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 612552f7358SJed Brown } 613552f7358SJed Brown PetscFunctionReturn(0); 614552f7358SJed Brown } 615552f7358SJed Brown 616552f7358SJed Brown #undef __FUNCT__ 617552f7358SJed Brown #define __FUNCT__ "DMPlexCreateBoxMesh" 6186b75479cSMatthew G Knepley /*@ 61926492d91SMatthew G. Knepley DMPlexCreateBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using simplices. 6206b75479cSMatthew G Knepley 6216b75479cSMatthew G Knepley Collective on MPI_Comm 6226b75479cSMatthew G Knepley 6236b75479cSMatthew G Knepley Input Parameters: 6246b75479cSMatthew G Knepley + comm - The communicator for the DM object 6256b75479cSMatthew G Knepley . dim - The spatial dimension 6266b75479cSMatthew G Knepley - interpolate - Flag to create intermediate mesh pieces (edges, faces) 6276b75479cSMatthew G Knepley 6286b75479cSMatthew G Knepley Output Parameter: 6296b75479cSMatthew G Knepley . dm - The DM object 6306b75479cSMatthew G Knepley 6316b75479cSMatthew G Knepley Level: beginner 6326b75479cSMatthew G Knepley 6336b75479cSMatthew G Knepley .keywords: DM, create 63426492d91SMatthew G. Knepley .seealso: DMPlexCreateHexBoxMesh(), DMSetType(), DMCreate() 6356b75479cSMatthew G Knepley @*/ 6360adebc6cSBarry Smith PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool interpolate, DM *dm) 6370adebc6cSBarry Smith { 638552f7358SJed Brown DM boundary; 639552f7358SJed Brown PetscErrorCode ierr; 640552f7358SJed Brown 641552f7358SJed Brown PetscFunctionBegin; 642552f7358SJed Brown PetscValidPointer(dm, 4); 643552f7358SJed Brown ierr = DMCreate(comm, &boundary);CHKERRQ(ierr); 644552f7358SJed Brown PetscValidLogicalCollectiveInt(boundary,dim,2); 645552f7358SJed Brown ierr = DMSetType(boundary, DMPLEX);CHKERRQ(ierr); 646552f7358SJed Brown ierr = DMPlexSetDimension(boundary, dim-1);CHKERRQ(ierr); 647552f7358SJed Brown switch (dim) { 648552f7358SJed Brown case 2: 649552f7358SJed Brown { 650552f7358SJed Brown PetscReal lower[2] = {0.0, 0.0}; 651552f7358SJed Brown PetscReal upper[2] = {1.0, 1.0}; 652552f7358SJed Brown PetscInt edges[2] = {2, 2}; 653552f7358SJed Brown 654552f7358SJed Brown ierr = DMPlexCreateSquareBoundary(boundary, lower, upper, edges);CHKERRQ(ierr); 655552f7358SJed Brown break; 656552f7358SJed Brown } 657552f7358SJed Brown case 3: 658552f7358SJed Brown { 659552f7358SJed Brown PetscReal lower[3] = {0.0, 0.0, 0.0}; 660552f7358SJed Brown PetscReal upper[3] = {1.0, 1.0, 1.0}; 661552f7358SJed Brown PetscInt faces[3] = {1, 1, 1}; 662552f7358SJed Brown 663552f7358SJed Brown ierr = DMPlexCreateCubeBoundary(boundary, lower, upper, faces);CHKERRQ(ierr); 664552f7358SJed Brown break; 665552f7358SJed Brown } 666552f7358SJed Brown default: 667552f7358SJed Brown SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 668552f7358SJed Brown } 6690298fd71SBarry Smith ierr = DMPlexGenerate(boundary, NULL, interpolate, dm);CHKERRQ(ierr); 670552f7358SJed Brown ierr = DMDestroy(&boundary);CHKERRQ(ierr); 671552f7358SJed Brown PetscFunctionReturn(0); 672552f7358SJed Brown } 673552f7358SJed Brown 674552f7358SJed Brown #undef __FUNCT__ 675552f7358SJed Brown #define __FUNCT__ "DMPlexCreateHexBoxMesh" 67626492d91SMatthew G. Knepley /*@ 67726492d91SMatthew G. Knepley DMPlexCreateHexBoxMesh - Creates a mesh on the tensor product of unit intervals (box) using hexahedra. 67826492d91SMatthew G. Knepley 67926492d91SMatthew G. Knepley Collective on MPI_Comm 68026492d91SMatthew G. Knepley 68126492d91SMatthew G. Knepley Input Parameters: 68226492d91SMatthew G. Knepley + comm - The communicator for the DM object 68326492d91SMatthew G. Knepley . dim - The spatial dimension 684fbf5b483SMatthew G. Knepley . periodicX - The boundary type for the X direction 685fbf5b483SMatthew G. Knepley . periodicY - The boundary type for the Y direction 686fbf5b483SMatthew G. Knepley . periodicZ - The boundary type for the Z direction 68726492d91SMatthew G. Knepley - cells - The number of cells in each direction 68826492d91SMatthew G. Knepley 68926492d91SMatthew G. Knepley Output Parameter: 69026492d91SMatthew G. Knepley . dm - The DM object 69126492d91SMatthew G. Knepley 69226492d91SMatthew G. Knepley Level: beginner 69326492d91SMatthew G. Knepley 69426492d91SMatthew G. Knepley .keywords: DM, create 69526492d91SMatthew G. Knepley .seealso: DMPlexCreateBoxMesh(), DMSetType(), DMCreate() 69626492d91SMatthew G. Knepley @*/ 697bff4a2f0SMatthew G. Knepley PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm comm, PetscInt dim, const PetscInt cells[], DMBoundaryType periodicX, DMBoundaryType periodicY, DMBoundaryType periodicZ, DM *dm) 698a6dfd86eSKarl Rupp { 699552f7358SJed Brown PetscErrorCode ierr; 700552f7358SJed Brown 701552f7358SJed Brown PetscFunctionBegin; 702552f7358SJed Brown PetscValidPointer(dm, 4); 703552f7358SJed Brown ierr = DMCreate(comm, dm);CHKERRQ(ierr); 704552f7358SJed Brown PetscValidLogicalCollectiveInt(*dm,dim,2); 705552f7358SJed Brown ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 706552f7358SJed Brown ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 707552f7358SJed Brown switch (dim) { 708552f7358SJed Brown case 2: 709552f7358SJed Brown { 710552f7358SJed Brown PetscReal lower[2] = {0.0, 0.0}; 711552f7358SJed Brown PetscReal upper[2] = {1.0, 1.0}; 712552f7358SJed Brown 713becd5721SMatthew G. Knepley ierr = DMPlexCreateSquareMesh(*dm, lower, upper, cells, periodicX, periodicY);CHKERRQ(ierr); 714552f7358SJed Brown break; 715552f7358SJed Brown } 716552f7358SJed Brown #if 0 717552f7358SJed Brown case 3: 718552f7358SJed Brown { 719552f7358SJed Brown PetscReal lower[3] = {0.0, 0.0, 0.0}; 720552f7358SJed Brown PetscReal upper[3] = {1.0, 1.0, 1.0}; 721552f7358SJed Brown 722becd5721SMatthew G. Knepley ierr = DMPlexCreateCubeMesh(boundary, lower, upper, cells, periodicX, periodicY, periodicZ);CHKERRQ(ierr); 723552f7358SJed Brown break; 724552f7358SJed Brown } 725552f7358SJed Brown #endif 726552f7358SJed Brown default: 727552f7358SJed Brown SETERRQ1(comm, PETSC_ERR_SUP, "Dimension not supported: %d", dim); 728552f7358SJed Brown } 729552f7358SJed Brown PetscFunctionReturn(0); 730552f7358SJed Brown } 731552f7358SJed Brown 732552f7358SJed Brown /* External function declarations here */ 733552f7358SJed Brown extern PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling); 734bceba477SMatthew G. Knepley extern PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, VecScatter *ctx); 735fd59a867SMatthew G. Knepley extern PetscErrorCode DMCreateDefaultSection_Plex(DM dm); 736b412c318SBarry Smith extern PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J); 737552f7358SJed Brown extern PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm); 738552f7358SJed Brown extern PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *dmRefined); 7390a6ba040SMatthew G. Knepley extern PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened); 7400a6ba040SMatthew G. Knepley extern PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM dmRefined[]); 7419a922efaSMatthew G. Knepley extern PetscErrorCode DMClone_Plex(DM dm, DM *newdm); 742552f7358SJed Brown extern PetscErrorCode DMSetUp_Plex(DM dm); 743552f7358SJed Brown extern PetscErrorCode DMDestroy_Plex(DM dm); 744552f7358SJed Brown extern PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer); 745552f7358SJed Brown extern PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm); 746552f7358SJed Brown extern PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS); 747552f7358SJed Brown 748552f7358SJed Brown #undef __FUNCT__ 7490a6ba040SMatthew G. Knepley #define __FUNCT__ "DMPlexReplace_Static" 7500a6ba040SMatthew G. Knepley /* Replace dm with the contents of dmNew 7510a6ba040SMatthew G. Knepley - Share the DM_Plex structure 7520a6ba040SMatthew G. Knepley - Share the coordinates 7530a6ba040SMatthew G. Knepley */ 7540a6ba040SMatthew G. Knepley static PetscErrorCode DMPlexReplace_Static(DM dm, DM dmNew) 7550a6ba040SMatthew G. Knepley { 7560a6ba040SMatthew G. Knepley PetscSection coordSection; 7570a6ba040SMatthew G. Knepley Vec coords; 7580a6ba040SMatthew G. Knepley PetscErrorCode ierr; 7590a6ba040SMatthew G. Knepley 7600a6ba040SMatthew G. Knepley PetscFunctionBegin; 7610a6ba040SMatthew G. Knepley ierr = DMGetCoordinateSection(dmNew, &coordSection);CHKERRQ(ierr); 7620a6ba040SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dmNew, &coords);CHKERRQ(ierr); 7630a6ba040SMatthew G. Knepley ierr = DMSetCoordinateSection(dm, coordSection);CHKERRQ(ierr); 7640a6ba040SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr); 7650a6ba040SMatthew G. Knepley ierr = DMDestroy_Plex(dm);CHKERRQ(ierr); 7660a6ba040SMatthew G. Knepley dm->data = dmNew->data; 7670a6ba040SMatthew G. Knepley ((DM_Plex *) dmNew->data)->refct++; 7680a6ba040SMatthew G. Knepley PetscFunctionReturn(0); 7690a6ba040SMatthew G. Knepley } 7700a6ba040SMatthew G. Knepley 7710a6ba040SMatthew G. Knepley #undef __FUNCT__ 7720a6ba040SMatthew G. Knepley #define __FUNCT__ "DMPlexSwap_Static" 7730a6ba040SMatthew G. Knepley /* Swap dm with the contents of dmNew 7740a6ba040SMatthew G. Knepley - Swap the DM_Plex structure 7750a6ba040SMatthew G. Knepley - Swap the coordinates 7760a6ba040SMatthew G. Knepley */ 7770a6ba040SMatthew G. Knepley static PetscErrorCode DMPlexSwap_Static(DM dmA, DM dmB) 7780a6ba040SMatthew G. Knepley { 7790a6ba040SMatthew G. Knepley DM coordDMA, coordDMB; 7800a6ba040SMatthew G. Knepley Vec coordsA, coordsB; 7810a6ba040SMatthew G. Knepley void *tmp; 7820a6ba040SMatthew G. Knepley PetscErrorCode ierr; 7830a6ba040SMatthew G. Knepley 7840a6ba040SMatthew G. Knepley PetscFunctionBegin; 7850a6ba040SMatthew G. Knepley ierr = DMGetCoordinateDM(dmA, &coordDMA);CHKERRQ(ierr); 7860a6ba040SMatthew G. Knepley ierr = DMGetCoordinateDM(dmB, &coordDMB);CHKERRQ(ierr); 7870a6ba040SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) coordDMA);CHKERRQ(ierr); 7880a6ba040SMatthew G. Knepley ierr = DMSetCoordinateDM(dmA, coordDMB);CHKERRQ(ierr); 7890a6ba040SMatthew G. Knepley ierr = DMSetCoordinateDM(dmB, coordDMA);CHKERRQ(ierr); 7900a6ba040SMatthew G. Knepley ierr = PetscObjectDereference((PetscObject) coordDMA);CHKERRQ(ierr); 7910a6ba040SMatthew G. Knepley 7920a6ba040SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dmA, &coordsA);CHKERRQ(ierr); 7930a6ba040SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dmB, &coordsB);CHKERRQ(ierr); 7940a6ba040SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) coordsA);CHKERRQ(ierr); 7950a6ba040SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dmA, coordsB);CHKERRQ(ierr); 7960a6ba040SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dmB, coordsA);CHKERRQ(ierr); 7970a6ba040SMatthew G. Knepley ierr = PetscObjectDereference((PetscObject) coordsA);CHKERRQ(ierr); 7980a6ba040SMatthew G. Knepley tmp = dmA->data; 7990a6ba040SMatthew G. Knepley dmA->data = dmB->data; 8000a6ba040SMatthew G. Knepley dmB->data = tmp; 8010a6ba040SMatthew G. Knepley PetscFunctionReturn(0); 8020a6ba040SMatthew G. Knepley } 8030a6ba040SMatthew G. Knepley 8040a6ba040SMatthew G. Knepley #undef __FUNCT__ 80568d4fef7SMatthew G. Knepley #define __FUNCT__ "DMSetFromOptions_NonRefinement_Plex" 80668d4fef7SMatthew G. Knepley PetscErrorCode DMSetFromOptions_NonRefinement_Plex(DM dm) 8070a6ba040SMatthew G. Knepley { 8080a6ba040SMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 8096c73c22cSMatthew G. Knepley DMBoundary b; 8100a6ba040SMatthew G. Knepley PetscErrorCode ierr; 8110a6ba040SMatthew G. Knepley 8120a6ba040SMatthew G. Knepley PetscFunctionBegin; 8136c73c22cSMatthew G. Knepley /* Handle boundary conditions */ 8146c73c22cSMatthew G. Knepley ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) dm), NULL, "Boundary condition options", "");CHKERRQ(ierr); 8156c73c22cSMatthew G. Knepley for (b = mesh->boundary; b; b = b->next) { 8166c73c22cSMatthew G. Knepley char optname[1024]; 8176c73c22cSMatthew G. Knepley PetscInt ids[1024], len = 1024, i; 8186c73c22cSMatthew G. Knepley PetscBool flg; 8196c73c22cSMatthew G. Knepley 8206c73c22cSMatthew G. Knepley ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr); 8216c73c22cSMatthew G. Knepley ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr); 8226c73c22cSMatthew G. Knepley ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr); 8236c73c22cSMatthew G. Knepley if (flg) { 8246c73c22cSMatthew G. Knepley DMLabel label; 8256c73c22cSMatthew G. Knepley 8266c73c22cSMatthew G. Knepley ierr = DMPlexGetLabel(dm, b->name, &label);CHKERRQ(ierr); 8276c73c22cSMatthew G. Knepley for (i = 0; i < len; ++i) { 8286c73c22cSMatthew G. Knepley PetscBool has; 8296c73c22cSMatthew G. Knepley 8306c73c22cSMatthew G. Knepley ierr = DMLabelHasValue(label, ids[i], &has);CHKERRQ(ierr); 8316c73c22cSMatthew G. Knepley if (!has) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Boundary id %D is not present in the label %s", ids[i], b->name); 8326c73c22cSMatthew G. Knepley } 8336c73c22cSMatthew G. Knepley b->numids = len; 8346c73c22cSMatthew G. Knepley ierr = PetscFree(b->ids);CHKERRQ(ierr); 8356c73c22cSMatthew G. Knepley ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr); 8366c73c22cSMatthew G. Knepley ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr); 8376c73c22cSMatthew G. Knepley } 8386c73c22cSMatthew G. Knepley } 8396c73c22cSMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 8400a6ba040SMatthew G. Knepley /* Handle viewing */ 8410a6ba040SMatthew G. Knepley ierr = PetscOptionsBool("-dm_plex_print_set_values", "Output all set values info", "DMView", PETSC_FALSE, &mesh->printSetValues, NULL);CHKERRQ(ierr); 8420a6ba040SMatthew G. Knepley ierr = PetscOptionsInt("-dm_plex_print_fem", "Debug output level all fem computations", "DMView", 0, &mesh->printFEM, NULL);CHKERRQ(ierr); 8430a6ba040SMatthew G. Knepley ierr = PetscOptionsReal("-dm_plex_print_tol", "Tolerance for FEM output", "DMView", PETSC_FALSE, &mesh->printTol, NULL);CHKERRQ(ierr); 84468d4fef7SMatthew G. Knepley PetscFunctionReturn(0); 84568d4fef7SMatthew G. Knepley } 84668d4fef7SMatthew G. Knepley 84768d4fef7SMatthew G. Knepley #undef __FUNCT__ 84868d4fef7SMatthew G. Knepley #define __FUNCT__ "DMSetFromOptions_Plex" 84968d4fef7SMatthew G. Knepley PetscErrorCode DMSetFromOptions_Plex(DM dm) 85068d4fef7SMatthew G. Knepley { 85168d4fef7SMatthew G. Knepley PetscInt refine = 0, r; 85268d4fef7SMatthew G. Knepley PetscBool isHierarchy; 85368d4fef7SMatthew G. Knepley PetscErrorCode ierr; 85468d4fef7SMatthew G. Knepley 85568d4fef7SMatthew G. Knepley PetscFunctionBegin; 85668d4fef7SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 85768d4fef7SMatthew G. Knepley ierr = PetscOptionsHead("DMPlex Options");CHKERRQ(ierr); 85868d4fef7SMatthew G. Knepley /* Handle DMPlex refinement */ 85968d4fef7SMatthew G. Knepley ierr = PetscOptionsInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL);CHKERRQ(ierr); 86068d4fef7SMatthew G. Knepley ierr = PetscOptionsInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy);CHKERRQ(ierr); 86168d4fef7SMatthew G. Knepley ierr = DMPlexSetRefinementUniform(dm, refine ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 86268d4fef7SMatthew G. Knepley if (refine && isHierarchy) { 86368d4fef7SMatthew G. Knepley DM *dms; 86468d4fef7SMatthew G. Knepley 86568d4fef7SMatthew G. Knepley ierr = PetscMalloc1(refine,&dms);CHKERRQ(ierr); 86668d4fef7SMatthew G. Knepley ierr = DMRefineHierarchy(dm, refine, dms);CHKERRQ(ierr); 86768d4fef7SMatthew G. Knepley /* Total hack since we do not pass in a pointer */ 86868d4fef7SMatthew G. Knepley ierr = DMPlexSwap_Static(dm, dms[refine-1]);CHKERRQ(ierr); 86968d4fef7SMatthew G. Knepley if (refine == 1) { 87068d4fef7SMatthew G. Knepley ierr = DMPlexSetCoarseDM(dm, dms[0]);CHKERRQ(ierr); 87168d4fef7SMatthew G. Knepley } else { 87268d4fef7SMatthew G. Knepley ierr = DMPlexSetCoarseDM(dm, dms[refine-2]);CHKERRQ(ierr); 87368d4fef7SMatthew G. Knepley ierr = DMPlexSetCoarseDM(dms[0], dms[refine-1]);CHKERRQ(ierr); 87468d4fef7SMatthew G. Knepley } 87568d4fef7SMatthew G. Knepley /* Free DMs */ 87668d4fef7SMatthew G. Knepley for (r = 0; r < refine; ++r) { 87768d4fef7SMatthew G. Knepley ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr); 87868d4fef7SMatthew G. Knepley ierr = DMDestroy(&dms[r]);CHKERRQ(ierr); 87968d4fef7SMatthew G. Knepley } 88068d4fef7SMatthew G. Knepley ierr = PetscFree(dms);CHKERRQ(ierr); 88168d4fef7SMatthew G. Knepley } else { 88268d4fef7SMatthew G. Knepley for (r = 0; r < refine; ++r) { 88368d4fef7SMatthew G. Knepley DM refinedMesh; 88468d4fef7SMatthew G. Knepley 88568d4fef7SMatthew G. Knepley ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr); 88668d4fef7SMatthew G. Knepley ierr = DMRefine(dm, PetscObjectComm((PetscObject) dm), &refinedMesh);CHKERRQ(ierr); 88768d4fef7SMatthew G. Knepley /* Total hack since we do not pass in a pointer */ 88868d4fef7SMatthew G. Knepley ierr = DMPlexReplace_Static(dm, refinedMesh);CHKERRQ(ierr); 88968d4fef7SMatthew G. Knepley ierr = DMDestroy(&refinedMesh);CHKERRQ(ierr); 89068d4fef7SMatthew G. Knepley } 89168d4fef7SMatthew G. Knepley } 89268d4fef7SMatthew G. Knepley ierr = DMSetFromOptions_NonRefinement_Plex(dm);CHKERRQ(ierr); 8930a6ba040SMatthew G. Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 8940a6ba040SMatthew G. Knepley PetscFunctionReturn(0); 8950a6ba040SMatthew G. Knepley } 8960a6ba040SMatthew G. Knepley 8970a6ba040SMatthew G. Knepley #undef __FUNCT__ 898552f7358SJed Brown #define __FUNCT__ "DMCreateGlobalVector_Plex" 899552f7358SJed Brown static PetscErrorCode DMCreateGlobalVector_Plex(DM dm,Vec *vec) 900552f7358SJed Brown { 901552f7358SJed Brown PetscErrorCode ierr; 902552f7358SJed Brown 903552f7358SJed Brown PetscFunctionBegin; 904552f7358SJed Brown ierr = DMCreateGlobalVector_Section_Private(dm,vec);CHKERRQ(ierr); 905552f7358SJed Brown /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 906552f7358SJed Brown ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Plex);CHKERRQ(ierr); 907552f7358SJed Brown PetscFunctionReturn(0); 908552f7358SJed Brown } 909552f7358SJed Brown 910552f7358SJed Brown #undef __FUNCT__ 911552f7358SJed Brown #define __FUNCT__ "DMCreateLocalVector_Plex" 912552f7358SJed Brown static PetscErrorCode DMCreateLocalVector_Plex(DM dm,Vec *vec) 913552f7358SJed Brown { 914552f7358SJed Brown PetscErrorCode ierr; 915552f7358SJed Brown 916552f7358SJed Brown PetscFunctionBegin; 917552f7358SJed Brown ierr = DMCreateLocalVector_Section_Private(dm,vec);CHKERRQ(ierr); 918552f7358SJed Brown ierr = VecSetOperation(*vec, VECOP_VIEW, (void(*)(void)) VecView_Plex_Local);CHKERRQ(ierr); 919552f7358SJed Brown PetscFunctionReturn(0); 920552f7358SJed Brown } 921552f7358SJed Brown 92238221697SMatthew G. Knepley #undef __FUNCT__ 923552f7358SJed Brown #define __FUNCT__ "DMInitialize_Plex" 924552f7358SJed Brown PetscErrorCode DMInitialize_Plex(DM dm) 925552f7358SJed Brown { 926552f7358SJed Brown PetscFunctionBegin; 927552f7358SJed Brown dm->ops->view = DMView_Plex; 928552f7358SJed Brown dm->ops->setfromoptions = DMSetFromOptions_Plex; 92938221697SMatthew G. Knepley dm->ops->clone = DMClone_Plex; 930552f7358SJed Brown dm->ops->setup = DMSetUp_Plex; 931fd59a867SMatthew G. Knepley dm->ops->createdefaultsection = DMCreateDefaultSection_Plex; 932552f7358SJed Brown dm->ops->createglobalvector = DMCreateGlobalVector_Plex; 933552f7358SJed Brown dm->ops->createlocalvector = DMCreateLocalVector_Plex; 934184d77edSJed Brown dm->ops->getlocaltoglobalmapping = NULL; 935184d77edSJed Brown dm->ops->getlocaltoglobalmappingblock = NULL; 9360298fd71SBarry Smith dm->ops->createfieldis = NULL; 937552f7358SJed Brown dm->ops->createcoordinatedm = DMCreateCoordinateDM_Plex; 9380a6ba040SMatthew G. Knepley dm->ops->getcoloring = NULL; 939552f7358SJed Brown dm->ops->creatematrix = DMCreateMatrix_Plex; 940bceba477SMatthew G. Knepley dm->ops->createinterpolation = DMCreateInterpolation_Plex; 941bceba477SMatthew G. Knepley dm->ops->getaggregates = NULL; 942bceba477SMatthew G. Knepley dm->ops->getinjection = DMCreateInjection_Plex; 943552f7358SJed Brown dm->ops->refine = DMRefine_Plex; 9440a6ba040SMatthew G. Knepley dm->ops->coarsen = DMCoarsen_Plex; 9450a6ba040SMatthew G. Knepley dm->ops->refinehierarchy = DMRefineHierarchy_Plex; 9460a6ba040SMatthew G. Knepley dm->ops->coarsenhierarchy = NULL; 9470298fd71SBarry Smith dm->ops->globaltolocalbegin = NULL; 9480298fd71SBarry Smith dm->ops->globaltolocalend = NULL; 9490298fd71SBarry Smith dm->ops->localtoglobalbegin = NULL; 9500298fd71SBarry Smith dm->ops->localtoglobalend = NULL; 951552f7358SJed Brown dm->ops->destroy = DMDestroy_Plex; 952552f7358SJed Brown dm->ops->createsubdm = DMCreateSubDM_Plex; 953552f7358SJed Brown dm->ops->locatepoints = DMLocatePoints_Plex; 954552f7358SJed Brown PetscFunctionReturn(0); 955552f7358SJed Brown } 956552f7358SJed Brown 95763a16f15SMatthew G. Knepley #undef __FUNCT__ 95863a16f15SMatthew G. Knepley #define __FUNCT__ "DMClone_Plex" 95963a16f15SMatthew G. Knepley PetscErrorCode DMClone_Plex(DM dm, DM *newdm) 96063a16f15SMatthew G. Knepley { 96163a16f15SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 96263a16f15SMatthew G. Knepley PetscErrorCode ierr; 96363a16f15SMatthew G. Knepley 96463a16f15SMatthew G. Knepley PetscFunctionBegin; 96563a16f15SMatthew G. Knepley mesh->refct++; 96663a16f15SMatthew G. Knepley (*newdm)->data = mesh; 96763a16f15SMatthew G. Knepley ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMPLEX);CHKERRQ(ierr); 96863a16f15SMatthew G. Knepley ierr = DMInitialize_Plex(*newdm);CHKERRQ(ierr); 96963a16f15SMatthew G. Knepley PetscFunctionReturn(0); 97063a16f15SMatthew G. Knepley } 97163a16f15SMatthew G. Knepley 9728818961aSMatthew G Knepley /*MC 9738818961aSMatthew G Knepley DMPLEX = "plex" - A DM object that encapsulates an unstructured mesh, or CW Complex, which can be expressed using a Hasse Diagram. 9748818961aSMatthew G Knepley In the local representation, Vecs contain all unknowns in the interior and shared boundary. This is 9758818961aSMatthew G Knepley specified by a PetscSection object. Ownership in the global representation is determined by 9768818961aSMatthew G Knepley ownership of the underlying DMPlex points. This is specified by another PetscSection object. 9778818961aSMatthew G Knepley 9788818961aSMatthew G Knepley Level: intermediate 9798818961aSMatthew G Knepley 9808818961aSMatthew G Knepley .seealso: DMType, DMPlexCreate(), DMCreate(), DMSetType() 9818818961aSMatthew G Knepley M*/ 9828818961aSMatthew G Knepley 983552f7358SJed Brown #undef __FUNCT__ 984552f7358SJed Brown #define __FUNCT__ "DMCreate_Plex" 9858cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Plex(DM dm) 986552f7358SJed Brown { 987552f7358SJed Brown DM_Plex *mesh; 988770b213bSMatthew G Knepley PetscInt unit, d; 989552f7358SJed Brown PetscErrorCode ierr; 990552f7358SJed Brown 991552f7358SJed Brown PetscFunctionBegin; 992552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 993b00a9115SJed Brown ierr = PetscNewLog(dm,&mesh);CHKERRQ(ierr); 994552f7358SJed Brown dm->data = mesh; 995552f7358SJed Brown 996552f7358SJed Brown mesh->refct = 1; 997552f7358SJed Brown mesh->dim = 0; 99882f516ccSBarry Smith ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->coneSection);CHKERRQ(ierr); 999552f7358SJed Brown mesh->maxConeSize = 0; 10000298fd71SBarry Smith mesh->cones = NULL; 10010298fd71SBarry Smith mesh->coneOrientations = NULL; 100282f516ccSBarry Smith ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &mesh->supportSection);CHKERRQ(ierr); 1003552f7358SJed Brown mesh->maxSupportSize = 0; 10040298fd71SBarry Smith mesh->supports = NULL; 1005552f7358SJed Brown mesh->refinementUniform = PETSC_TRUE; 1006552f7358SJed Brown mesh->refinementLimit = -1.0; 1007552f7358SJed Brown 10080298fd71SBarry Smith mesh->facesTmp = NULL; 1009552f7358SJed Brown 10100298fd71SBarry Smith mesh->subpointMap = NULL; 1011552f7358SJed Brown 10128865f1eaSKarl Rupp for (unit = 0; unit < NUM_PETSC_UNITS; ++unit) mesh->scale[unit] = 1.0; 1013552f7358SJed Brown 10140298fd71SBarry Smith mesh->labels = NULL; 10158af19771SMatthew G. Knepley mesh->depthLabel = NULL; 10160298fd71SBarry Smith mesh->globalVertexNumbers = NULL; 10170298fd71SBarry Smith mesh->globalCellNumbers = NULL; 10188865f1eaSKarl Rupp for (d = 0; d < 8; ++d) mesh->hybridPointMax[d] = PETSC_DETERMINE; 1019552f7358SJed Brown mesh->vtkCellHeight = 0; 1020552f7358SJed Brown mesh->preallocCenterDim = -1; 1021552f7358SJed Brown 1022552f7358SJed Brown mesh->printSetValues = PETSC_FALSE; 1023552f7358SJed Brown mesh->printFEM = 0; 10246113b454SMatthew G. Knepley mesh->printTol = 1.0e-10; 1025552f7358SJed Brown 1026552f7358SJed Brown ierr = DMInitialize_Plex(dm);CHKERRQ(ierr); 1027552f7358SJed Brown PetscFunctionReturn(0); 1028552f7358SJed Brown } 1029552f7358SJed Brown 1030552f7358SJed Brown #undef __FUNCT__ 1031552f7358SJed Brown #define __FUNCT__ "DMPlexCreate" 1032552f7358SJed Brown /*@ 1033552f7358SJed Brown DMPlexCreate - Creates a DMPlex object, which encapsulates an unstructured mesh, or CW complex, which can be expressed using a Hasse Diagram. 1034552f7358SJed Brown 1035552f7358SJed Brown Collective on MPI_Comm 1036552f7358SJed Brown 1037552f7358SJed Brown Input Parameter: 1038552f7358SJed Brown . comm - The communicator for the DMPlex object 1039552f7358SJed Brown 1040552f7358SJed Brown Output Parameter: 1041552f7358SJed Brown . mesh - The DMPlex object 1042552f7358SJed Brown 1043552f7358SJed Brown Level: beginner 1044552f7358SJed Brown 1045552f7358SJed Brown .keywords: DMPlex, create 1046552f7358SJed Brown @*/ 1047552f7358SJed Brown PetscErrorCode DMPlexCreate(MPI_Comm comm, DM *mesh) 1048552f7358SJed Brown { 1049552f7358SJed Brown PetscErrorCode ierr; 1050552f7358SJed Brown 1051552f7358SJed Brown PetscFunctionBegin; 1052552f7358SJed Brown PetscValidPointer(mesh,2); 1053552f7358SJed Brown ierr = DMCreate(comm, mesh);CHKERRQ(ierr); 1054552f7358SJed Brown ierr = DMSetType(*mesh, DMPLEX);CHKERRQ(ierr); 1055552f7358SJed Brown PetscFunctionReturn(0); 1056552f7358SJed Brown } 1057552f7358SJed Brown 1058552f7358SJed Brown #undef __FUNCT__ 10599298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildFromCellList_Private" 10609298eaa6SMatthew G Knepley /* 10619298eaa6SMatthew G Knepley This takes as input the common mesh generator output, a list of the vertices for each cell 10629298eaa6SMatthew G Knepley */ 10639298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildFromCellList_Private(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, const int cells[]) 10649298eaa6SMatthew G Knepley { 10659298eaa6SMatthew G Knepley PetscInt *cone, c, p; 10669298eaa6SMatthew G Knepley PetscErrorCode ierr; 10679298eaa6SMatthew G Knepley 10689298eaa6SMatthew G Knepley PetscFunctionBegin; 10699298eaa6SMatthew G Knepley ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr); 10709298eaa6SMatthew G Knepley for (c = 0; c < numCells; ++c) { 10719298eaa6SMatthew G Knepley ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr); 10729298eaa6SMatthew G Knepley } 10739298eaa6SMatthew G Knepley ierr = DMSetUp(dm);CHKERRQ(ierr); 10749298eaa6SMatthew G Knepley ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 10759298eaa6SMatthew G Knepley for (c = 0; c < numCells; ++c) { 10769298eaa6SMatthew G Knepley for (p = 0; p < numCorners; ++p) { 10779298eaa6SMatthew G Knepley cone[p] = cells[c*numCorners+p]+numCells; 10789298eaa6SMatthew G Knepley } 10799298eaa6SMatthew G Knepley ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr); 10809298eaa6SMatthew G Knepley } 10819298eaa6SMatthew G Knepley ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr); 10829298eaa6SMatthew G Knepley ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 10839298eaa6SMatthew G Knepley ierr = DMPlexStratify(dm);CHKERRQ(ierr); 10849298eaa6SMatthew G Knepley PetscFunctionReturn(0); 10859298eaa6SMatthew G Knepley } 10869298eaa6SMatthew G Knepley 10879298eaa6SMatthew G Knepley #undef __FUNCT__ 10889298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexBuildCoordinates_Private" 10899298eaa6SMatthew G Knepley /* 10909298eaa6SMatthew G Knepley This takes as input the coordinates for each vertex 10919298eaa6SMatthew G Knepley */ 10929298eaa6SMatthew G Knepley PetscErrorCode DMPlexBuildCoordinates_Private(DM dm, PetscInt spaceDim, PetscInt numCells, PetscInt numVertices, const double vertexCoords[]) 10939298eaa6SMatthew G Knepley { 10949298eaa6SMatthew G Knepley PetscSection coordSection; 10959298eaa6SMatthew G Knepley Vec coordinates; 10969298eaa6SMatthew G Knepley PetscScalar *coords; 10979298eaa6SMatthew G Knepley PetscInt coordSize, v, d; 10989298eaa6SMatthew G Knepley PetscErrorCode ierr; 10999298eaa6SMatthew G Knepley 11009298eaa6SMatthew G Knepley PetscFunctionBegin; 1101c2166f76SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 11029298eaa6SMatthew G Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 11039298eaa6SMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, spaceDim);CHKERRQ(ierr); 11049298eaa6SMatthew G Knepley ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 11059298eaa6SMatthew G Knepley for (v = numCells; v < numCells+numVertices; ++v) { 11069298eaa6SMatthew G Knepley ierr = PetscSectionSetDof(coordSection, v, spaceDim);CHKERRQ(ierr); 11079298eaa6SMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, spaceDim);CHKERRQ(ierr); 11089298eaa6SMatthew G Knepley } 11099298eaa6SMatthew G Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 11109298eaa6SMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 11119298eaa6SMatthew G Knepley ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 11129298eaa6SMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 11139298eaa6SMatthew G Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 11142eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 11159298eaa6SMatthew G Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 11169298eaa6SMatthew G Knepley for (v = 0; v < numVertices; ++v) { 11179298eaa6SMatthew G Knepley for (d = 0; d < spaceDim; ++d) { 11189298eaa6SMatthew G Knepley coords[v*spaceDim+d] = vertexCoords[v*spaceDim+d]; 11199298eaa6SMatthew G Knepley } 11209298eaa6SMatthew G Knepley } 11219298eaa6SMatthew G Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 11229298eaa6SMatthew G Knepley ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 11239298eaa6SMatthew G Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 11249298eaa6SMatthew G Knepley PetscFunctionReturn(0); 11259298eaa6SMatthew G Knepley } 11269298eaa6SMatthew G Knepley 11279298eaa6SMatthew G Knepley #undef __FUNCT__ 11289298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromCellList" 11299298eaa6SMatthew G Knepley /*@C 11309298eaa6SMatthew G Knepley DMPlexCreateFromCellList - This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM 11319298eaa6SMatthew G Knepley 11329298eaa6SMatthew G Knepley Input Parameters: 11339298eaa6SMatthew G Knepley + comm - The communicator 11349298eaa6SMatthew G Knepley . dim - The topological dimension of the mesh 11359298eaa6SMatthew G Knepley . numCells - The number of cells 11369298eaa6SMatthew G Knepley . numVertices - The number of vertices 11379298eaa6SMatthew G Knepley . numCorners - The number of vertices for each cell 11389298eaa6SMatthew G Knepley . interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically 11399298eaa6SMatthew G Knepley . cells - An array of numCells*numCorners numbers, the vertices for each cell 11409298eaa6SMatthew G Knepley . spaceDim - The spatial dimension used for coordinates 11419298eaa6SMatthew G Knepley - vertexCoords - An array of numVertices*spaceDim numbers, the coordinates of each vertex 11429298eaa6SMatthew G Knepley 11439298eaa6SMatthew G Knepley Output Parameter: 11449298eaa6SMatthew G Knepley . dm - The DM 11459298eaa6SMatthew G Knepley 11469298eaa6SMatthew G Knepley Note: Two triangles sharing a face 11479298eaa6SMatthew G Knepley $ 11489298eaa6SMatthew G Knepley $ 2 11499298eaa6SMatthew G Knepley $ / | \ 11509298eaa6SMatthew G Knepley $ / | \ 11519298eaa6SMatthew G Knepley $ / | \ 11529298eaa6SMatthew G Knepley $ 0 0 | 1 3 11539298eaa6SMatthew G Knepley $ \ | / 11549298eaa6SMatthew G Knepley $ \ | / 11559298eaa6SMatthew G Knepley $ \ | / 11569298eaa6SMatthew G Knepley $ 1 11579298eaa6SMatthew G Knepley would have input 11589298eaa6SMatthew G Knepley $ numCells = 2, numVertices = 4 11599298eaa6SMatthew G Knepley $ cells = [0 1 2 1 3 2] 11609298eaa6SMatthew G Knepley $ 11619298eaa6SMatthew G Knepley which would result in the DMPlex 11629298eaa6SMatthew G Knepley $ 11639298eaa6SMatthew G Knepley $ 4 11649298eaa6SMatthew G Knepley $ / | \ 11659298eaa6SMatthew G Knepley $ / | \ 11669298eaa6SMatthew G Knepley $ / | \ 11679298eaa6SMatthew G Knepley $ 2 0 | 1 5 11689298eaa6SMatthew G Knepley $ \ | / 11699298eaa6SMatthew G Knepley $ \ | / 11709298eaa6SMatthew G Knepley $ \ | / 11719298eaa6SMatthew G Knepley $ 3 11729298eaa6SMatthew G Knepley 11739298eaa6SMatthew G Knepley Level: beginner 11749298eaa6SMatthew G Knepley 1175939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromDAG(), DMPlexCreate() 11769298eaa6SMatthew G Knepley @*/ 11779298eaa6SMatthew 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) 11789298eaa6SMatthew G Knepley { 11799298eaa6SMatthew G Knepley PetscErrorCode ierr; 11809298eaa6SMatthew G Knepley 11819298eaa6SMatthew G Knepley PetscFunctionBegin; 11829298eaa6SMatthew G Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 11839298eaa6SMatthew G Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 11849298eaa6SMatthew G Knepley ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 11859298eaa6SMatthew G Knepley ierr = DMPlexBuildFromCellList_Private(*dm, numCells, numVertices, numCorners, cells);CHKERRQ(ierr); 11869298eaa6SMatthew G Knepley if (interpolate) { 11879298eaa6SMatthew G Knepley DM idm; 11889298eaa6SMatthew G Knepley 11899298eaa6SMatthew G Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 11909298eaa6SMatthew G Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 11919298eaa6SMatthew G Knepley *dm = idm; 11929298eaa6SMatthew G Knepley } 11939298eaa6SMatthew G Knepley ierr = DMPlexBuildCoordinates_Private(*dm, spaceDim, numCells, numVertices, vertexCoords);CHKERRQ(ierr); 11949298eaa6SMatthew G Knepley PetscFunctionReturn(0); 11959298eaa6SMatthew G Knepley } 11969298eaa6SMatthew G Knepley 11979298eaa6SMatthew G Knepley #undef __FUNCT__ 11989298eaa6SMatthew G Knepley #define __FUNCT__ "DMPlexCreateFromDAG" 1199939f6067SMatthew G. Knepley /*@ 1200939f6067SMatthew 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 1201939f6067SMatthew G. Knepley 1202939f6067SMatthew G. Knepley Input Parameters: 1203939f6067SMatthew G. Knepley + dm - The empty DM object, usually from DMCreate() and DMPlexSetDimension() 1204939f6067SMatthew G. Knepley . depth - The depth of the DAG 1205939f6067SMatthew G. Knepley . numPoints - The number of points at each depth 1206939f6067SMatthew G. Knepley . coneSize - The cone size of each point 1207939f6067SMatthew G. Knepley . cones - The concatenation of the cone points for each point, the cone list must be oriented correctly for each point 1208939f6067SMatthew G. Knepley . coneOrientations - The orientation of each cone point 1209939f6067SMatthew G. Knepley - vertexCoords - An array of numVertices*dim numbers, the coordinates of each vertex 1210939f6067SMatthew G. Knepley 1211939f6067SMatthew G. Knepley Output Parameter: 1212939f6067SMatthew G. Knepley . dm - The DM 1213939f6067SMatthew G. Knepley 1214939f6067SMatthew G. Knepley Note: Two triangles sharing a face would have input 1215939f6067SMatthew G. Knepley $ depth = 1, numPoints = [4 2], coneSize = [3 3 0 0 0 0] 1216939f6067SMatthew G. Knepley $ cones = [2 3 4 3 5 4], coneOrientations = [0 0 0 0 0 0] 1217939f6067SMatthew G. Knepley $ vertexCoords = [-1.0 0.0 0.0 -1.0 0.0 1.0 1.0 0.0] 1218939f6067SMatthew G. Knepley $ 1219939f6067SMatthew G. Knepley which would result in the DMPlex 1220939f6067SMatthew G. Knepley $ 1221939f6067SMatthew G. Knepley $ 4 1222939f6067SMatthew G. Knepley $ / | \ 1223939f6067SMatthew G. Knepley $ / | \ 1224939f6067SMatthew G. Knepley $ / | \ 1225939f6067SMatthew G. Knepley $ 2 0 | 1 5 1226939f6067SMatthew G. Knepley $ \ | / 1227939f6067SMatthew G. Knepley $ \ | / 1228939f6067SMatthew G. Knepley $ \ | / 1229939f6067SMatthew G. Knepley $ 3 1230939f6067SMatthew G. Knepley $ 1231939f6067SMatthew G. Knepley $ Notice that all points are numbered consecutively, unlikely DMPlexCreateFromCellList() 1232939f6067SMatthew G. Knepley 1233939f6067SMatthew G. Knepley Level: advanced 1234939f6067SMatthew G. Knepley 1235939f6067SMatthew G. Knepley .seealso: DMPlexCreateFromCellList(), DMPlexCreate() 1236939f6067SMatthew G. Knepley @*/ 12379298eaa6SMatthew G Knepley PetscErrorCode DMPlexCreateFromDAG(DM dm, PetscInt depth, const PetscInt numPoints[], const PetscInt coneSize[], const PetscInt cones[], const PetscInt coneOrientations[], const PetscScalar vertexCoords[]) 12389298eaa6SMatthew G Knepley { 12399298eaa6SMatthew G Knepley Vec coordinates; 12409298eaa6SMatthew G Knepley PetscSection coordSection; 12419298eaa6SMatthew G Knepley PetscScalar *coords; 12429298eaa6SMatthew G Knepley PetscInt coordSize, firstVertex = numPoints[depth], pStart = 0, pEnd = 0, p, v, dim, d, off; 12439298eaa6SMatthew G Knepley PetscErrorCode ierr; 12449298eaa6SMatthew G Knepley 12459298eaa6SMatthew G Knepley PetscFunctionBegin; 12469298eaa6SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 12479298eaa6SMatthew G Knepley for (d = 0; d <= depth; ++d) pEnd += numPoints[d]; 12489298eaa6SMatthew G Knepley ierr = DMPlexSetChart(dm, pStart, pEnd);CHKERRQ(ierr); 12499298eaa6SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 12509298eaa6SMatthew G Knepley ierr = DMPlexSetConeSize(dm, p, coneSize[p-pStart]);CHKERRQ(ierr); 12519298eaa6SMatthew G Knepley } 12529298eaa6SMatthew G Knepley ierr = DMSetUp(dm);CHKERRQ(ierr); /* Allocate space for cones */ 12539298eaa6SMatthew G Knepley for (p = pStart, off = 0; p < pEnd; off += coneSize[p-pStart], ++p) { 12549298eaa6SMatthew G Knepley ierr = DMPlexSetCone(dm, p, &cones[off]);CHKERRQ(ierr); 12559298eaa6SMatthew G Knepley ierr = DMPlexSetConeOrientation(dm, p, &coneOrientations[off]);CHKERRQ(ierr); 12569298eaa6SMatthew G Knepley } 12579298eaa6SMatthew G Knepley ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr); 12589298eaa6SMatthew G Knepley ierr = DMPlexStratify(dm);CHKERRQ(ierr); 12599298eaa6SMatthew G Knepley /* Build coordinates */ 1260c2166f76SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 12619298eaa6SMatthew G Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 12629298eaa6SMatthew G Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 12639298eaa6SMatthew G Knepley ierr = PetscSectionSetChart(coordSection, firstVertex, firstVertex+numPoints[0]);CHKERRQ(ierr); 12649298eaa6SMatthew G Knepley for (v = firstVertex; v < firstVertex+numPoints[0]; ++v) { 12659298eaa6SMatthew G Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 12666f8cbbeeSMatthew G Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 12679298eaa6SMatthew G Knepley } 12689298eaa6SMatthew G Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 12699298eaa6SMatthew G Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 12709298eaa6SMatthew G Knepley ierr = VecCreate(PetscObjectComm((PetscObject)dm), &coordinates);CHKERRQ(ierr); 12716f8cbbeeSMatthew G Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 12729298eaa6SMatthew G Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 12732eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 12749298eaa6SMatthew G Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 12759298eaa6SMatthew G Knepley for (v = 0; v < numPoints[0]; ++v) { 12769298eaa6SMatthew G Knepley PetscInt off; 12779298eaa6SMatthew G Knepley 12789298eaa6SMatthew G Knepley ierr = PetscSectionGetOffset(coordSection, v+firstVertex, &off);CHKERRQ(ierr); 12799298eaa6SMatthew G Knepley for (d = 0; d < dim; ++d) { 12809298eaa6SMatthew G Knepley coords[off+d] = vertexCoords[v*dim+d]; 12819298eaa6SMatthew G Knepley } 12829298eaa6SMatthew G Knepley } 12839298eaa6SMatthew G Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 12849298eaa6SMatthew G Knepley ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 12859298eaa6SMatthew G Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 12869298eaa6SMatthew G Knepley PetscFunctionReturn(0); 12879298eaa6SMatthew G Knepley } 1288