<div dir="ltr">Hi Matt<br><br>I think I've gotten it just about there. I'm just having an issue with the VecISCopy. I have an IS built that matches size correctly to map from the full state to the filtered state. The core issue I think, is should the expanded IS the ownership range of the vector subtracted out. Looking at the implementation, it looks like VecISCopy takes care of that for me. (Line 573 in src/vec/vec/utils/projection.c) But I could be mistaken. <br><br>The error I am getting is: <br>[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br>[0]PETSC ERROR: No support for this operation for this object type<br>[0]PETSC ERROR: Only owned values supported<br><br><br>Here is what I am currently doing.<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(235,187,255)">call</span> DMPlexFilter(dmplex_full, iBlankLabel, <span style="color:rgb(255,197,143)">1</span>, dmplex_filtered,ierr)</div><div><span style="color:rgb(235,187,255)">call</span> DMPlexGetSubpointIS(dmplex_filtered, subpointsIS,ierr)</div><br><div><span style="color:rgb(114,133,183)">! adds section to dmplex_filtered and allocates vec_filtered using DMCreateGlobalVector</span></div><div><span style="color:rgb(235,187,255)">call</span> addSectionToDMPlex(dmplex_filtered,vec_filtered) </div><br><div><span style="color:rgb(114,133,183)">! Get Sections for dmplex_filtered and dmplex_full</span></div><div><span style="color:rgb(235,187,255)">call</span> DMGetGlobalSection(dmplex_filtered,filteredfieldSection,ierr)</div><div><span style="color:rgb(235,187,255)">call</span> DMGetGlobalSection(dmplex_full,fullfieldSection,ierr)</div><br><br><br><div><span style="color:rgb(235,187,255)">call</span> ISGetIndicesF90(subpointsIS, subPointKey,ierr)</div><div>ExpandedIndexSize <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(255,197,143)">0</span></div><div><span style="color:rgb(235,187,255)">do</span> i <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(255,197,143)">1</span>, <span style="color:rgb(235,187,255)">size</span>(subPointKey)</div><div> <span style="color:rgb(235,187,255)">call</span> PetscSectionGetDof(fullfieldSection, subPointKey(i), dof,ierr)</div><div> ExpandedIndexSize <span style="color:rgb(153,255,255)">=</span> ExpandedIndexSize <span style="color:rgb(153,255,255)">+</span> dof</div><div><span style="color:rgb(235,187,255)">enddo</span></div><br><br><div><span style="color:rgb(114,133,183)">!Create expandedIS from offset sections of full and filtered sections</span></div><div><span style="color:rgb(235,187,255)">allocate</span>(ExpandedIndex(ExpandedIndexSize))</div><div><span style="color:rgb(235,187,255)">call</span> VecGetOwnershipRange(vec_full,oStart,oEnd,ierr)</div><div><span style="color:rgb(235,187,255)">do</span> i <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(255,197,143)">1</span>, <span style="color:rgb(235,187,255)">size</span>(subPointKey)</div><div> <span style="color:rgb(235,187,255)">call</span> PetscSectionGetOffset(fullfieldSection, subPointKey(i), offset,ierr)</div><div> <span style="color:rgb(235,187,255)">call</span> PetscSectionGetDof(fullfieldSection, subPointKey(i), dof,ierr)</div><div> <span style="color:rgb(114,133,183)">!offset=offset-oStart !looking at VecIScopy it takes care of this subtraction (not sure)</span></div><div> <span style="color:rgb(235,187,255)">do</span> j <span style="color:rgb(153,255,255)">=</span> <span style="color:rgb(255,197,143)">1</span>, (dof)</div><div> ExpandedIndex((i<span style="color:rgb(255,197,143)">-1</span>)<span style="color:rgb(153,255,255)">*</span>dof<span style="color:rgb(153,255,255)">+</span>j) <span style="color:rgb(153,255,255)">=</span> offset<span style="color:rgb(153,255,255)">+</span>j</div><div> <span style="color:rgb(235,187,255)">end do</span></div><div><span style="color:rgb(235,187,255)">enddo</span></div><br><div><span style="color:rgb(235,187,255)">call</span> ISCreateGeneral(PETSC_COMM_WORLD, ExpandedIndexSize, ExpandedIndex, PETSC_COPY_VALUES, expandedIS,ierr)</div><div><span style="color:rgb(235,187,255)">call</span> ISRestoreIndicesF90(subpointsIS, subPointKey,ierr)</div><div><span style="color:rgb(235,187,255)">deallocate</span>(ExpandedIndex)</div><br><br><div><span style="color:rgb(235,187,255)">call</span> VecGetLocalSize(vec_full,sizeVec,ierr)</div><div><span style="color:rgb(235,187,255)">write</span>(<span style="color:rgb(153,255,255)">*</span>,<span style="color:rgb(153,255,255)">*</span>) sizeVec</div><div><span style="color:rgb(235,187,255)">call</span> VecGetLocalSize(vec_filtered,sizeVec,ierr)</div><div><span style="color:rgb(235,187,255)">write</span>(<span style="color:rgb(153,255,255)">*</span>,<span style="color:rgb(153,255,255)">*</span>) sizeVec</div><div><span style="color:rgb(235,187,255)">call</span> ISGetLocalSize(expandedIS,sizeVec,ierr)</div><div><span style="color:rgb(235,187,255)">write</span>(<span style="color:rgb(153,255,255)">*</span>,<span style="color:rgb(153,255,255)">*</span>) sizeVec</div><div><span style="color:rgb(235,187,255)">call</span> PetscSynchronizedFlush(PETSC_COMM_WORLD,ierr)</div><br><br><div><span style="color:rgb(235,187,255)">call</span> VecISCopy(vec_full,expandedIS,SCATTER_REVERSE,vec_filtered,ierr)</div></div><br><br>Thanks again for the great help.<br><br>Sincerely<br>Nicholas<br><br></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Wed, Dec 7, 2022 at 9:29 PM 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 Wed, Dec 7, 2022 at 9:21 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">Thank you for the help.<br><br>I think the last piece of the puzzle is how do I create the "expanded IS" from the subpoint IS using the section? <br></div></blockquote><div><br></div><div>Loop over the points in the IS. For each point, get the dof and offset from the Section. Make a new IS that has all the</div><div>dogs, namely each run [offset, offset+dof).</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">Sincerely<br>Nicholas<br><br></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Wed, Dec 7, 2022 at 7:06 AM Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">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 Wed, Dec 7, 2022 at 6:51 AM 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">Hi <br><br>Thank you so much for your patience. One thing to note: I don't have any need to go back from the filtered distributed mapping back to the full but it is good to know. <br><br>One aside question.<br>1) Is natural and global ordering the same in this context? <br></div></blockquote><div><br></div><div>No.</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">As far as implementing what you have described.<br><br>When I call ISView on the generated SubpointIS, I get an unusual error which I'm not sure how to interpret. (this case is running on 2 ranks and the filter label has points located on both ranks of the original DM. However, if I manually get the indices (the commented lines), it seems to not have any issues. <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> DMPlexFilter(dmplex_full, iBlankLabel, <span style="color:rgb(255,197,143)">1</span>, dmplex_filtered,ierr)</div><div><span style="color:rgb(235,187,255)">call</span> DMPlexGetSubpointIS(dmplex_filtered, subpointsIS,ierr)</div><div><span style="color:rgb(114,133,183)">!call ISGetIndicesF90(subpointsIS, subPointKey,ierr)</span></div><div><span style="color:rgb(114,133,183)">!write(*,*) subPointKey</span></div><div><span style="color:rgb(114,133,183)">!call ISRestoreIndicesF90(subpointsIS, subPointKey,ierr)</span></div><div><span style="color:rgb(235,187,255)">call</span> ISView(subpointsIS,PETSC_VIEWER_STDOUT_WORLD,ierr)</div></div><br>[1]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br>[1]PETSC ERROR: Arguments must have same communicators<br>[1]PETSC ERROR: Different communicators in the two objects: Argument # 1 and 2 flag 3<br>[1]PETSC ERROR: See <a href="https://petsc.org/release/faq/" target="_blank">https://petsc.org/release/faq/</a> for trouble shooting.<br>[1]PETSC ERROR: Petsc Development GIT revision: v3.18.1-320-g7810d690132 GIT Date: 2022-11-20 20:25:41 -0600<br>[1]PETSC ERROR: Configure options with-fc=mpiifort with-mpi-f90=mpiifort --download-triangle --download-parmetis --download-metis --with-debugging=1 --download-hdf5 --prefix=/home/narnoldm/packages/petsc_install<br>[1]PETSC ERROR: #1 ISView() at /home/narnoldm/packages/petsc/src/vec/is/is/interface/index.c:1629<br></div></blockquote><div><br></div><div>The problem here is the subpointsIS is a _serial_ object, and you are using a parallel viewer. You can use PETSC_VIEWER_STDOUT_SELF,</div><div>or you can pull out the singleton viewer from STDOUT_WORLD if you want them all to print in order.</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">As far as the overall process you have described my question on first glance is do I have to allocate/create the vector that is output by VecISCopy before calling it, or does it create the vector automatically?</div></blockquote><div><br></div><div>You create both vectors. I would do it using DMCreateGlobalVector() from both DMs.</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"> I think I would need to create it first using a section and Setting the Vec in the filtered DM?</div></blockquote><div><br></div><div>Setting the Section in the filtered 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"> And I presume in this case I would be using the scatter reverse option to go from the full set to the reduced set? <br></div></blockquote><div><br></div><div>Yes</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">Sincerely<br>Nicholas<br><br><br><br><br><br><br>Sincerely <br>Nick</div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Wed, Dec 7, 2022 at 6:00 AM Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">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 Wed, Dec 7, 2022 at 3:35 AM 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">Hi Matthew<br><br>Thank you for the help. This clarified a great deal.<br><br>I have a follow-up question related to DMPlexFilter. It may be better to describe what I'm trying to achieve. <br><br>I have a general mesh I am solving which has a section with cell center finite volume states, as described in my initial email. After calculating some metrics, I tag a bunch of cells with an identifying Label and use DMFilter to generate a new DM which is only that subset of cells. Generally, this leads to a pretty unbalanced DM so I then plan to use DMPlexDIstribute to balance that DM across the processors. The coordinates pass along fine, but the state(or I should say Section) does not at least as far as I can tell.<br> <br>Assuming I can get a filtered DM I then distribute the DM and state using the method you described above and it seems to be working ok. <br><br>The last connection I have to make is the transfer of information from the full mesh to the "sampled" filtered mesh. From what I can gather I would need to get the mapping of points using DMPlexGetSubpointIS and then manually copy the values from the full DM section to the filtered DM? I have the process from full->filtered->distributed all working for the coordinates so its just a matter of transferring the section correctly. <br><br>I appreciate all the help you have provided. <br></div></blockquote><div><br></div><div>Let's do this in two steps, which makes it easier to debug. First, do not redistribute the submesh. Just use DMPlexGetSubpointIS()</div><div>to get the mapping of filtered points to points in the original mesh. Then create an expanded IS using the Section which makes</div><div>dofs in the filtered mesh to dofs in the original mesh. >From this use</div><div><br></div><div> <a href="https://petsc.org/main/docs/manualpages/Vec/VecISCopy/" target="_blank">https://petsc.org/main/docs/manualpages/Vec/VecISCopy/</a></div><div><br></div><div>to move values between the original vector and the filtered vector.</div><div><br></div><div>Once that works, you can try redistributing the filtered mesh. Before calling DMPlexDistribute() on the filtered mesh, you need to call</div><div><br></div><div> <a href="https://petsc.org/main/docs/manualpages/DM/DMSetUseNatural/" target="_blank">https://petsc.org/main/docs/manualpages/DM/DMSetUseNatural/</a></div><div><br></div><div>When you redistribute, it will compute a mapping back to the original layout. Now when you want to transfer values, you</div><div><br></div><div> 1) Create a natural vector with DMCreateNaturalVec()</div><div><br></div><div> 2) Use DMGlobalToNaturalBegin/End() to move values from the filtered vector to the natural vector</div><div><br></div><div> 3) Use VecISCopy() to move values from the natural vector to the original vector</div><div><br></div><div>Let me know if you have any problems.</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">Sincerely<br>Nicholas<br><br><br></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Mon, Nov 28, 2022 at 6:19 AM Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">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 Sun, Nov 27, 2022 at 10:22 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">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"><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"><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>
</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"><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>
</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"><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>
</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>