<div dir="ltr"><div dir="ltr">On Sun, Nov 27, 2022 at 10:22 PM Nicholas Arnold-Medabalimi <<a href="mailto:narnoldm@umich.edu">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">Hi Petsc Users<br><br>I have a question about properly using PetscSection to assign state variables to a DM. I have an existing DMPlex mesh distributed on 2 processors. My goal is to have state variables set to the cell centers. I then want to call DMPlexDistribute, which I hope will balance the mesh elements and hopefully transport the state variables to the hosting processors as the cells are distributed to a different processor count or simply just redistributing after doing mesh adaption.<br><br>Looking at the DMPlex User guide, I should be able to achieve this with a single field section using SetDof and assigning the DOF to the points corresponding to cells. <br></div></blockquote><div><br></div><div>Note that if you want several different fields, you can clone the DM first for this field</div><div><br></div><div> call DMClone(dm,dmState,ierr)</div><div><br></div><div>and use dmState in your calls below.</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"><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-wrap"><div><span style="color:rgb(235,187,255)">call</span> DMPlexGetHeightStratum(dm,<span style="color:rgb(255,197,143)">0</span>,c0,c1,ierr)</div><div><span style="color:rgb(235,187,255)">call</span> DMPlexGetChart(dm,p0,p1,ierr)</div><div><span style="color:rgb(235,187,255)">call</span> PetscSectionCreate(PETSC_COMM_WORLD,section,ierr)</div><div><span style="color:rgb(235,187,255)">call</span> PetscSectionSetNumFields(section,<span style="color:rgb(255,197,143)">1</span>,ierr)
<span style="color:rgb(235,187,255)">call</span> PetscSectionSetChart(section,p0,p1,ierr)</div><div><div style="line-height:19px"><div><span style="color:rgb(235,187,255)">do</span> i <span style="color:rgb(153,255,255)">=</span> c0, (c1<span style="color:rgb(255,197,143)">-1</span>)</div><div> <span style="color:rgb(235,187,255)">call</span> PetscSectionSetDof(section,i,nvar,ierr)</div><div><span style="color:rgb(235,187,255)">end do</span></div></div><div style="line-height:19px"><div><span style="color:rgb(235,187,255)">call</span> PetscSectionSetup(section,ierr)</div><div><span style="color:rgb(235,187,255)">call</span> DMSetLocalSection(dm,section,ierr)</div></div></div></div></div></blockquote><div><br></div><div>In the loop, I would add a call to</div><div><br></div><div> call PetscSectionSetFieldDof(section,i,0,nvar,ierr)</div><div><br></div><div>This also puts in the field breakdown. It is not essential, but nicer.</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"><div>From here, it looks like I can access and set the state vars using <br><br><div style="background-color:rgb(0,36,81);font-family:Consolas,"Courier New",monospace;font-size:14px;line-height:19px;white-space:pre-wrap"><div style="color:rgb(255,255,255)"><span style="color:rgb(235,187,255)">call</span> DMGetGlobalVector(dmplex,state,ierr) </div><div><span style="color:rgb(235,187,255)">call</span><font color="#ffffff"> DMGetGlobalSection(dmplex,section,ierr)
</font><div style="color:rgb(255,255,255);line-height:19px"><div><span style="color:rgb(235,187,255)">call</span> VecGetArrayF90(state,stateVec,ierr)</div></div><div><span style="color:rgb(235,187,255)">do</span><font color="#ffffff"> i </font><span style="color:rgb(153,255,255)">=</span><font color="#ffffff"> c0, (c1</font><span style="color:rgb(255,197,143)">-1</span><font color="#ffffff">)
</font><div style="color:rgb(255,255,255);line-height:19px"><div><span style="color:rgb(235,187,255)"> call</span> PetscSectionGetOffset(section,i,offset,ierr)</div></div><div style="line-height:19px"><div><font color="#ffffff"> stateVec(offset</font><font color="#99ffff">:(offset+nvar)</font><font color="#ffffff">)=state_i(:) !simplified assignment</font></div></div></div><div style="color:rgb(255,255,255)"><span style="color:rgb(235,187,255)">end do</span>
<div style="line-height:19px"><div><span style="color:rgb(235,187,255)">call</span> VecRestoreArrayF90(state,stateVec,ierr)</div></div><div style="line-height:19px"><div><span style="color:rgb(235,187,255)">call</span> DMRestoreGlobalVector(dmplex,state,ierr)</div></div></div></div></div></div><br>To my understanding, I should be using Global vector since this is a pure assignment operation and I don't need the ghost cells. <br></div></blockquote><div><br></div><div>Yes. </div><div><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">But the behavior I am seeing isn't exactly what I'd expect.<br><br>To be honest, I'm somewhat unclear on a few things<br><br>1) Should be using nvar fields with 1 DOF each or 1 field with nvar DOFs or what the distinction between the two methods are?<br></div></blockquote><div><br></div><div>We have two divisions in a Section. A field can have a number of components. This is intended to model a vector or tensor field.</div><div>Then a Section can have a number of fields, such as velocity and pressure for a Stokes problem. The division is mainly to help the</div><div>user, so I would use the most natural one.</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) Adding a print statement after the offset assignment I get (on rank 0 of 2) <br>cell 1 offset 0<br> cell 2 offset 18<br> cell 3 offset 36<br>which is expected and works but on rank 1 I get<br> cell 1 offset 9000<br> cell 2 offset 9018<br> cell 3 offset 9036<br><br>which isn't exactly what I would expect. Shouldn't the offsets reset at 0 for the next rank? <br></div></blockquote><div><br></div><div>The local and global sections hold different information. This is the source of the confusion. The local section does describe a local</div><div>vector, and thus includes overlap or "ghost" dofs. The global section describes a global vector. However, it is intended to deliver</div><div>global indices, and thus the offsets give back global indices. When you use VecGetArray*() you are getting out the local array, and</div><div>thus you have to subtract the first index on this process. You can get that from</div><div><br></div><div> VecGetOwnershipRange(v, &rstart, &rEnd);</div><div><br></div><div>This is the same whether you are using DMDA or DMPlex or any other DM.</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) Does calling DMPlexDistribute also distribute the section data associated with the DOF, based on the description in DMPlexDistribute it looks like it should?</div></blockquote><div><br></div><div>No. By default, DMPlexDistribute() only distributes coordinate data. I you want to distribute your field, it would look something like this:</div><div><br></div><div>DMPlexDistribute(dm, 0, &sfDist, &dmDist);</div>VecCreate(comm, &stateDist);</div><div class="gmail_quote">VecSetDM(sateDist, dmDist);<br>PetscSectionCreate(comm §ionDist);</div><div class="gmail_quote">DMSetLocalSection(dmDist, sectionDist);<br><div>DMPlexDistributeField(dmDist, sfDist, section, state, sectionDist, stateDist);<br></div><div><br></div><div>We do this in src/dm/impls/plex/tests/ex36.c</div><div><br></div><div> THanks,</div><div><br></div><div> Matt</div><div><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">I'd appreciate any insight into the specifics of this usage. I expect I have a misconception on the local vs global section. Thank you. <br><br>Sincerely<br>Nicholas<br><br>-- <br><div dir="ltr"><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></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_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>