xref: /petsc/src/dm/impls/plex/plexinterpolate.c (revision 07b0a1fcbb200e7259840b2fcdb7d753db089615)
19f074e33SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
29f074e33SMatthew G Knepley #include <../src/sys/utils/hash.h>
39f074e33SMatthew G Knepley 
49f074e33SMatthew G Knepley #undef __FUNCT__
59f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexGetFaces"
69f074e33SMatthew G Knepley /*
79f074e33SMatthew G Knepley   DMPlexGetFaces -
89f074e33SMatthew G Knepley 
99f074e33SMatthew G Knepley   Note: This will only work for cell-vertex meshes.
109f074e33SMatthew G Knepley */
119f074e33SMatthew G Knepley static PetscErrorCode DMPlexGetFaces(DM dm, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[])
129f074e33SMatthew G Knepley {
139f074e33SMatthew G Knepley   DM_Plex        *mesh = (DM_Plex*) dm->data;
149f074e33SMatthew G Knepley   const PetscInt *cone = NULL;
159f074e33SMatthew G Knepley   PetscInt        depth = 0, dim, coneSize;
169f074e33SMatthew G Knepley   PetscErrorCode  ierr;
179f074e33SMatthew G Knepley 
189f074e33SMatthew G Knepley   PetscFunctionBegin;
199f074e33SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
209f074e33SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
219f074e33SMatthew G Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
229f074e33SMatthew G Knepley   if (depth > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Faces can only be returned for cell-vertex meshes.");
239f074e33SMatthew G Knepley   if (!mesh->facesTmp) {ierr = PetscMalloc(PetscSqr(PetscMax(mesh->maxConeSize, mesh->maxSupportSize)) * sizeof(PetscInt), &mesh->facesTmp);CHKERRQ(ierr);}
249f074e33SMatthew G Knepley   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
259f074e33SMatthew G Knepley   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
269f074e33SMatthew G Knepley   switch (dim) {
279f074e33SMatthew G Knepley   case 2:
289f074e33SMatthew G Knepley     switch (coneSize) {
299f074e33SMatthew G Knepley     case 3:
309f074e33SMatthew G Knepley       mesh->facesTmp[0] = cone[0]; mesh->facesTmp[1] = cone[1];
319f074e33SMatthew G Knepley       mesh->facesTmp[2] = cone[1]; mesh->facesTmp[3] = cone[2];
329f074e33SMatthew G Knepley       mesh->facesTmp[4] = cone[2]; mesh->facesTmp[5] = cone[0];
339f074e33SMatthew G Knepley       *numFaces         = 3;
349f074e33SMatthew G Knepley       *faceSize         = 2;
359f074e33SMatthew G Knepley       *faces            = mesh->facesTmp;
369f074e33SMatthew G Knepley       break;
379f074e33SMatthew G Knepley     case 4:
389f074e33SMatthew G Knepley       mesh->facesTmp[0] = cone[0]; mesh->facesTmp[1] = cone[1];
399f074e33SMatthew G Knepley       mesh->facesTmp[2] = cone[1]; mesh->facesTmp[3] = cone[2];
409f074e33SMatthew G Knepley       mesh->facesTmp[4] = cone[2]; mesh->facesTmp[5] = cone[3];
419f074e33SMatthew G Knepley       mesh->facesTmp[6] = cone[3]; mesh->facesTmp[7] = cone[0];
429f074e33SMatthew G Knepley       *numFaces         = 4;
439f074e33SMatthew G Knepley       *faceSize         = 2;
449f074e33SMatthew G Knepley       *faces            = mesh->facesTmp;
459f074e33SMatthew G Knepley       break;
469f074e33SMatthew G Knepley     default:
479f074e33SMatthew G Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
489f074e33SMatthew G Knepley     }
499f074e33SMatthew G Knepley     break;
509f074e33SMatthew G Knepley   case 3:
519f074e33SMatthew G Knepley     switch (coneSize) {
529f074e33SMatthew G Knepley     case 3:
539f074e33SMatthew G Knepley       mesh->facesTmp[0] = cone[0]; mesh->facesTmp[1] = cone[1];
549f074e33SMatthew G Knepley       mesh->facesTmp[2] = cone[1]; mesh->facesTmp[3] = cone[2];
559f074e33SMatthew G Knepley       mesh->facesTmp[4] = cone[2]; mesh->facesTmp[5] = cone[0];
569f074e33SMatthew G Knepley       *numFaces         = 3;
579f074e33SMatthew G Knepley       *faceSize         = 2;
589f074e33SMatthew G Knepley       *faces            = mesh->facesTmp;
599f074e33SMatthew G Knepley       break;
609f074e33SMatthew G Knepley     case 4:
619f074e33SMatthew G Knepley       mesh->facesTmp[0] = cone[0]; mesh->facesTmp[1]  = cone[1]; mesh->facesTmp[2]  = cone[2];
629f074e33SMatthew G Knepley       mesh->facesTmp[3] = cone[0]; mesh->facesTmp[4]  = cone[2]; mesh->facesTmp[5]  = cone[3];
639f074e33SMatthew G Knepley       mesh->facesTmp[6] = cone[0]; mesh->facesTmp[7]  = cone[3]; mesh->facesTmp[8]  = cone[1];
649f074e33SMatthew G Knepley       mesh->facesTmp[9] = cone[1]; mesh->facesTmp[10] = cone[3]; mesh->facesTmp[11] = cone[2];
659f074e33SMatthew G Knepley       *numFaces         = 4;
669f074e33SMatthew G Knepley       *faceSize         = 3;
679f074e33SMatthew G Knepley       *faces            = mesh->facesTmp;
689f074e33SMatthew G Knepley       break;
699f074e33SMatthew G Knepley     default:
709f074e33SMatthew G Knepley       SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim);
719f074e33SMatthew G Knepley     }
729f074e33SMatthew G Knepley     break;
739f074e33SMatthew G Knepley   default:
749f074e33SMatthew G Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
759f074e33SMatthew G Knepley   }
769f074e33SMatthew G Knepley   PetscFunctionReturn(0);
779f074e33SMatthew G Knepley }
789f074e33SMatthew G Knepley 
799f074e33SMatthew G Knepley #undef __FUNCT__
809f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolate_2D"
819f074e33SMatthew G Knepley static PetscErrorCode DMPlexInterpolate_2D(DM dm, DM *dmInt)
829f074e33SMatthew G Knepley {
839f074e33SMatthew G Knepley   DM             idm;
849f074e33SMatthew G Knepley   DM_Plex       *mesh;
859f074e33SMatthew G Knepley   PetscHashIJ    edgeTable;
869f074e33SMatthew G Knepley   PetscInt      *off;
879f074e33SMatthew G Knepley   PetscInt       dim, numCells, cStart, cEnd, c, numVertices, vStart, vEnd;
889f074e33SMatthew G Knepley   PetscInt       numEdges, firstEdge, edge, e;
899f074e33SMatthew G Knepley   PetscErrorCode ierr;
909f074e33SMatthew G Knepley 
919f074e33SMatthew G Knepley   PetscFunctionBegin;
929f074e33SMatthew G Knepley   ierr        = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
939f074e33SMatthew G Knepley   ierr        = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
949f074e33SMatthew G Knepley   ierr        = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
959f074e33SMatthew G Knepley   numCells    = cEnd - cStart;
969f074e33SMatthew G Knepley   numVertices = vEnd - vStart;
979f074e33SMatthew G Knepley   firstEdge   = numCells + numVertices;
989f074e33SMatthew G Knepley   numEdges    = 0;
999f074e33SMatthew G Knepley   /* Count edges using algorithm from CreateNeighborCSR */
1009f074e33SMatthew G Knepley   ierr = DMPlexCreateNeighborCSR(dm, NULL, &off, NULL);CHKERRQ(ierr);
1019f074e33SMatthew G Knepley   if (off) {
1029f074e33SMatthew G Knepley     PetscInt numCorners = 0;
1039f074e33SMatthew G Knepley 
1049f074e33SMatthew G Knepley     numEdges = off[numCells]/2;
1059f074e33SMatthew G Knepley #if 0
1069f074e33SMatthew G Knepley     /* Account for boundary edges: \sum_c 3 - neighbors = 3*numCells - totalNeighbors */
1079f074e33SMatthew G Knepley     numEdges += 3*numCells - off[numCells];
1089f074e33SMatthew G Knepley #else
1099f074e33SMatthew G Knepley     /* Account for boundary edges: \sum_c #faces - #neighbors = \sum_c #cellVertices - #neighbors = totalCorners - totalNeighbors */
1109f074e33SMatthew G Knepley     for (c = cStart; c < cEnd; ++c) {
1119f074e33SMatthew G Knepley       PetscInt coneSize;
1129f074e33SMatthew G Knepley 
1139f074e33SMatthew G Knepley       ierr        = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
1149f074e33SMatthew G Knepley       numCorners += coneSize;
1159f074e33SMatthew G Knepley     }
1169f074e33SMatthew G Knepley     numEdges += numCorners - off[numCells];
1179f074e33SMatthew G Knepley #endif
1189f074e33SMatthew G Knepley   }
1199f074e33SMatthew G Knepley #if 0
1209f074e33SMatthew G Knepley   /* Check Euler characteristic V - E + F = 1 */
1219f074e33SMatthew G Knepley   if (numVertices && (numVertices-numEdges+numCells != 1)) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Euler characteristic of mesh is %d  != 1", numVertices-numEdges+numCells);
1229f074e33SMatthew G Knepley #endif
1239f074e33SMatthew G Knepley   /* Create interpolated mesh */
1249f074e33SMatthew G Knepley   ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);
1259f074e33SMatthew G Knepley   ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
1269f074e33SMatthew G Knepley   ierr = DMPlexSetDimension(idm, dim);CHKERRQ(ierr);
1279f074e33SMatthew G Knepley   ierr = DMPlexSetChart(idm, 0, numCells+numVertices+numEdges);CHKERRQ(ierr);
1289f074e33SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
1299f074e33SMatthew G Knepley     PetscInt numCorners;
1309f074e33SMatthew G Knepley 
1319f074e33SMatthew G Knepley     ierr = DMPlexGetConeSize(dm, c, &numCorners);CHKERRQ(ierr);
1329f074e33SMatthew G Knepley     ierr = DMPlexSetConeSize(idm, c, numCorners);CHKERRQ(ierr);
1339f074e33SMatthew G Knepley   }
1349f074e33SMatthew G Knepley   for (e = firstEdge; e < firstEdge+numEdges; ++e) {
1359f074e33SMatthew G Knepley     ierr = DMPlexSetConeSize(idm, e, 2);CHKERRQ(ierr);
1369f074e33SMatthew G Knepley   }
1379f074e33SMatthew G Knepley   ierr = DMSetUp(idm);CHKERRQ(ierr);
1389f074e33SMatthew G Knepley   /* Get edge cones from subsets of cell vertices */
1399f074e33SMatthew G Knepley   ierr = PetscHashIJCreate(&edgeTable);CHKERRQ(ierr);
1409f074e33SMatthew G Knepley   ierr = PetscHashIJSetMultivalued(edgeTable, PETSC_FALSE);CHKERRQ(ierr);
1419f074e33SMatthew G Knepley 
1429f074e33SMatthew G Knepley   for (c = 0, edge = firstEdge; c < numCells; ++c) {
1439f074e33SMatthew G Knepley     const PetscInt *cellFaces;
1449f074e33SMatthew G Knepley     PetscInt        numCellFaces, faceSize, cf;
1459f074e33SMatthew G Knepley 
1469f074e33SMatthew G Knepley     ierr = DMPlexGetFaces(dm, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
1479f074e33SMatthew G Knepley     if (faceSize != 2) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Triangles cannot have face of size %D", faceSize);
1489f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
1499f074e33SMatthew G Knepley #if 1
1509f074e33SMatthew G Knepley       PetscHashIJKey key;
1519f074e33SMatthew G Knepley 
1529f074e33SMatthew G Knepley       key.i = PetscMin(cellFaces[cf*faceSize+0], cellFaces[cf*faceSize+1]);
1539f074e33SMatthew G Knepley       key.j = PetscMax(cellFaces[cf*faceSize+0], cellFaces[cf*faceSize+1]);
1549f074e33SMatthew G Knepley       ierr  = PetscHashIJGet(edgeTable, key, &e);CHKERRQ(ierr);
1559f074e33SMatthew G Knepley       if (e < 0) {
1569f074e33SMatthew G Knepley         ierr = DMPlexSetCone(idm, edge, &cellFaces[cf*faceSize]);CHKERRQ(ierr);
1579f074e33SMatthew G Knepley         ierr = PetscHashIJAdd(edgeTable, key, edge);CHKERRQ(ierr);
1589f074e33SMatthew G Knepley         e    = edge++;
1599f074e33SMatthew G Knepley       }
1609f074e33SMatthew G Knepley #else
1619f074e33SMatthew G Knepley       PetscBool found = PETSC_FALSE;
1629f074e33SMatthew G Knepley 
1639f074e33SMatthew G Knepley       /* TODO Need join of vertices to check for existence of edges, which needs support (could set edge support), so just brute force for now */
1649f074e33SMatthew G Knepley       for (e = firstEdge; e < edge; ++e) {
1659f074e33SMatthew G Knepley         const PetscInt *cone;
1669f074e33SMatthew G Knepley 
1679f074e33SMatthew G Knepley         ierr = DMPlexGetCone(idm, e, &cone);CHKERRQ(ierr);
1689f074e33SMatthew G Knepley         if (((cellFaces[cf*faceSize+0] == cone[0]) && (cellFaces[cf*faceSize+1] == cone[1])) ||
1699f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[1]) && (cellFaces[cf*faceSize+1] == cone[0]))) {
1709f074e33SMatthew G Knepley           found = PETSC_TRUE;
1719f074e33SMatthew G Knepley           break;
1729f074e33SMatthew G Knepley         }
1739f074e33SMatthew G Knepley       }
1749f074e33SMatthew G Knepley       if (!found) {
1759f074e33SMatthew G Knepley         ierr = DMPlexSetCone(idm, edge, &cellFaces[cf*faceSize]);CHKERRQ(ierr);
1769f074e33SMatthew G Knepley         ++edge;
1779f074e33SMatthew G Knepley       }
1789f074e33SMatthew G Knepley #endif
1799f074e33SMatthew G Knepley       ierr = DMPlexInsertCone(idm, c, cf, e);CHKERRQ(ierr);
1809f074e33SMatthew G Knepley     }
1819f074e33SMatthew G Knepley   }
1829f074e33SMatthew G Knepley   if (edge != firstEdge+numEdges) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid number of edges %D should be %D", edge-firstEdge, numEdges);
1839f074e33SMatthew G Knepley   ierr = PetscHashIJDestroy(&edgeTable);CHKERRQ(ierr);
1849f074e33SMatthew G Knepley   ierr = PetscFree(off);CHKERRQ(ierr);
1859f074e33SMatthew G Knepley   ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr);
1869f074e33SMatthew G Knepley   ierr = DMPlexStratify(idm);CHKERRQ(ierr);
1879f074e33SMatthew G Knepley   mesh = (DM_Plex*) (idm)->data;
1889f074e33SMatthew G Knepley   /* Orient edges */
1899f074e33SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
1909f074e33SMatthew G Knepley     const PetscInt *cone = NULL, *cellFaces;
1919f074e33SMatthew G Knepley     PetscInt        coneSize, coff, numCellFaces, faceSize, cf;
1929f074e33SMatthew G Knepley 
1939f074e33SMatthew G Knepley     ierr = DMPlexGetConeSize(idm, c, &coneSize);CHKERRQ(ierr);
1949f074e33SMatthew G Knepley     ierr = DMPlexGetCone(idm, c, &cone);CHKERRQ(ierr);
1959f074e33SMatthew G Knepley     ierr = PetscSectionGetOffset(mesh->coneSection, c, &coff);CHKERRQ(ierr);
1969f074e33SMatthew G Knepley     ierr = DMPlexGetFaces(dm, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
1979f074e33SMatthew G Knepley     if (coneSize != numCellFaces) SETERRQ3(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid number of edges %D for cell %D should be %D", coneSize, c, numCellFaces);
1989f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
1999f074e33SMatthew G Knepley       const PetscInt *econe = NULL;
2009f074e33SMatthew G Knepley       PetscInt        esize;
2019f074e33SMatthew G Knepley 
2029f074e33SMatthew G Knepley       ierr = DMPlexGetConeSize(idm, cone[cf], &esize);CHKERRQ(ierr);
2039f074e33SMatthew G Knepley       ierr = DMPlexGetCone(idm, cone[cf], &econe);CHKERRQ(ierr);
2049f074e33SMatthew G Knepley       if (esize != 2) SETERRQ2(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid number of edge endpoints %D for edge %D should be 2", esize, cone[cf]);
2059f074e33SMatthew G Knepley       if ((cellFaces[cf*faceSize+0] == econe[0]) && (cellFaces[cf*faceSize+1] == econe[1])) {
2069f074e33SMatthew G Knepley         /* Correctly oriented */
2079f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = 0;
2089f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == econe[1]) && (cellFaces[cf*faceSize+1] == econe[0])) {
2099f074e33SMatthew G Knepley         /* Start at index 1, and reverse orientation */
2109f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = -(1+1);
2119f074e33SMatthew G Knepley       }
2129f074e33SMatthew G Knepley     }
2139f074e33SMatthew G Knepley   }
2149f074e33SMatthew G Knepley   *dmInt = idm;
2159f074e33SMatthew G Knepley   PetscFunctionReturn(0);
2169f074e33SMatthew G Knepley }
2179f074e33SMatthew G Knepley 
2189f074e33SMatthew G Knepley #undef __FUNCT__
2199f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolate_3D"
2209f074e33SMatthew G Knepley static PetscErrorCode DMPlexInterpolate_3D(DM dm, DM *dmInt)
2219f074e33SMatthew G Knepley {
2229f074e33SMatthew G Knepley   DM             idm, fdm;
2239f074e33SMatthew G Knepley   DM_Plex       *mesh;
2249f074e33SMatthew G Knepley   PetscInt      *off;
2259f074e33SMatthew G Knepley   const PetscInt numCorners = 4;
2269f074e33SMatthew G Knepley   PetscInt       dim, numCells, cStart, cEnd, c, numVertices, vStart, vEnd;
2279f074e33SMatthew G Knepley   PetscInt       numFaces, firstFace, face, f, numEdges, firstEdge, edge, e;
2289f074e33SMatthew G Knepley   PetscErrorCode ierr;
2299f074e33SMatthew G Knepley 
2309f074e33SMatthew G Knepley   PetscFunctionBegin;
2319f074e33SMatthew G Knepley   {
2329f074e33SMatthew G Knepley     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr);
2339f074e33SMatthew G Knepley     ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
2349f074e33SMatthew G Knepley     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
2359f074e33SMatthew G Knepley   }
2369f074e33SMatthew G Knepley   ierr        = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2379f074e33SMatthew G Knepley   ierr        = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2389f074e33SMatthew G Knepley   ierr        = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
2399f074e33SMatthew G Knepley   numCells    = cEnd - cStart;
2409f074e33SMatthew G Knepley   numVertices = vEnd - vStart;
2419f074e33SMatthew G Knepley   firstFace   = numCells + numVertices;
2429f074e33SMatthew G Knepley   numFaces    = 0;
2439f074e33SMatthew G Knepley   /* Count faces using algorithm from CreateNeighborCSR */
2449f074e33SMatthew G Knepley   ierr = DMPlexCreateNeighborCSR(dm, NULL, &off, NULL);CHKERRQ(ierr);
2459f074e33SMatthew G Knepley   if (off) {
2469f074e33SMatthew G Knepley     numFaces = off[numCells]/2;
2479f074e33SMatthew G Knepley     /* Account for boundary faces: \sum_c 4 - neighbors = 4*numCells - totalNeighbors */
2489f074e33SMatthew G Knepley     numFaces += 4*numCells - off[numCells];
2499f074e33SMatthew G Knepley   }
2509f074e33SMatthew G Knepley   /* Use Euler characteristic to get edges V - E + F - C = 1 */
2519f074e33SMatthew G Knepley   firstEdge = firstFace + numFaces;
2529f074e33SMatthew G Knepley   numEdges  = numVertices + numFaces - numCells - 1;
2539f074e33SMatthew G Knepley   /* Create interpolated mesh */
2549f074e33SMatthew G Knepley   ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);
2559f074e33SMatthew G Knepley   ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr);
2569f074e33SMatthew G Knepley   ierr = DMPlexSetDimension(idm, dim);CHKERRQ(ierr);
2579f074e33SMatthew G Knepley   ierr = DMPlexSetChart(idm, 0, numCells+numVertices+numFaces+numEdges);CHKERRQ(ierr);
2589f074e33SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
2599f074e33SMatthew G Knepley     ierr = DMPlexSetConeSize(idm, c, numCorners);CHKERRQ(ierr);
2609f074e33SMatthew G Knepley   }
2619f074e33SMatthew G Knepley   for (f = firstFace; f < firstFace+numFaces; ++f) {
2629f074e33SMatthew G Knepley     ierr = DMPlexSetConeSize(idm, f, 3);CHKERRQ(ierr);
2639f074e33SMatthew G Knepley   }
2649f074e33SMatthew G Knepley   for (e = firstEdge; e < firstEdge+numEdges; ++e) {
2659f074e33SMatthew G Knepley     ierr = DMPlexSetConeSize(idm, e, 2);CHKERRQ(ierr);
2669f074e33SMatthew G Knepley   }
2679f074e33SMatthew G Knepley   ierr = DMSetUp(idm);CHKERRQ(ierr);
2689f074e33SMatthew G Knepley   /* Get face cones from subsets of cell vertices */
2699f074e33SMatthew G Knepley   ierr = DMCreate(PetscObjectComm((PetscObject)dm), &fdm);CHKERRQ(ierr);
2709f074e33SMatthew G Knepley   ierr = DMSetType(fdm, DMPLEX);CHKERRQ(ierr);
2719f074e33SMatthew G Knepley   ierr = DMPlexSetDimension(fdm, dim);CHKERRQ(ierr);
2729f074e33SMatthew G Knepley   ierr = DMPlexSetChart(fdm, numCells, firstFace+numFaces);CHKERRQ(ierr);
2739f074e33SMatthew G Knepley   for (f = firstFace; f < firstFace+numFaces; ++f) {
2749f074e33SMatthew G Knepley     ierr = DMPlexSetConeSize(fdm, f, 3);CHKERRQ(ierr);
2759f074e33SMatthew G Knepley   }
2769f074e33SMatthew G Knepley   ierr = DMSetUp(fdm);CHKERRQ(ierr);
2779f074e33SMatthew G Knepley   for (c = 0, face = firstFace; c < numCells; ++c) {
2789f074e33SMatthew G Knepley     const PetscInt *cellFaces;
2799f074e33SMatthew G Knepley     PetscInt        numCellFaces, faceSize, cf;
2809f074e33SMatthew G Knepley 
2819f074e33SMatthew G Knepley     ierr = DMPlexGetFaces(dm, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
2829f074e33SMatthew G Knepley     if (faceSize != 3) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Tetrahedra cannot have face of size %D", faceSize);
2839f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
2849f074e33SMatthew G Knepley       PetscBool found = PETSC_FALSE;
2859f074e33SMatthew G Knepley 
2869f074e33SMatthew G Knepley       /* TODO Need join of vertices to check for existence of edges, which needs support (could set edge support), so just brute force for now */
2879f074e33SMatthew G Knepley       for (f = firstFace; f < face; ++f) {
2889f074e33SMatthew G Knepley         const PetscInt *cone = NULL;
2899f074e33SMatthew G Knepley 
2909f074e33SMatthew G Knepley         ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
2919f074e33SMatthew G Knepley         if (((cellFaces[cf*faceSize+0] == cone[0]) && (cellFaces[cf*faceSize+1] == cone[1]) && (cellFaces[cf*faceSize+2] == cone[2])) ||
2929f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[1]) && (cellFaces[cf*faceSize+1] == cone[2]) && (cellFaces[cf*faceSize+2] == cone[0])) ||
2939f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[2]) && (cellFaces[cf*faceSize+1] == cone[0]) && (cellFaces[cf*faceSize+2] == cone[1])) ||
2949f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[0]) && (cellFaces[cf*faceSize+1] == cone[2]) && (cellFaces[cf*faceSize+2] == cone[1])) ||
2959f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[2]) && (cellFaces[cf*faceSize+1] == cone[1]) && (cellFaces[cf*faceSize+2] == cone[0])) ||
2969f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[1]) && (cellFaces[cf*faceSize+1] == cone[0]) && (cellFaces[cf*faceSize+2] == cone[2]))) {
2979f074e33SMatthew G Knepley           found = PETSC_TRUE;
2989f074e33SMatthew G Knepley           break;
2999f074e33SMatthew G Knepley         }
3009f074e33SMatthew G Knepley       }
3019f074e33SMatthew G Knepley       if (!found) {
3029f074e33SMatthew G Knepley         ierr = DMPlexSetCone(idm, face, &cellFaces[cf*faceSize]);CHKERRQ(ierr);
3039f074e33SMatthew G Knepley         /* Save the vertices for orientation calculation */
3049f074e33SMatthew G Knepley         ierr = DMPlexSetCone(fdm, face, &cellFaces[cf*faceSize]);CHKERRQ(ierr);
3059f074e33SMatthew G Knepley         ++face;
3069f074e33SMatthew G Knepley       }
3079f074e33SMatthew G Knepley       ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr);
3089f074e33SMatthew G Knepley     }
3099f074e33SMatthew G Knepley   }
3109f074e33SMatthew G Knepley   if (face != firstFace+numFaces) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid number of faces %D should be %D", face-firstFace, numFaces);
3119f074e33SMatthew G Knepley   /* Get edge cones from subsets of face vertices */
3129f074e33SMatthew G Knepley   for (f = firstFace, edge = firstEdge; f < firstFace+numFaces; ++f) {
3139f074e33SMatthew G Knepley     const PetscInt *cellFaces;
3149f074e33SMatthew G Knepley     PetscInt        numCellFaces, faceSize, cf;
3159f074e33SMatthew G Knepley 
3169f074e33SMatthew G Knepley     ierr = DMPlexGetFaces(idm, f, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
3179f074e33SMatthew G Knepley     if (faceSize != 2) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Triangles cannot have face of size %D", faceSize);
3189f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
3199f074e33SMatthew G Knepley       PetscBool found = PETSC_FALSE;
3209f074e33SMatthew G Knepley 
3219f074e33SMatthew G Knepley       /* TODO Need join of vertices to check for existence of edges, which needs support (could set edge support), so just brute force for now */
3229f074e33SMatthew G Knepley       for (e = firstEdge; e < edge; ++e) {
3239f074e33SMatthew G Knepley         const PetscInt *cone = NULL;
3249f074e33SMatthew G Knepley 
3259f074e33SMatthew G Knepley         ierr = DMPlexGetCone(idm, e, &cone);CHKERRQ(ierr);
3269f074e33SMatthew G Knepley         if (((cellFaces[cf*faceSize+0] == cone[0]) && (cellFaces[cf*faceSize+1] == cone[1])) ||
3279f074e33SMatthew G Knepley             ((cellFaces[cf*faceSize+0] == cone[1]) && (cellFaces[cf*faceSize+1] == cone[0]))) {
3289f074e33SMatthew G Knepley           found = PETSC_TRUE;
3299f074e33SMatthew G Knepley           break;
3309f074e33SMatthew G Knepley         }
3319f074e33SMatthew G Knepley       }
3329f074e33SMatthew G Knepley       if (!found) {
3339f074e33SMatthew G Knepley         ierr = DMPlexSetCone(idm, edge, &cellFaces[cf*faceSize]);CHKERRQ(ierr);
3349f074e33SMatthew G Knepley         ++edge;
3359f074e33SMatthew G Knepley       }
3369f074e33SMatthew G Knepley       ierr = DMPlexInsertCone(idm, f, cf, e);CHKERRQ(ierr);
3379f074e33SMatthew G Knepley     }
3389f074e33SMatthew G Knepley   }
3399f074e33SMatthew G Knepley   if (edge != firstEdge+numEdges) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid number of edges %D should be %D", edge-firstEdge, numEdges);
3409f074e33SMatthew G Knepley   ierr = PetscFree(off);CHKERRQ(ierr);
3419f074e33SMatthew G Knepley   ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr);
3429f074e33SMatthew G Knepley   ierr = DMPlexStratify(idm);CHKERRQ(ierr);
3439f074e33SMatthew G Knepley   mesh = (DM_Plex*) (idm)->data;
3449f074e33SMatthew G Knepley   /* Orient edges */
3459f074e33SMatthew G Knepley   for (f = firstFace; f < firstFace+numFaces; ++f) {
3469f074e33SMatthew G Knepley     const PetscInt *cone, *cellFaces;
3479f074e33SMatthew G Knepley     PetscInt        coneSize, coff, numCellFaces, faceSize, cf;
3489f074e33SMatthew G Knepley 
3499f074e33SMatthew G Knepley     ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr);
3509f074e33SMatthew G Knepley     ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr);
3519f074e33SMatthew G Knepley     ierr = PetscSectionGetOffset(mesh->coneSection, f, &coff);CHKERRQ(ierr);
3529f074e33SMatthew G Knepley     ierr = DMPlexGetFaces(fdm, f, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
3539f074e33SMatthew G Knepley     if (coneSize != numCellFaces) SETERRQ3(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid number of edges %D for face %D should be %D", coneSize, f, numCellFaces);
3549f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
3559f074e33SMatthew G Knepley       const PetscInt *econe;
3569f074e33SMatthew G Knepley       PetscInt        esize;
3579f074e33SMatthew G Knepley 
3589f074e33SMatthew G Knepley       ierr = DMPlexGetConeSize(idm, cone[cf], &esize);CHKERRQ(ierr);
3599f074e33SMatthew G Knepley       ierr = DMPlexGetCone(idm, cone[cf], &econe);CHKERRQ(ierr);
3609f074e33SMatthew G Knepley       if (esize != 2) SETERRQ2(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid number of edge endpoints %D for edge %D should be 2", esize, cone[cf]);
3619f074e33SMatthew G Knepley       if ((cellFaces[cf*faceSize+0] == econe[0]) && (cellFaces[cf*faceSize+1] == econe[1])) {
3629f074e33SMatthew G Knepley         /* Correctly oriented */
3639f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = 0;
3649f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == econe[1]) && (cellFaces[cf*faceSize+1] == econe[0])) {
3659f074e33SMatthew G Knepley         /* Start at index 1, and reverse orientation */
3669f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = -(1+1);
3679f074e33SMatthew G Knepley       }
3689f074e33SMatthew G Knepley     }
3699f074e33SMatthew G Knepley   }
3709f074e33SMatthew G Knepley   ierr = DMDestroy(&fdm);CHKERRQ(ierr);
3719f074e33SMatthew G Knepley   /* Orient faces */
3729f074e33SMatthew G Knepley   for (c = 0; c < numCells; ++c) {
3739f074e33SMatthew G Knepley     const PetscInt *cone, *cellFaces;
3749f074e33SMatthew G Knepley     PetscInt        coneSize, coff, numCellFaces, faceSize, cf;
3759f074e33SMatthew G Knepley 
3769f074e33SMatthew G Knepley     ierr = DMPlexGetConeSize(idm, c, &coneSize);CHKERRQ(ierr);
3779f074e33SMatthew G Knepley     ierr = DMPlexGetCone(idm, c, &cone);CHKERRQ(ierr);
3789f074e33SMatthew G Knepley     ierr = PetscSectionGetOffset(mesh->coneSection, c, &coff);CHKERRQ(ierr);
3799f074e33SMatthew G Knepley     ierr = DMPlexGetFaces(dm, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr);
3809f074e33SMatthew G Knepley     if (coneSize != numCellFaces) SETERRQ3(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid number of edges %D for cell %D should be %D", coneSize, c, numCellFaces);
3819f074e33SMatthew G Knepley     for (cf = 0; cf < numCellFaces; ++cf) {
3829f074e33SMatthew G Knepley       PetscInt *origClosure = NULL, *closure;
3839f074e33SMatthew G Knepley       PetscInt closureSize, i;
3849f074e33SMatthew G Knepley 
3859f074e33SMatthew G Knepley       ierr = DMPlexGetTransitiveClosure(idm, cone[cf], PETSC_TRUE, &closureSize, &origClosure);CHKERRQ(ierr);
3869f074e33SMatthew G Knepley       if (closureSize != 7) SETERRQ2(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid closure size %D for face %D should be 7", closureSize, cone[cf]);
3879f074e33SMatthew G Knepley       for (i = 4; i < 7; ++i) {
3889f074e33SMatthew G Knepley         if ((origClosure[i*2] < vStart) || (origClosure[i*2] >= vEnd)) SETERRQ3(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Invalid closure point %D should be a vertex in [%D, %D)", origClosure[i*2], vStart, vEnd);
3899f074e33SMatthew G Knepley       }
3909f074e33SMatthew G Knepley       closure = &origClosure[4*2];
3919f074e33SMatthew G Knepley       /* Remember that this is the orientation for edges, not vertices */
3929f074e33SMatthew G Knepley       if        ((cellFaces[cf*faceSize+0] == closure[0*2]) && (cellFaces[cf*faceSize+1] == closure[1*2]) && (cellFaces[cf*faceSize+2] == closure[2*2])) {
3939f074e33SMatthew G Knepley         /* Correctly oriented */
3949f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = 0;
3959f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == closure[1*2]) && (cellFaces[cf*faceSize+1] == closure[2*2]) && (cellFaces[cf*faceSize+2] == closure[0*2])) {
3969f074e33SMatthew G Knepley         /* Shifted by 1 */
3979f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = 1;
3989f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == closure[2*2]) && (cellFaces[cf*faceSize+1] == closure[0*2]) && (cellFaces[cf*faceSize+2] == closure[1*2])) {
3999f074e33SMatthew G Knepley         /* Shifted by 2 */
4009f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = 2;
4019f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == closure[2*2]) && (cellFaces[cf*faceSize+1] == closure[1*2]) && (cellFaces[cf*faceSize+2] == closure[0*2])) {
4029f074e33SMatthew G Knepley         /* Start at edge 1, and reverse orientation */
4039f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = -(1+1);
4049f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == closure[1*2]) && (cellFaces[cf*faceSize+1] == closure[0*2]) && (cellFaces[cf*faceSize+2] == closure[2*2])) {
4059f074e33SMatthew G Knepley         /* Start at index 0, and reverse orientation */
4069f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = -(0+1);
4079f074e33SMatthew G Knepley       } else if ((cellFaces[cf*faceSize+0] == closure[0*2]) && (cellFaces[cf*faceSize+1] == closure[2*2]) && (cellFaces[cf*faceSize+2] == closure[1*2])) {
4089f074e33SMatthew G Knepley         /* Start at index 2, and reverse orientation */
4099f074e33SMatthew G Knepley         mesh->coneOrientations[coff+cf] = -(2+1);
4109f074e33SMatthew G Knepley       } else SETERRQ3(PetscObjectComm((PetscObject)idm), PETSC_ERR_PLIB, "Face %D did not match local face %D in cell %D for any orientation", cone[cf], cf, c);
4119f074e33SMatthew G Knepley       ierr = DMPlexRestoreTransitiveClosure(idm, cone[cf], PETSC_TRUE, &closureSize, &origClosure);CHKERRQ(ierr);
4129f074e33SMatthew G Knepley     }
4139f074e33SMatthew G Knepley   }
4149f074e33SMatthew G Knepley   {
4159f074e33SMatthew G Knepley     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr);
4169f074e33SMatthew G Knepley     ierr = DMView(idm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
4179f074e33SMatthew G Knepley     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
4189f074e33SMatthew G Knepley   }
4199f074e33SMatthew G Knepley   *dmInt = idm;
4209f074e33SMatthew G Knepley   PetscFunctionReturn(0);
4219f074e33SMatthew G Knepley }
4229f074e33SMatthew G Knepley 
4239f074e33SMatthew G Knepley #undef __FUNCT__
4249f074e33SMatthew G Knepley #define __FUNCT__ "DMPlexInterpolate"
4259f074e33SMatthew G Knepley PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
4269f074e33SMatthew G Knepley {
4279f074e33SMatthew G Knepley   PetscInt       dim;
4289f074e33SMatthew G Knepley   PetscErrorCode ierr;
4299f074e33SMatthew G Knepley 
4309f074e33SMatthew G Knepley   PetscFunctionBegin;
4319f074e33SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
4329f074e33SMatthew G Knepley   switch (dim) {
4339f074e33SMatthew G Knepley   case 2:
4349f074e33SMatthew G Knepley     ierr = DMPlexInterpolate_2D(dm, dmInt);CHKERRQ(ierr);break;
4359f074e33SMatthew G Knepley   case 3:
4369f074e33SMatthew G Knepley     ierr = DMPlexInterpolate_3D(dm, dmInt);CHKERRQ(ierr);break;
4379f074e33SMatthew G Knepley   default:
4389f074e33SMatthew G Knepley     SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No mesh interpolation support for dimension %D", dim);
4399f074e33SMatthew G Knepley   }
4409f074e33SMatthew G Knepley   PetscFunctionReturn(0);
4419f074e33SMatthew G Knepley }
442*07b0a1fcSMatthew G Knepley 
443*07b0a1fcSMatthew G Knepley #undef __FUNCT__
444*07b0a1fcSMatthew G Knepley #define __FUNCT__ "DMPlexCopyCoordinates"
445*07b0a1fcSMatthew G Knepley PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB)
446*07b0a1fcSMatthew G Knepley {
447*07b0a1fcSMatthew G Knepley   Vec            coordinatesA, coordinatesB;
448*07b0a1fcSMatthew G Knepley   PetscSection   coordSectionA, coordSectionB;
449*07b0a1fcSMatthew G Knepley   PetscScalar   *coordsA, *coordsB;
450*07b0a1fcSMatthew G Knepley   PetscInt       spaceDim, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d;
451*07b0a1fcSMatthew G Knepley   PetscErrorCode ierr;
452*07b0a1fcSMatthew G Knepley 
453*07b0a1fcSMatthew G Knepley   PetscFunctionBegin;
454*07b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr);
455*07b0a1fcSMatthew G Knepley   ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr);
456*07b0a1fcSMatthew G Knepley   if ((vEndA-vStartA) != (vEndB-vStartB)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of vertices in first DM %d != %d in the second DM", vEndA-vStartA, vEndB-vStartB);
457*07b0a1fcSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr);
458*07b0a1fcSMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr);
459*07b0a1fcSMatthew G Knepley   ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr);
460*07b0a1fcSMatthew G Knepley   ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr);
461*07b0a1fcSMatthew G Knepley   ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr);
462*07b0a1fcSMatthew G Knepley   ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr);
463*07b0a1fcSMatthew G Knepley   for (v = vStartB; v < vEndB; ++v) {
464*07b0a1fcSMatthew G Knepley     ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr);
465*07b0a1fcSMatthew G Knepley     ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr);
466*07b0a1fcSMatthew G Knepley   }
467*07b0a1fcSMatthew G Knepley   ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr);
468*07b0a1fcSMatthew G Knepley   ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr);
469*07b0a1fcSMatthew G Knepley   ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr);
470*07b0a1fcSMatthew G Knepley   ierr = VecCreate(PetscObjectComm((PetscObject) dmB), &coordinatesB);CHKERRQ(ierr);
471*07b0a1fcSMatthew G Knepley   ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr);
472*07b0a1fcSMatthew G Knepley   ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr);
473*07b0a1fcSMatthew G Knepley   ierr = VecSetFromOptions(coordinatesB);CHKERRQ(ierr);
474*07b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr);
475*07b0a1fcSMatthew G Knepley   ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr);
476*07b0a1fcSMatthew G Knepley   for (v = 0; v < vEndB-vStartB; ++v) {
477*07b0a1fcSMatthew G Knepley     for (d = 0; d < spaceDim; ++d) {
478*07b0a1fcSMatthew G Knepley       coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d];
479*07b0a1fcSMatthew G Knepley     }
480*07b0a1fcSMatthew G Knepley   }
481*07b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr);
482*07b0a1fcSMatthew G Knepley   ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr);
483*07b0a1fcSMatthew G Knepley   ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr);
484*07b0a1fcSMatthew G Knepley   ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr);
485*07b0a1fcSMatthew G Knepley   PetscFunctionReturn(0);
486*07b0a1fcSMatthew G Knepley }
487