<html><head><style type='text/css'>p { margin: 0; }</style></head><body><div style='font-family: times new roman,new york,times,serif; font-size: 12pt; color: #000000'><br>Hello,<br>sorry, we talked about it few months ago.<br>For the moment, it is assumed in DMPlexCreateFromCellList function that cells have all the same number of vertices.<br>In other words, it means that meshes containing one and only one type of cell (from a topological viewpoint) can be handled by this function only.<br>However, several mesh generators are able to discretize geometries by means of different kinds of cells.<br>Resulting meshes are said to be mixed or hybrid because they contain different types of cells.<br>To create a DMPlex structure from such a mesh, one needs to slightly modify DMPlexCreateFromCellList prototype.<br>To my mind, the simplest way to do that is to redefine numCorners as an integer array instead of a simple integer value.<br>Its length should be numCells and it should contain the number of vertices for each cell.<br>Note that the length of "cells" parameter should be modified as well.<br>But, of course, if you have a better idea about the new prototype, it's ok for me.<br>Best regards,<br>Cédric Doucet<br><br><br><br><hr id="zwchr"><blockquote style="border-left:2px solid #1010FF;margin-left:5px;padding-left:5px;color:#000;font-weight:normal;font-style:normal;text-decoration:none;font-family:Helvetica,Arial,sans-serif;font-size:12pt;"><b>De: </b>"Matthew Knepley" <knepley@gmail.com><br><b>À: </b>"Cedric Doucet" <cedric.doucet@inria.fr><br><b>Cc: </b>petsc-users@mcs.anl.gov<br><b>Envoyé: </b>Jeudi 2 Janvier 2014 17:01:24<br><b>Objet: </b>Re: [petsc-users] DMPlexCreateFromCellList and hybrid meshes<br><br><div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Thu, Jan 2, 2014 at 8:08 AM, Cedric Doucet <span dir="ltr"><<a href="mailto:cedric.doucet@inria.fr" target="_blank">cedric.doucet@inria.fr</a>></span> wrote:<br>
<blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div><div style="font-size:12pt;font-family:times new roman,new york,times,serif"><br>Hello,<br><br>indeed, I need to have a generalized version of DMPlexCreateFromCellList function which can take an hybrid mesh as an input parameter.<br>
</div></div></blockquote><div><br></div><div>What interface?</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div><div style="font-size:12pt;font-family:times new roman,new york,times,serif">
Modern mesh generators often provide mixed generation functionalities (e.g. gmsh which can generate pyramids and prisms).<br>As a matter of fact, cells are never numbered independently of vertices in hybrid meshes because different generators are used (structured and unstructured ones).<br>
>From a practical point of view, it is far much easier to let each generator number its own cells without renumbering shared vertices.<br>Therefore, cells of hybrid meshes are sorted in a canonical way: cells of generator 1, cells of generator 2, etc.<br>
</div></div></blockquote><div><br></div><div>I don't understand what you want here. Its imprecise.</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><div style="font-size:12pt;font-family:times new roman,new york,times,serif">I deeply think that this new version is very useful for real-life applications.<br>Please let me know when this new version will be available.<br>
<br>Thank you very much,<br><br>Cédric Doucet<br><br><br><br><br><br><hr><blockquote style="padding-left:5px;font-size:12pt;font-style:normal;margin-left:5px;font-family:Helvetica,Arial,sans-serif;text-decoration:none;font-weight:normal;border-left:2px solid #1010ff">
<b>De: </b>"Matthew Knepley" <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br><b>À: </b>"Cedric Doucet" <<a href="mailto:cedric.doucet@inria.fr" target="_blank">cedric.doucet@inria.fr</a>><br>
<b>Cc: </b><a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a><br><b>Envoyé: </b>Jeudi 31 Octobre 2013 14:13:14<br><b>Objet: </b>Re: [petsc-users] DMPlexCreateFromCellList and hybrid meshes<br>
<br><div dir="ltr">On Thu, Oct 31, 2013 at 5:20 AM, Cedric Doucet <span dir="ltr"><<a href="mailto:cedric.doucet@inria.fr" target="_blank">cedric.doucet@inria.fr</a>></span> wrote:<br><div class="gmail_extra"><div class="gmail_quote">

<blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div><div style="font-size:12pt;font-family:times new roman,new york,times,serif"><br>Thank you very much for your help.<br>

