<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Mon, May 5, 2014 at 9:35 AM, Lawrence Mitchell <span dir="ltr"><<a href="mailto:lawrence.mitchell@imperial.ac.uk" target="_blank">lawrence.mitchell@imperial.ac.uk</a>></span> wrote:<br>
<blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">If the plex is an embedded manifold, its topological dimension will not<br>
match the dimension of the coordinate section.  Refine coordinates by<br>
inspecting the latter (rather than the former) when carrying out uniform<br>
refinement.<br>
---<br>
 I sent this a while ago, but I'm not sure if it got lost in transit<br>
 (hence reposting to petsc-dev for comment).  It's possible that<br>
 something similar needs to be done for non-uniform refinement, but<br>
 I'm not sure how to drive those interfaces.  We need this, for<br>
 example, when refining a mesh of the surface of the earth.<br></blockquote><div><br></div><div>Yes, I lost it. WIll commit to next tonight.</div><div><br></div><div>  Thanks,</div><div><br></div><div>    Matt</div><div>
 </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
Cheers,<br>
Lawrence<br>
<br>
 src/dm/impls/plex/plexrefine.c | 32 ++++++++++++++++----------------<br>
 1 file changed, 16 insertions(+), 16 deletions(-)<br>
<br>
diff --git a/src/dm/impls/plex/plexrefine.c b/src/dm/impls/plex/plexrefine.c<br>
index fe4d71d..d6cbc74 100644<br>
--- a/src/dm/impls/plex/plexrefine.c<br>
+++ b/src/dm/impls/plex/plexrefine.c<br>
@@ -5339,11 +5339,10 @@ static PetscErrorCode CellRefinerSetCoordinates(CellRefiner refiner, DM dm, Pets<br>
   Vec            coordinates, coordinatesNew;<br>
   PetscScalar   *coords, *coordsNew;<br>
   const PetscInt numVertices = depthSize ? depthSize[0] : 0;<br>
-  PetscInt       dim, depth, bs, coordSizeNew, cStart, cEnd, cMax, c, vStart, vStartNew, vEnd, v, eStart, eEnd, eMax, e, fStart, fEnd, fMax, f;<br>
+  PetscInt       spaceDim, depth, bs, coordSizeNew, cStart, cEnd, cMax, c, vStart, vStartNew, vEnd, v, eStart, eEnd, eMax, e, fStart, fEnd, fMax, f;<br>
   PetscErrorCode ierr;<br>
<br>
   PetscFunctionBegin;<br>
-  ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);<br>
   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);<br>
   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);<br>
   ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);<br>
@@ -5353,17 +5352,18 @@ static PetscErrorCode CellRefinerSetCoordinates(CellRefiner refiner, DM dm, Pets<br>
   if (refiner) {ierr = GetDepthStart_Private(depth, depthSize, NULL, NULL, NULL, &vStartNew);CHKERRQ(ierr);}<br>
   ierr = GetDepthStart_Private(depth, depthSize, NULL, NULL, NULL, &vStartNew);CHKERRQ(ierr);<br>
   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);<br>
+  ierr = PetscSectionGetFieldComponents(coordSection, 0, &spaceDim);CHKERRQ(ierr);<br>
   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &coordSectionNew);CHKERRQ(ierr);<br>
   ierr = PetscSectionSetNumFields(coordSectionNew, 1);CHKERRQ(ierr);<br>
-  ierr = PetscSectionSetFieldComponents(coordSectionNew, 0, dim);CHKERRQ(ierr);<br>
+  ierr = PetscSectionSetFieldComponents(coordSectionNew, 0, spaceDim);CHKERRQ(ierr);<br>
   ierr = PetscSectionSetChart(coordSectionNew, vStartNew, vStartNew+numVertices);CHKERRQ(ierr);<br>
   if (cMax < 0) cMax = cEnd;<br>
   if (fMax < 0) fMax = fEnd;<br>
   if (eMax < 0) eMax = eEnd;<br>
-  /* All vertices have the dim coordinates */<br>
+  /* All vertices have spaceDim coordinates */<br>
   for (v = vStartNew; v < vStartNew+numVertices; ++v) {<br>
-    ierr = PetscSectionSetDof(coordSectionNew, v, dim);CHKERRQ(ierr);<br>
-    ierr = PetscSectionSetFieldDof(coordSectionNew, v, 0, dim);CHKERRQ(ierr);<br>
+    ierr = PetscSectionSetDof(coordSectionNew, v, spaceDim);CHKERRQ(ierr);<br>
+    ierr = PetscSectionSetFieldDof(coordSectionNew, v, 0, spaceDim);CHKERRQ(ierr);<br>
   }<br>
   ierr = PetscSectionSetUp(coordSectionNew);CHKERRQ(ierr);<br>
   ierr = DMSetCoordinateSection(rdm, coordSectionNew);CHKERRQ(ierr);<br>
@@ -5396,16 +5396,16 @@ static PetscErrorCode CellRefinerSetCoordinates(CellRefiner refiner, DM dm, Pets<br>
         ierr = PetscSectionGetOffset(coordSection, cone[v], &off[v]);CHKERRQ(ierr);<br>
       }<br>
       ierr = PetscSectionGetOffset(coordSectionNew, newv, &offnew);CHKERRQ(ierr);<br>
