<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
</head>
<body style="word-wrap: break-word; -webkit-nbsp-mode: space; -webkit-line-break: after-white-space;">
<div>
<div style="color: rgb(0, 0, 0); font-family: Calibri, sans-serif; font-size: 14px;">
Hi Matt</div>
<div style="color: rgb(0, 0, 0); font-family: Calibri, sans-serif; font-size: 14px;">
<br>
</div>
<div style="color: rgb(0, 0, 0); font-family: Calibri, sans-serif; font-size: 14px;">
I don’t think the problem is within Petsc - rather somewhere in my code. When I dump the DMPlex using DMView (ascii-info–detail) the ghost mapping seems to be setup correctly.</div>
<div style="color: rgb(0, 0, 0); font-family: Calibri, sans-serif; font-size: 14px;">
<br>
</div>
<div style="color: rgb(0, 0, 0); font-family: Calibri, sans-serif; font-size: 14px;">
Is there a better way to determine if a local point is a ghost point? </div>
<div style="color: rgb(0, 0, 0); font-family: Calibri, sans-serif; font-size: 14px;">
<br>
</div>
<div style="color: rgb(0, 0, 0); font-family: Calibri, sans-serif; font-size: 14px;">
The way I iterate the DMPlex is like this:</div>
<div>
<pre style="background-color: rgb(255, 255, 255);"><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">void </span><font face="Menlo"><span style="font-size: 9pt;">iterateDMPlex(DM dm){<br>    Vec coordinates;<br>    DMGetCoordinatesLocal(dm, &coordinates);<br>    PetscSection defaultSection;<br>    DMGetDefaultSection(dm, &defaultSection);<br>    PetscSection coordSection;<br>    DMGetCoordinateSection(dm, &coordSection);<br>    PetscScalar *coords;<br>    VecGetArray(coordinates, &coords);<br>    DM cdm;<br>    DMGetCoordinateDM(dm, &cdm);<br><br>    </span></font><span style="color: rgb(128, 128, 128); font-family: Menlo; font-size: 9pt; font-style: italic;">// iterate (local) mesh<br></span><span style="color: rgb(128, 128, 128); font-family: Menlo; font-size: 9pt; font-style: italic;">    </span><font face="Menlo"><span style="font-size: 9pt;">PetscInt cellsFrom, cellsTo;<br>    std::string s = </span></font><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">""</span><font face="Menlo"><span style="font-size: 9pt;">;<br>    DMPlexGetHeightStratum(dm, </span></font><span style="color: rgb(0, 0, 255); font-family: Menlo; font-size: 9pt;">0</span><font face="Menlo"><span style="font-size: 9pt;">, &cellsFrom, &cellsTo);<br>    </span></font><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">for </span><font face="Menlo"><span style="font-size: 9pt;">(PetscInt i=cellsFrom;i<cellsTo;i++) {<br>        PetscInt edgesSize;<br>        </span></font><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">const </span><font face="Menlo"><span style="font-size: 9pt;">PetscInt *edgeIndices;<br>        DMPlexGetConeSize(dm, i, &edgesSize);<br>        DMPlexGetCone(dm, i, &edgeIndices);<br>        s = s + </span></font><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">"Element: "</span><font face="Menlo"><span style="font-size: 9pt;">+std::to_string(i)+</span></font><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">"</span><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">\n</span><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">"</span><font face="Menlo"><span style="font-size: 9pt;">;<br><br>        </span></font><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">for </span><font face="Menlo"><span style="font-size: 9pt;">(</span></font><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">int </span><font face="Menlo"><span style="font-size: 9pt;">edgeId = </span></font><span style="color: rgb(0, 0, 255); font-family: Menlo; font-size: 9pt;">0</span><font face="Menlo"><span style="font-size: 9pt;">;edgeId <edgesSize;edgeId ++){ </span></font><span style="color: rgb(128, 128, 128); font-family: Menlo; font-size: 9pt; font-style: italic;">// ignore edge orientation<br></span><span style="color: rgb(128, 128, 128); font-family: Menlo; font-size: 9pt; font-style: italic;">            </span><font face="Menlo"><span style="font-size: 9pt;">PetscInt points = edgeIndices[edgeId];<br>            PetscInt edgePoint = edgeIndices[edgeId];<br>            s = s + </span></font><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">"</span><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">\t</span><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">Edge: "</span><font face="Menlo"><span style="font-size: 9pt;">+std::to_string(edgePoint)+</span></font><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">"</span><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">\n</span><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">"</span><font face="Menlo"><span style="font-size: 9pt;">;<br>            PetscInt vertexSize;<br>            </span></font><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">const </span><font face="Menlo"><span style="font-size: 9pt;">PetscInt *vertexIndices;<br>            DMPlexGetConeSize(dm, edgePoint, &vertexSize);<br>            DMPlexGetCone(dm, edgePoint, &vertexIndices);<br><br>            </span></font><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">for </span><font face="Menlo"><span style="font-size: 9pt;">(</span></font><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">int </span><font face="Menlo"><span style="font-size: 9pt;">j=</span></font><span style="color: rgb(0, 0, 255); font-family: Menlo; font-size: 9pt;">0</span><font face="Menlo"><span style="font-size: 9pt;">;j<vertexSize;j++){<br>                s = s + </span></font><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">"</span><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">\t\t</span><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">Vertex: "</span><font face="Menlo"><span style="font-size: 9pt;">+std::to_string(vertexIndices[j]);<br><br>                s = s + </span></font><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">" coords "</span><font face="Menlo"><span style="font-size: 9pt;">;<br>                PetscScalar* values;<br>                VecGetValuesSection(coordinates, coordSection,vertexIndices[j],&values);<br><br>                </span></font><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">int </span><font face="Menlo"><span style="font-size: 9pt;">dim = </span></font><span style="color: rgb(0, 0, 255); font-family: Menlo; font-size: 9pt;">2</span><font face="Menlo"><span style="font-size: 9pt;">;<br>                </span></font><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">for </span><font face="Menlo"><span style="font-size: 9pt;">(</span></font><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">int </span><font face="Menlo"><span style="font-size: 9pt;">k=</span></font><span style="color: rgb(0, 0, 255); font-family: Menlo; font-size: 9pt;">0</span><font face="Menlo"><span style="font-size: 9pt;">;k<dim;k++){<br>                    </span></font><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">double </span><font face="Menlo"><span style="font-size: 9pt;">coordinate = values[k];<br><br>                    s = s +std::to_string(coordinate)+</span></font><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">" "</span><font face="Menlo"><span style="font-size: 9pt;">;<br>                }<br><br>                s = s + (isGhost(cdm, vertexIndices[j])?</span></font><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">" ghost"</span><font face="Menlo"><span style="font-size: 9pt;">:</span></font><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">""</span><font face="Menlo"><span style="font-size: 9pt;">);<br><br>                s = s + </span></font><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">"</span><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">\n</span><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">"</span><font face="Menlo"><span style="font-size: 9pt;">;<br><br>            }<br>        }<br>    }<br><br>    VecRestoreArray(coordinates, &coords);<br><br>    </span></font><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">int </span><font face="Menlo"><span style="font-size: 9pt;">rank;<br>    MPI_Comm_rank (PETSC_COMM_WORLD, &rank);   </span></font><span style="color: rgb(128, 128, 128); font-family: Menlo; font-size: 9pt; font-style: italic;">// get current process id<br></span><span style="color: rgb(128, 128, 128); font-family: Menlo; font-size: 9pt; font-style: italic;"><br></span><span style="color: rgb(128, 128, 128); font-family: Menlo; font-size: 9pt; font-style: italic;">    </span><font face="Menlo"><span style="font-size: 9pt;">PetscSynchronizedPrintf(PETSC_COMM_WORLD,</span></font><font color="#008000" face="Menlo"><span style="font-size: 12px;"><b>”</b></span></font><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">dmplex iteration rank %d </span><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">\n</span><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;"> %s</span><span style="color: rgb(0, 0, 128); font-family: Menlo; font-size: 9pt; font-weight: bold;">\n</span><span style="color: rgb(0, 128, 0); font-family: Menlo; font-size: 9pt; font-weight: bold;">"</span><font face="Menlo"><span style="font-size: 9pt;">,rank, s.c_str());<br>    PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);<br>}</span></font></pre>
</div>
<div style="color: rgb(0, 0, 0); font-family: Calibri, sans-serif; font-size: 14px;">
<div id="MAC_OUTLOOK_SIGNATURE"></div>
</div>
</div>
<div style="color: rgb(0, 0, 0); font-family: Calibri, sans-serif; font-size: 14px;">
<br>
</div>
<span id="OLK_SRC_BODY_SECTION" style="color: rgb(0, 0, 0); font-family: Calibri, sans-serif; font-size: 14px;">
<div style="font-family:Calibri; font-size:12pt; text-align:left; color:black; BORDER-BOTTOM: medium none; BORDER-LEFT: medium none; PADDING-BOTTOM: 0in; PADDING-LEFT: 0in; PADDING-RIGHT: 0in; BORDER-TOP: #b5c4df 1pt solid; BORDER-RIGHT: medium none; PADDING-TOP: 3pt">
<span style="font-weight:bold">From: </span>Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>><br>
<span style="font-weight:bold">Date: </span>Monday 30 November 2015 at 14:08 <br>
<span style="font-weight:bold">To: </span>Morten Nobel-Jørgensen <<a href="mailto:mono@mek.dtu.dk">mono@mek.dtu.dk</a>><br>
<span style="font-weight:bold">Cc: </span>"<a href="mailto:petsc-users@mcs.anl.gov">petsc-users@mcs.anl.gov</a>" <<a href="mailto:petsc-users@mcs.anl.gov">petsc-users@mcs.anl.gov</a>><br>
<span style="font-weight:bold">Subject: </span>Re: [petsc-users] DMPlex: Ghost points after DMRefine<br>
</div>
<div><br>
</div>
<div>
<div>
<div dir="ltr">
<div class="gmail_extra">
<div class="gmail_quote">On Mon, Nov 30, 2015 at 7:01 AM, Morten Nobel-Jørgensen <span dir="ltr">
<<a href="mailto:mono@mek.dtu.dk" target="_blank">mono@mek.dtu.dk</a>></span> wrote:<br>
<blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
<div style="word-wrap:break-word;color:rgb(0,0,0);font-size:14px">
<div style="font-family:Calibri,sans-serif">I have a very simple unstructured mesh composed of two triangles (four vertices) with one shared edge using a DMPlex:</div>
<div><font face="Courier"><br>
</font></div>
<div><font face="Courier"> /|\</font></div>
<div><font face="Courier">/ | \</font></div>
<div><font face="Courier">\ | /</font></div>
<div><font face="Courier"> \|/</font></div>
<div style="font-family:Calibri,sans-serif"><br>
</div>
<div style="font-family:Calibri,sans-serif">After distributing this mesh to two processes, each process owns a triangle. However one process owns tree vertices, while the last vertex is owned by the other process.</div>
<div style="font-family:Calibri,sans-serif"><br>
</div>
<div style="font-family:Calibri,sans-serif">The problem occurs after uniformly refining the dm. The mesh now looks like this: </div>
<div style="font-family:Calibri,sans-serif"><br>
</div>
<div style="font-family:Calibri,sans-serif">
<div style="font-family:-webkit-standard"><font face="Courier"> /|\</font></div>
<div style="font-family:-webkit-standard"><font face="Courier">/\|/\</font></div>
<div style="font-family:-webkit-standard"><font face="Courier">\/|\/</font></div>
<div style="font-family:-webkit-standard"><font face="Courier"> \|/</font></div>
<div style="font-family:-webkit-standard"><font face="Courier"><br>
</font></div>
<div style="font-family:-webkit-standard"><span style="font-family:Calibri,sans-serif">The new center vertex is now not listed as a ghost vertex but instead exists as two individual points.</span></div>
<div style="font-family:-webkit-standard"><span style="font-family:Calibri,sans-serif"><br>
</span></div>
<div style="font-family:-webkit-standard"><span style="font-family:Calibri,sans-serif">Is there any way that this new center vertex could be created as a ghost vertex during refinement?</span></div>
</div>
</div>
</blockquote>
<div><br>
</div>
<div>This could be a bug with the l2g mapping. I do not recreate it when refining, only the SF defining the mapping.</div>
<div>Here is an experiment: do not retrieve the mapping until after the refinement. Do you get what you want? If so,</div>
<div>I can easily fix this by destroying the map when I refine.</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 style="word-wrap:break-word;color:rgb(0,0,0);font-size:14px">
<div style="font-family:Calibri,sans-serif">
<div style="font-family:-webkit-standard"><span style="font-family:Calibri,sans-serif">Kind regards,</span></div>
<div style="font-family:-webkit-standard"><span style="font-family:Calibri,sans-serif">Morten</span></div>
<div style="font-family:-webkit-standard"><span style="font-family:Calibri,sans-serif"><br>
</span></div>
<div style="font-family:-webkit-standard"><span style="font-family:Calibri,sans-serif">Ps. Here are some code snippets for getting global point index and test of point is a ghost point:</span></div>
<div style="font-family:-webkit-standard"><span style="font-family:Calibri,sans-serif"><br>
</span></div>
<div style="font-family:-webkit-standard">
<pre style="background-color:rgb(255,255,255);font-family:Menlo;font-size:9pt"><span style="color:#000080;font-weight:bold">int </span>localToGlobal(DM dm, PetscInt point){<br>    <span style="color:#000080;font-weight:bold">const </span>PetscInt* array;<br>    ISLocalToGlobalMapping ltogm;<br>    DMGetLocalToGlobalMapping(dm,&ltogm);<br>    ISLocalToGlobalMappingGetIndices(ltogm, &array);<br>    PetscInt res = array[point];<br>    <span style="color:#000080;font-weight:bold">if </span>(res < <span style="color:#0000ff">0</span>){ <span style="color:#808080;font-style:italic">// if ghost<br></span><span style="color:#808080;font-style:italic">        </span>res = -res +<span style="color:#0000ff">1</span>;<br>    }<br>    <span style="color:#000080;font-weight:bold">return </span>res;<br>}<br><br><span style="color:#000080;font-weight:bold">bool </span>isGhost(DM dm, PetscInt point){<br>    <span style="color:#000080;font-weight:bold">const </span>PetscInt* array;<br>    ISLocalToGlobalMapping ltogm;<br>    DMGetLocalToGlobalMapping(dm,&ltogm);<br>    ISLocalToGlobalMappingGetIndices(ltogm, &array);<br>    <span style="color:#000080;font-weight:bold">return </span>array[point]<<span style="color:#0000ff">0</span>;<br>}</pre>
</div>
</div>
<div style="font-family:Calibri,sans-serif">
<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>
</div>
</div>
</span>
</body>
</html>