Do you plan to add this code to Petsc?<br>It would be better for us to call an official Petsc function than maintaining a Petsc-style code into our C++ software.<br></div></div></blockquote><div><br></div><div>If this is really what you need, then I will put it in. However, I intended this to be used primarily for mesh generators, since they</div>

<div>tend to number cells and vertices independently and have the int and double types hardwired.</div><div><br></div><div>If you are writing your own code, please consider DMPlexCreateFromDAG() which has the flexibility for cell topology that you</div>

<div>want. The main differences are that it uses PETSc types and a single contiguous numbering. Would this work for you?</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><div style="font-size:12pt;font-family:times new roman,new york,times,serif">Best regards,<br><br>Cédric Doucet<br><br><br><br><hr><blockquote style="padding-left:5px;font-size:12pt;font-style:normal;margin-left:5px;font-family:Helvetica,Arial,sans-serif;text-decoration:none;font-weight:normal;border-left:2px solid #1010ff">

<b>De: </b>"Matthew Knepley" <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br><b>À: </b>"Cedric Doucet" <<a href="mailto:cedric.doucet@inria.fr" target="_blank">cedric.doucet@inria.fr</a>><br>

<b>Cc: </b><a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a><br><b>Envoyé: </b>Mardi 29 Octobre 2013 17:24:03<br><b>Objet: </b>Re: [petsc-users] DMPlexCreateFromCellList and hybrid meshes<br>

<br><div dir="ltr">On Tue, Oct 29, 2013 at 10:24 AM, Cedric Doucet <span dir="ltr"><<a href="mailto:cedric.doucet@inria.fr" target="_blank">cedric.doucet@inria.fr</a>></span> wrote:<br><div class="gmail_extra"><div class="gmail_quote">


<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><div style="font-size:12pt;font-family:'times new roman','new york',times,serif">


Hello,<br><br>I have a short question about the code I need to modify.<br>As far as I understand, numCorner should be replaced by an array numCorners[numCells] containing the number of vertices of each cell.<br>The body of DMPlexBuildFromCellList_Private function becomes<br>


---------------------------------------------------------------------------------------------------------------------------------------<br>ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr);</div></div></blockquote>


<div>    maxCorners = 0; </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><div style="font-size:12pt;font-family:'times new roman','new york',times,serif">



<div>  for (c = 0; c < numCells; ++c) {</div><div>    ierr = DMPlexSetConeSize(dm, c, numCorners[c]);CHKERRQ(ierr);</div></div></div></blockquote><div>         maxCorners = PetscMax(maxCorners, numCorners[c]); </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 style="font-size:12pt;font-family:'times new roman','new york',times,serif"><div>  }</div><div>  ierr = DMSetUp(dm);CHKERRQ(ierr);</div><div></div></div></blockquote><div> </div><div><span style="font-family:'times new roman','new york',times,serif;font-size:15.555556297302246px">    ierr = DMGetWorkArray(dm, maxCorners, PETSC_INT, &cone);CHKERRQ(ierr);</span><span style="font-family:'times new roman','new york',times,serif;font-size:12pt">  </span></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 style="font-size:12pt;font-family:'times new roman','new york',times,serif">



<div>  for (c = 0, off = 0; c < numCells; ++c) {<br></div><div>    for (p = 0; p < numCorners[c]; ++p) {</div><div>      cone[p] = cells[off+p]+numCells;</div><div>    }<br></div><div>    ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);</div>


</div></blockquote><div>         off += numCorners[c]; </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 style="font-size:12pt;font-family:'times new roman','new york',times,serif">
<div>  }</div></div></blockquote><div><span style="font-family:'times new roman','new york',times,serif;font-size:15.555556297302246px">    ierr = DMRestoreWorkArray(dm, maxCorners, PETSC_INT, &cone);CHKERRQ(ierr);</span><span style="font-family:'times new roman','new york',times,serif;font-size:12pt">  </span></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 style="font-size:12pt;font-family:'times new roman','new york',times,serif">


</div></blockquote><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 style="font-size:12pt;font-family:'times new roman','new york',times,serif">


<div> </div><div>  ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);</div><div>  ierr = DMPlexStratify(dm);CHKERRQ(ierr);<br></div></div></blockquote><div><br></div><div>   Matt</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 style="font-size:12pt;font-family:'times new roman','new york',times,serif"><div>------------------------------------------------------------------------------------------------------------------------------------------<br>