-      for (d = 0; d < dim; ++d) coordsNew[offnew+d] = 0.0;<br>
-      for (v = 0; v < coneSize; ++v) {ierr = DMPlexLocalizeAddCoordinate_Internal(dm, dim, &coords[off[0]], &coords[off[v]], &coordsNew[offnew]);CHKERRQ(ierr);}<br>
-      for (d = 0; d < dim; ++d) coordsNew[offnew+d] /= coneSize;<br>
+      for (d = 0; d < spaceDim; ++d) coordsNew[offnew+d] = 0.0;<br>
+      for (v = 0; v < coneSize; ++v) {ierr = DMPlexLocalizeAddCoordinate_Internal(dm, spaceDim, &coords[off[0]], &coords[off[v]], &coordsNew[offnew]);CHKERRQ(ierr);}<br>
+      for (d = 0; d < spaceDim; ++d) coordsNew[offnew+d] /= coneSize;<br>
       ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &cone);CHKERRQ(ierr);<br>
     }<br>
   case 2: /* Hex 2D */<br>
   case 4: /* Hybrid Hex 2D */<br>
     /* Cell vertices have the average of corner coordinates */<br>
     for (c = cStart; c < cMax; ++c) {<br>
-      const PetscInt newv = vStartNew + (vEnd - vStart) + (eMax - eStart) + (c - cStart) + (dim > 2 ? (fMax - fStart) : 0);<br>
+      const PetscInt newv = vStartNew + (vEnd - vStart) + (eMax - eStart) + (c - cStart) + (spaceDim > 2 ? (fMax - fStart) : 0);<br>
       PetscInt      *cone = NULL;<br>
       PetscInt       closureSize, coneSize = 0, off[8], offnew, p, d;<br>
<br>
@@ -5418,9 +5418,9 @@ static PetscErrorCode CellRefinerSetCoordinates(CellRefiner refiner, DM dm, Pets<br>
         ierr = PetscSectionGetOffset(coordSection, cone[v], &off[v]);CHKERRQ(ierr);<br>
       }<br>
       ierr = PetscSectionGetOffset(coordSectionNew, newv, &offnew);CHKERRQ(ierr);<br>
-      for (d = 0; d < dim; ++d) coordsNew[offnew+d] = 0.0;<br>
-      for (v = 0; v < coneSize; ++v) {ierr = DMPlexLocalizeAddCoordinate_Internal(dm, dim, &coords[off[0]], &coords[off[v]], &coordsNew[offnew]);CHKERRQ(ierr);}<br>
-      for (d = 0; d < dim; ++d) coordsNew[offnew+d] /= coneSize;<br>
+      for (d = 0; d < spaceDim; ++d) coordsNew[offnew+d] = 0.0;<br>
+      for (v = 0; v < coneSize; ++v) {ierr = DMPlexLocalizeAddCoordinate_Internal(dm, spaceDim, &coords[off[0]], &coords[off[v]], &coordsNew[offnew]);CHKERRQ(ierr);}<br>
+      for (d = 0; d < spaceDim; ++d) coordsNew[offnew+d] /= coneSize;<br>
       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &cone);CHKERRQ(ierr);<br>
     }<br>
   case 1: /* Simplicial 2D */<br>
@@ -5439,8 +5439,8 @@ static PetscErrorCode CellRefinerSetCoordinates(CellRefiner refiner, DM dm, Pets<br>
       ierr = PetscSectionGetOffset(coordSection, cone[0], &offA);CHKERRQ(ierr);<br>
       ierr = PetscSectionGetOffset(coordSection, cone[1], &offB);CHKERRQ(ierr);<br>
       ierr = PetscSectionGetOffset(coordSectionNew, newv, &offnew);CHKERRQ(ierr);<br>
-      ierr = DMPlexLocalizeCoordinate_Internal(dm, dim, &coords[offA], &coords[offB], &coordsNew[offnew]);CHKERRQ(ierr);<br>
-      for (d = 0; d < dim; ++d) {<br>
+      ierr = DMPlexLocalizeCoordinate_Internal(dm, spaceDim, &coords[offA], &coords[offB], &coordsNew[offnew]);CHKERRQ(ierr);<br>
+      for (d = 0; d < spaceDim; ++d) {<br>
         coordsNew[offnew+d] = 0.5*(coords[offA+d] + coordsNew[offnew+d]);<br>
       }<br>
     }<br>
@@ -5451,7 +5451,7 @@ static PetscErrorCode CellRefinerSetCoordinates(CellRefiner refiner, DM dm, Pets<br>
<br>
       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);<br>
       ierr = PetscSectionGetOffset(coordSectionNew, newv, &offnew);CHKERRQ(ierr);<br>
-      for (d = 0; d < dim; ++d) {<br>
+      for (d = 0; d < spaceDim; ++d) {<br>
         coordsNew[offnew+d] = coords[off+d];<br>
       }<br>
     }<br>
<span class="HOEnZb"><font color="#888888">--<br>
1.8.4.474.g128a96c<br>
<br>
</font></span></blockquote></div><br><br clear="all"><div><br></div>-- <br>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>
-- Norbert Wiener
</div></div>