<div dir="ltr">Hi <br><br>Thanks for the response. I still have a bit of confusion here. To make sure I follow I have setup a code block (at end) with options for overlap and CreateGhostCells. <br><br>Using a 2x2x2 Box mesh I've seen the following behavior. <br><br>When I use overlap=1 I get approximately what I would expect. <br>If I run on 1 processor there are only 8cells and everything is normal <br>If I run on 2 processors each processor has 8cells (4 from the split domain and an extra 4 extending into the neighboring cell) <br>My question here is how do I then differentiate the overlap cells from each processor so that I can synchronize the extra cells as I iterate?<br>Both GetHeightStratum GetSimplex will identify the 8cells with no differentiation between the partitions vs the overlap. <br>I believe this is the strategy that will get me the result I want where I can synch the overlaps so that when I iterate through the "owned cells" the stencil can reach into the overlap cells to get the information they need?<br><br><br>When I use ghostCellGenerate I also get what I would expect<br>If I run on 1 processor each of the outside faces of the cube gets an extra cell attached so I get Box Cells 0 through 8 and then Ghost Cells 8 through 32. <br>If I run on 2 processors only the outside faces of the cub get the extra cells so I get Box Cells 0 through 4 and then Ghost Cells 4 through 16 <br>So I don't think this is what I want for my intended use?<br><br>Thanks Again<br><br>Sincerely<br>Nicholas <br><br><div style="color:rgb(255,255,255);background-color:rgb(0,36,81);font-family:Consolas,"Courier New",monospace;font-size:14px;line-height:19px;white-space:pre"><div> <span style="color:rgb(255,238,173)">PetscInt</span> <span style="color:rgb(255,157,164)">facecount</span> <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(255,197,143)">10</span>;</div><div> PetscInt i, dim <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(255,197,143)">3</span>, overlap <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(255,197,143)">0</span>;</div><div> <span style="color:rgb(255,157,164)">ierr</span> <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(187,218,255)">PetscOptionsGetInt</span>(<span style="color:rgb(187,218,255)">NULL</span>, <span style="color:rgb(187,218,255)">NULL</span>, <span style="color:rgb(209,241,169)">"-facecount"</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">facecount</span>, <span style="color:rgb(187,218,255)">NULL</span>);</div><div> <span style="color:rgb(187,218,255)">CHKERRQ</span>(<span style="color:rgb(255,157,164)">ierr</span>);</div><div> <span style="color:rgb(235,187,255)">for</span> (<span style="color:rgb(255,157,164)">i</span> <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(255,197,143)">0</span>; <span style="color:rgb(255,157,164)">i</span> <span style="color:rgb(153,255,255)"><</span> <span style="color:rgb(255,157,164)">dim</span>; <span style="color:rgb(255,157,164)">i</span><span style="color:rgb(153,255,255)">++</span>)</div><div> <span style="color:rgb(255,157,164)">faces</span>[<span style="color:rgb(255,157,164)">i</span>] <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(255,157,164)">facecount</span>;</div><div> <span style="color:rgb(187,218,255)">PetscViewerASCIIPrintf</span>(<span style="color:rgb(255,157,164)">viewerstdout</span>, <span style="color:rgb(209,241,169)">"facecount is </span>%d<span style="color:rgb(255,197,143)">\n</span><span style="color:rgb(209,241,169)">"</span>, <span style="color:rgb(255,157,164)">facecount</span>);</div><div> <span style="color:rgb(255,238,173)">PetscReal</span> <span style="color:rgb(255,157,164)">lower</span>[<span style="color:rgb(255,197,143)">3</span>] <span style="color:rgb(153,255,255)">=</span> {<span style="color:rgb(255,197,143)">0.0</span>, <span style="color:rgb(255,197,143)">0.0</span>, <span style="color:rgb(255,197,143)">0.0</span>};</div><div> <span style="color:rgb(255,238,173)">PetscReal</span> <span style="color:rgb(255,157,164)">upper</span>[<span style="color:rgb(255,197,143)">3</span>] <span style="color:rgb(153,255,255)">=</span> {<span style="color:rgb(255,197,143)">2.0</span>, <span style="color:rgb(255,197,143)">2.0</span>, <span style="color:rgb(255,197,143)">2.0</span>};</div><div> <span style="color:rgb(255,157,164)">ierr</span> <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(255,157,164)">DMPlexCreateBoxMesh</span>(<span style="color:rgb(255,157,164)">comm</span>, <span style="color:rgb(255,157,164)">dim</span>, <span style="color:rgb(255,157,164)">simplex</span>, <span style="color:rgb(255,157,164)">faces</span>, <span style="color:rgb(255,157,164)">lower</span>, <span style="color:rgb(255,157,164)">upper</span>,<span style="color:rgb(114,133,183)"> /* periodicity */</span> <span style="color:rgb(187,218,255)">NULL</span>, <span style="color:rgb(255,157,164)">dmInterped</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">dm</span>);</div><div> <span style="color:rgb(187,218,255)">CHKERRQ</span>(<span style="color:rgb(255,157,164)">ierr</span>);</div><div> <span style="color:rgb(255,157,164)">ierr</span> <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(187,218,255)">DMPlexDistribute</span>(<span style="color:rgb(255,157,164)">dm</span>, <span style="color:rgb(255,157,164)">overlap</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">distributionSF</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">dmDist</span>);</div><div> <span style="color:rgb(187,218,255)">CHKERRQ</span>(<span style="color:rgb(255,157,164)">ierr</span>);</div><div> <span style="color:rgb(235,187,255)">if</span> (<span style="color:rgb(255,157,164)">dmDist</span>)</div><div> {</div><div> <span style="color:rgb(255,157,164)">ierr</span> <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(187,218,255)">DMDestroy</span>(<span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">dm</span>);</div><div> <span style="color:rgb(187,218,255)">CHKERRQ</span>(<span style="color:rgb(255,157,164)">ierr</span>);</div><div> <span style="color:rgb(255,157,164)">dm</span> <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(255,157,164)">dmDist</span>;</div><div> }</div><br><div> PetscInt ghostCells;</div><br><div> PetscBool use_ghostcells <span style="color:rgb(153,255,255)">=</span> PETSC_TRUE;</div><div> ierr <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(187,218,255)">PetscOptionsGetBool</span>(<span style="color:rgb(255,197,143)">NULL</span>, <span style="color:rgb(255,197,143)">NULL</span>, <span style="color:rgb(209,241,169)">"-use_ghostcells"</span>, <span style="color:rgb(153,255,255)">&</span>use_ghostcells, <span style="color:rgb(255,197,143)">NULL</span>);</div><div> DM dmGhost;</div><div> <span style="color:rgb(235,187,255)">if</span> (use_ghostcells)</div><div> {</div><div> ierr <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(187,218,255)">DMPlexConstructGhostCells</span>(dm, <span style="color:rgb(255,197,143)">NULL</span>, <span style="color:rgb(153,255,255)">&</span>ghostCells, <span style="color:rgb(153,255,255)">&</span>dmGhost);</div><div> <span style="color:rgb(187,218,255)">CHKERRQ</span>(ierr);</div><div> <span style="color:rgb(235,187,255)">if</span> (dmGhost)</div><div> {</div><div> <span style="color:rgb(187,218,255)">PetscViewerASCIIPrintf</span>(viewerstdout, <span style="color:rgb(209,241,169)">"</span>%d<span style="color:rgb(209,241,169)">ghost Cells generated</span><span style="color:rgb(255,197,143)">\n</span><span style="color:rgb(209,241,169)">"</span>, ghostCells);</div><div> ierr <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(187,218,255)">DMDestroy</span>(<span style="color:rgb(153,255,255)">&</span>dm);</div><div> <span style="color:rgb(187,218,255)">CHKERRQ</span>(ierr);</div><div> dm <span style="color:rgb(153,255,255)">=</span> dmGhost;</div><div> }</div><div> <span style="color:rgb(235,187,255)">else</span></div><div> <span style="color:rgb(187,218,255)">PetscViewerASCIIPrintf</span>(viewerstdout, <span style="color:rgb(209,241,169)">"ghost Cells not generated</span><span style="color:rgb(255,197,143)">\n</span><span style="color:rgb(209,241,169)">"</span>);</div><div> }</div><div> <span style="color:rgb(235,187,255)">else</span></div><div> <span style="color:rgb(187,218,255)">PetscViewerASCIIPrintf</span>(viewerstdout, <span style="color:rgb(209,241,169)">"ghost Cells not generated</span><span style="color:rgb(255,197,143)">\n</span><span style="color:rgb(209,241,169)">"</span>);</div><br><div> <span style="color:rgb(255,238,173)">PetscInt</span> <span style="color:rgb(255,157,164)">c0</span>, <span style="color:rgb(255,157,164)">c1</span>, <span style="color:rgb(255,157,164)">f0</span>, <span style="color:rgb(255,157,164)">f1</span>, <span style="color:rgb(255,157,164)">e0</span>, <span style="color:rgb(255,157,164)">e1</span>, <span style="color:rgb(255,157,164)">v0</span>, <span style="color:rgb(255,157,164)">v1</span>, <span style="color:rgb(255,157,164)">s0, s1, g</span>0, <span style="color:rgb(255,157,164)">g1</span>;</div><br><div> <span style="color:rgb(187,218,255)">DMPlexGetHeightStratum</span>(<span style="color:rgb(255,157,164)">dm</span>, <span style="color:rgb(255,197,143)">0</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">c0</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">c1</span>);</div><div> <span style="color:rgb(187,218,255)">DMPlexGetHeightStratum</span>(<span style="color:rgb(255,157,164)">dm</span>, <span style="color:rgb(255,197,143)">1</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">f0</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">f1</span>);</div><div> <span style="color:rgb(187,218,255)">DMPlexGetHeightStratum</span>(<span style="color:rgb(255,157,164)">dm</span>, <span style="color:rgb(255,197,143)">2</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">e0</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">e1</span>);</div><div> <span style="color:rgb(187,218,255)">DMPlexGetHeightStratum</span>(<span style="color:rgb(255,157,164)">dm</span>, <span style="color:rgb(255,197,143)">3</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">v0</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">v1</span>);</div><div> <span style="color:rgb(187,218,255)">DMPlexGetSimplexOrBoxCells</span>(<span style="color:rgb(255,157,164)">dm</span>, <span style="color:rgb(255,197,143)">0</span>, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">s</span>0, <span style="color:rgb(153,255,255)">&</span><span style="color:rgb(255,157,164)">s</span>1);</div><div> <span style="color:rgb(187,218,255)">DMPlexGetGhostCellStratum</span>(dm, <span style="color:rgb(153,255,255)">&</span>g0, <span style="color:rgb(153,255,255)">&</span>g1);</div><div> <span style="color:rgb(187,218,255)">PetscViewerASCIIPrintf</span>(<span style="color:rgb(255,157,164)">viewerstdout</span>, <span style="color:rgb(209,241,169)">"Cells: p0:</span>%d<span style="color:rgb(209,241,169)"> ,p1:</span>%d<span style="color:rgb(209,241,169)"> </span><span style="color:rgb(255,197,143)">\n</span><span style="color:rgb(209,241,169)">"</span>, <span style="color:rgb(255,157,164)">c0</span>, <span style="color:rgb(255,157,164)">c1</span>);</div><div> <span style="color:rgb(187,218,255)">PetscViewerASCIIPrintf</span>(<span style="color:rgb(255,157,164)">viewerstdout</span>, <span style="color:rgb(209,241,169)">"Faces: p0:</span>%d<span style="color:rgb(209,241,169)"> ,p1:</span>%d<span style="color:rgb(209,241,169)"> </span><span style="color:rgb(255,197,143)">\n</span><span style="color:rgb(209,241,169)">"</span>, <span style="color:rgb(255,157,164)">f0</span>, <span style="color:rgb(255,157,164)">f1</span>);</div><div> <span style="color:rgb(187,218,255)">PetscViewerASCIIPrintf</span>(<span style="color:rgb(255,157,164)">viewerstdout</span>, <span style="color:rgb(209,241,169)">"Edges: p0:</span>%d<span style="color:rgb(209,241,169)"> ,p1:</span>%d<span style="color:rgb(209,241,169)"> </span><span style="color:rgb(255,197,143)">\n</span><span style="color:rgb(209,241,169)">"</span>, <span style="color:rgb(255,157,164)">e0</span>, <span style="color:rgb(255,157,164)">e1</span>);</div><div> <span style="color:rgb(187,218,255)">PetscViewerASCIIPrintf</span>(<span style="color:rgb(255,157,164)">viewerstdout</span>, <span style="color:rgb(209,241,169)">"Vertices: p0:</span>%d<span style="color:rgb(209,241,169)"> ,p1:</span>%d<span style="color:rgb(209,241,169)"> </span><span style="color:rgb(255,197,143)">\n</span><span style="color:rgb(209,241,169)">"</span>, <span style="color:rgb(255,157,164)">v0</span>, <span style="color:rgb(255,157,164)">v1</span>);</div><div> <span style="color:rgb(187,218,255)">PetscViewerASCIIPrintf</span>(<span style="color:rgb(255,157,164)">viewerstdout</span>, <span style="color:rgb(209,241,169)">"Box Cells: p0:</span>%d<span style="color:rgb(209,241,169)"> ,p1:</span>%d<span style="color:rgb(209,241,169)"> </span><span style="color:rgb(255,197,143)">\n</span><span style="color:rgb(209,241,169)">"</span>, <span style="color:rgb(255,157,164)">s</span>0, <span style="color:rgb(255,157,164)">s</span>1);</div><div> <span style="color:rgb(187,218,255)">PetscViewerASCIIPrintf</span>(viewerstdout, <span style="color:rgb(209,241,169)">"Ghost Cells: p0:</span>%d<span style="color:rgb(209,241,169)"> ,p1:</span>%d<span style="color:rgb(209,241,169)"> </span><span style="color:rgb(255,197,143)">\n</span><span style="color:rgb(209,241,169)">"</span>, g0, g1);</div></div><br><br></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Sun, Jul 3, 2022 at 12:42 AM Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div dir="ltr">On Sat, Jul 2, 2022 at 6:20 PM Nicholas Arnold-Medabalimi <<a href="mailto:narnoldm@umich.edu" target="_blank">narnoldm@umich.edu</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">Dear PETSc users, <br><br>I have some general novice questions regarding how to properly access a distributed DMPLEX vector that has ghost cells for a finite volume use case.<br><br>My process for setup<br>1) generated a simple box mesh using DMPlexCreateBoxMesh <br>2) distribute the mesh using DMPlexDistribute (at this point everything looks fine when I View (as paraview file) and I can see the partitions)<br>3) I generate a section with 5 variables located at the cells. (I don't think it matters if I use 5 fields or 1 field with 5 components?) I can alternatively do this using the PetscFV variable with 5 components. <br>4) I then have the dm setup properly (with no ghost cells at this point) and I use DMCreateGlobalVector and DMCreateLocalVector to get the Vectors to work with<br>5) I'm just setting a quadratic IC. The overall loop over the cells is constructed via getting the cell stratum which has the Cell point list and then using the SectionOffset to assign each of the components at each cell using DMPlexComputeCellGeometryFVM to get the cell centroid.<br><br>As far as I can tell this is an adequate solution but now I'm moving into learning how to access neighboring cells for flux calculations which isn't just the local cell but requires access to neighboring cells which at some point will be located on neighboring processors and synched to be accessible. When I did this for a DMDA it was pretty straightforward insofar as I just needed to call DMDAGetCorners and the LocalVector indexing would allow the stencil to extend gracefully into the ghost values. <br><br>My general question is how do I achieve a similar loop structure for DMPlex as far as the indexing over the owned cells while being able to access the neighbor halo cells. I have been looking at example 52 and 11 in the ts tutorial examples but I'm struggling to extract exactly what I need. <br><br>To boil it down into 3 starting questions<br>1) What is the preferred way to setup the cells and then access the internal cells for say N cell halos for the boundaries between processors in a DMPlex?<br></div></blockquote><div><br></div><div>If you want a halo during distribution, you can give the "overlap" argument to DMPlexDistribute().</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">2) What is the distinction between the overlap argument in DMPlexDistribute and DMPlexConstructGhostCells (which I see used in the examples)? (at least based on my tests it seems GhostCells is generating cells at the overall boundaries not the internal partition boundaries?)<br></div></blockquote><div><br></div><div>The examples use both. DMPlexConstructGhostCells() puts a layer of cells around the domain boundary that are special. They have no cell boundary, just the face they are attached to. The overlap adds true cells to each domain across the parallel boundary.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">3) Once I do have a DMPlex object that has the appropriate halos, how would I get the access range for those cells? Based on my current assumption the Stratum range for the Cells would include the ghost cells as well which I would like to avoid?</div></blockquote><div><br></div><div>I use <a href="https://petsc.org/main/docs/manualpages/DMPLEX/DMPlexGetSimplexOrBoxCells" target="_blank">https://petsc.org/main/docs/manualpages/DMPLEX/DMPlexGetSimplexOrBoxCells</a> to avoid the FV ghost cells. However, you can also do it by hand by checking the type, or use <a href="https://petsc.org/main/docs/manualpages/DMPLEX/DMPlexGetGhostCellStratum" target="_blank">https://petsc.org/main/docs/manualpages/DMPLEX/DMPlexGetGhostCellStratum</a></div><div><br></div><div> Thanks,</div><div><br></div><div> Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"> (Ultimately in the examples I see cases where there is the synch step for Local to Global and vice versa but I'm trying to just get to the point where I have the DM setup properly to handle the internal mesh boundaries and have the correct indexing from the StratumGet commands.)<br><br>Any advice and corrections would be welcome. <br><br>Sincerely<br>Nicholas <br clear="all"><div><br></div>-- <br><div dir="ltr"><div dir="ltr"><div dir="ltr"><div style="font-family:arial;font-size:small"><font color="#000000">Nicholas Arnold-Medabalimi<br></font><br></div></div></div></div></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr"><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><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div style="font-family:arial;font-size:small"><font color="#000000">Nicholas Arnold-Medabalimi<br><br></font><span style="font-family:sans-serif;font-size:14px">Ph.D. Candidate</span><font color="#000000"><br>Computational Aeroscience Lab<br>University of Michigan</font></div></div></div></div></div>