<div dir="ltr"><div dir="ltr">On Tue, 5 Feb 2019 at 15:27, Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-style:solid;border-left-color:rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div dir="ltr">On Tue, Feb 5, 2019 at 9:47 AM Andrew Parker via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-style:solid;border-left-color:rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div dir="ltr"><div dir="ltr">Does anyone have a MWE for DMPlexCreateCGNS to use in parallel? Ideally, read parallel, distribute in parallel, construct ghost cells (for parallel comms and halos for physical boundaries)? It's for a cell-centered solver working with cgns meshes.  Is there any limitation on cell-types?</div></div></div></blockquote><div><br></div><div>1) CGNS is a terrible format, unfortunately. Are you sure about committing to it? Always best to question design decisions before lots of code has been written.</div></div></div></blockquote><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-style:solid;border-left-color:rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div class="gmail_quote"><div><br></div><div>2) DMPlexCreateCGNS() probably works, but it has not been tested exhaustively.</div></div></div></blockquote><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-style:solid;border-left-color:rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div class="gmail_quote"><div><br></div><div>3) It does not work in parallel, and is unlikely to in the near future. In order to read in parallel, you should be able to nicely select</div><div>    blocks from disk in the format, and also have local BC specification (rather than global numbers). MED is a nice format like this.</div><div>    It is the CASCADE format, and GMsh uses it internally. We can now read MED in parallel. Without too much work, we could also</div><div>    read ExodusII in parallel I think.</div></div></div></blockquote><div> </div><div>Do you have a MWE for either of those two formats then?  </div><div><br></div><div>Key for me is to be in parallel by the time the setup stage is finished, if there is some serial reading here and there at the initial phase that's fine.  For a cell-centered code I just want to go from file read to parallel distribute, parallel ghost-cell creation, physical halo creation, and boundary condition (faces, edge in 2D) marking re-using as much of the petsc DMPlex system as possible.  </div><div><br></div><div>My main requirement with any file format is the ability in a mesh/cad package to mark complex groups of edges/faces/cells together (on the boundary but also groupings of cells in the volume) to infer groups and therefore bcs, or material regions etc.  I am aware med/exo provide this: my reason for picking cgns is that is central to their format also.  I need to ensure for complex geometry and meshes, that if the file format does encode groupings/bcs then these are accessible after the parallel distribution phase in complete.  It would help me massively if there were an example that does all of this....I'll jump ship on file format to get that.</div><div><br></div><div>Thanks,</div><div>Andy   </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-style:solid;border-left-color:rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div class="gmail_quote"><div><br></div><div>4) Once read, you get a regular Plex, so you can redistribute, make ghost cells, etc. as you can for any Plex mesh.</div><div><br></div><div>5) As long as all you want Plex to do is manage topology and field data, then you can have whatever mix of cell types you want.</div><div><br></div><div>6) Limitations on cell types come from routines that calculate geometric quantities, integrate, etc.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div></div></div></blockquote><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-style:solid;border-left-color:rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div class="gmail_quote"><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-style:solid;border-left-color:rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div dir="ltr"><div>Thanks,</div><div>Andy</div><div><br></div></div></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail-m_-5642496311768580030gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><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><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>
</blockquote></div></div>