However, I am not quite sure that DMGetWorkArray and DMRestoreWorkArray should be used like this.<br>I guess that DMRestoreWorkArray behaves like free function in ansi c but I am not sure.<br>Do you think that calling DMGetWorkArray and DMRestoreWorkArray inside a loop is a good thing from the point of view of computational time and memory management?<br>


<br>Note: the best way to do this may be to first sort cells by numCorners to avoid many calls to DMGetWorkArray and DMRestoreWorkArray. This is actually what I have in my own code but I would like to maintain Petsc's philosophy.<br>


<br>Thanks,<br><br>Cédric<br></div><br><br><br><br><br><hr><blockquote style="padding-left:5px;font-size:12pt;font-style:normal;margin-left:5px;font-family:Helvetica,Arial,sans-serif;text-decoration:none;font-weight:normal;border-left-width:2px;border-left-style:solid;border-left-color:rgb(16,16,255)">


<b>De: </b>"Matthew Knepley" <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br><b>À: </b>"Cedric Doucet" <<a href="mailto:cedric.doucet@inria.fr" target="_blank">cedric.doucet@inria.fr</a>><br>


<b>Cc: </b><a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a><br><b>Envoyé: </b>Vendredi 25 Octobre 2013 16:31:59<br><b>Objet: </b>Re: [petsc-users] DMPlexCreateFromCellList and hybrid meshes<br>


<br><div dir="ltr">On Fri, Oct 25, 2013 at 7:23 AM, Cedric Doucet <span dir="ltr"><<a href="mailto:cedric.doucet@inria.fr" target="_blank">cedric.doucet@inria.fr</a>></span> wrote:<br><div class="gmail_extra"><div class="gmail_quote">



<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><div style="font-size:12pt;font-family:'times new roman','new york',times,serif">



Hello,<br><br>I've noticed that DMPlexCreateFromCellList assumes that cells have the same number of vertices (numcorners argument).<br>What should be done when one wants to create a DMPlex from a mesh containing different types of cells?<br>



Does one have to create several DMPlex structures and merge them?<br>Does one have to create a unique DMPlex by hand?<br></div></div></blockquote><div><br></div><div>The code is very short:</div><div><br></div><div><div>


  ierr = DMPlexSetChart(dm, 0, numCells+numVertices);CHKERRQ(ierr);</div>
<div>  for (c = 0; c < numCells; ++c) {</div><div>    ierr = DMPlexSetConeSize(dm, c, numCorners);CHKERRQ(ierr);</div><div>  }</div><div>  ierr = DMSetUp(dm);CHKERRQ(ierr);</div><div>  ierr = DMGetWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);</div>



<div>  for (c = 0; c < numCells; ++c) {</div><div>    for (p = 0; p < numCorners; ++p) {</div><div>      cone[p] = cells[c*numCorners+p]+numCells;</div><div>    }</div><div>    ierr = DMPlexSetCone(dm, c, cone);CHKERRQ(ierr);</div>



<div>  }</div><div>  ierr = DMRestoreWorkArray(dm, numCorners, PETSC_INT, &cone);CHKERRQ(ierr);</div><div>  ierr = DMPlexSymmetrize(dm);CHKERRQ(ierr);</div><div>  ierr = DMPlexStratify(dm);CHKERRQ(ierr);</div></div><div>



<br></div><div>This code is all in plexcreate.c. If you want different cells, you can change numCorners for each cell.</div><div>I could make a convenience form if necessary.</div><div><br></div><div>   Matt</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><div style="font-size:12pt;font-family:'times new roman','new york',times,serif">



<br>Thank you very much for your help.<br><br>Cédric<br><br><br><br></div></div></blockquote></div><br><br clear="all"><span><font color="#888888"><span><font color="#888888"><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
</font></span></font></span></div></div><span><font color="#888888">
</font></span></blockquote><span><font color="#888888"><br></font></span></div></blockquote></div><span><font color="#888888"><br><br clear="all"><span class="HOEnZb"><font color="#888888"><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
</font></span></font></span></div></div><span class="HOEnZb"><font color="#888888">
</font></span></blockquote><span class="HOEnZb"><font color="#888888"><br></font></span></div></div></blockquote></div><span class="HOEnZb"><font color="#888888"><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
</font></span></div></div>
</blockquote><br></div></div></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>
</blockquote><br></div></body></html>