<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Mon, Apr 20, 2015 at 6:09 AM, Justin Chang <span dir="ltr"><<a href="mailto:jchang27@uh.edu" target="_blank">jchang27@uh.edu</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div><div><div>Matt, <br><br></div>I have attempted to implement this, called UninterpolateSF(), and attached is the code for the function. However, I suspect that there are bugs because when I refine the mesh via -dm_refine my program crashes.<br></div></div></div></blockquote><div><br></div><div>I have finally fixed this. It is in next.</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"><div dir="ltr"><div><div></div><div>In my program I have the following lines for uninterpolating a distributed mesh:<br><br>DM idm=NULL;<br>ierr = DMPlexUninterpolate(*dm, &idm);CHKERRQ(ierr);<br>ierr = DMPlexCopyCoordinates(*dm, idm);CHKERRQ(ierr);<br>ierr = DMPlexCopyLabels(*dm, idm);CHKERRQ(ierr);<br>if (size>1) {<br>  ierr = UninterpolateSF(*dm,idm);CHKERRQ(ierr);<br>}<br>ierr = DMDestroy(dm);CHKERRQ(ierr);<br>*dm = idm;<br></div><div><br></div>Can you have a look into the source code I attached and see what may potentially be the issue? <br><br></div>Thanks,<br><div><div><br></div></div></div><div class="HOEnZb"><div class="h5"><div class="gmail_extra"><br><div class="gmail_quote">On Mon, Apr 13, 2015 at 4:28 PM, Matthew Knepley <span dir="ltr"><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div class="gmail_extra"><div class="gmail_quote"><span>On Mon, Apr 13, 2015 at 4:22 PM, Justin Chang <span dir="ltr"><<a href="mailto:jchang27@uh.edu" target="_blank">jchang27@uh.edu</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr"><div>Is there an example somewhere that does something similar to this?<br></div></div></blockquote><div><br></div></span><div>There is DMPlexShiftSF_Internal() in plexsubmesh.c</div><div><br></div><div>  Matt</div><div><div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr"><div></div>Thanks,<br></div><div><div><div class="gmail_extra"><br><div class="gmail_quote">On Mon, Apr 13, 2015 at 12:31 PM, Matthew Knepley <span dir="ltr"><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr"><div class="gmail_extra"><div class="gmail_quote"><span>On Sat, Apr 11, 2015 at 6:55 PM, Justin Chang <span dir="ltr"><<a href="mailto:jchang27@uh.edu" target="_blank">jchang27@uh.edu</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr"><div><div><div>Hello,<br><br></div>When I call DMPlexUninterpolate(...) on a distributed mesh (say 2 processors), it seems to overwrite the "ghost" points (i.e., the points not locally owned by the processor) and treats all points as if they are local to the processor.<br></div></div></div></blockquote><div><br></div></span><div>Yes, I wrote Uninterpolate() just for testing, and do not currently handle the SF. I put it on my TODO list.</div><div><br></div><div>Its not hard if you want to try. You just filter out any points that are not cells and vertices from the SF, so</div><div><br></div><div>  PetscSFGetGraph()</div><div>  for (leaves)</div><div>    if leaf not a cell or vertex, skip</div><div>  PetscSFSetGraph()</div><div><br></div><div>  Thanks,</div><div><br></div><div>    Matt</div><div><div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr"><div><div></div>Say I have this 2D cell-vertex mesh:<br><br></div>14----15-----16<div>| \  5  | \ 7   |<br>|   \   |   \   |<br>|  4  \ |  6  \ |</div><div>11----12-----13<br>| \  1  | \ 3   |<br>|   \   |   \   |<br>|  0  \ |  2  \ |</div><div>8------9------10<br><br></div><div>Which results in the following DM:<br><br>DM Object: 2 MPI processes<br>  type: plex<br>DM_0x84000004_0 in 2 dimensions:<br>  0-cells: 9 0<br>  2-cells: 8 0<br>Labels:<br>  marker: 1 strata of sizes (8)<br>  depth: 2 strata of sizes (9, 8)<br><br>I proceed by interpolating this DM:<br><br>DM Object: 2 MPI processes<br>  type: plex<br>DM_0x84000004_1 in 2 dimensions:<br>  0-cells: 9 0<br>  1-cells: 16 0<br>  2-cells: 8 0<br>Labels:<br>  marker: 1 strata of sizes (16)<br>  depth: 3 strata of sizes (9, 16, 8)<br><br></div><div>Then distributing across 2 processors:<br><br>DM Object:Parallel Mesh 2 MPI processes<br>  type: plex<br>Parallel Mesh in 2 dimensions:<br>  0-cells: 6 6<br>  1-cells: 9 9<br>  2-cells: 4 4<br>Labels:<br>  marker: 1 strata of sizes (9)<br>  depth: 3 strata of sizes (6, 9, 4)<br><br></div><div>I have the option of uniformly refining the mesh here but I choose not to for now. If my dofs are vertex based, then the global size of my DM vector is 9 and the local sizes for ranks 0 and 1 are 3 and 6 respectively. However, if I choose to uninterpolate the mesh by calling DMPlexUninterpolate(...), I get this:<br><br>DM Object: 2 MPI processes<br>  type: plex<br>DM_0x84000004_2 in 2 dimensions:<br>  0-cells: 6 6<br>  2-cells: 4 4<br>Labels:<br>  marker: 1 strata of sizes (5)<br>  depth: 2 strata of sizes (6, 4)<br><br></div><div>And the global size of my DM vector becomes 12 and the local size for both ranks is 6. It looks like the ghost points in rank 0 have been duplicated, which is not suppose to happen.<br><br></div><div>Is there a way to capture the ghost point information when uninterpolating the DM?<br><br></div><div>Thanks,<span><font color="#888888"><br></font></span></div><span><font color="#888888"><br>-- <br><div><div dir="ltr"><div><div><div>Justin Chang<br></div>PhD Candidate, Civil Engineering - Computational Sciences<br></div>University of Houston, Department of Civil and Environmental Engineering<br></div>Houston, TX 77004<br><a href="tel:%28512%29%20963-3262" value="+15129633262" target="_blank">(512) 963-3262</a><br></div></div></font></span></div>
</blockquote></div></div></div><span><font color="#888888"><br><br clear="all"><div><br></div>-- <br><div>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>
</font></span></div></div>
</blockquote></div><br></div><br clear="all"><br>-- <br><div><div dir="ltr"><div><div><div>Justin Chang<br></div>PhD Candidate, Civil Engineering - Computational Sciences<br></div>University of Houston, Department of Civil and Environmental Engineering<br></div>Houston, TX 77004<br><a href="tel:%28512%29%20963-3262" value="+15129633262" target="_blank">(512) 963-3262</a><br></div></div>
</div></div></blockquote></div></div></div><div><div><br><br clear="all"><div><br></div>-- <br><div>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></div></div></div>
</blockquote></div><br></div><br clear="all"><br>-- <br><div><div dir="ltr"><div><div><div>Justin Chang<br></div>PhD Candidate, Civil Engineering - Computational Sciences<br></div>University of Houston, Department of Civil and Environmental Engineering<br></div>Houston, TX 77004<br><a href="tel:%28512%29%20963-3262" value="+15129633262" target="_blank">(512) 963-3262</a><br></div></div>
</div></div></blockquote></div><br><br clear="all"><div><br></div>-- <br><div class="gmail_signature">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></div>