<html>
  <head>
    <meta http-equiv="Content-Type" content="text/html; charset=utf-8">
  </head>
  <body text="#000000" bgcolor="#FFFFFF">
    On 2018-04-27 04:11 AM, Matthew Knepley wrote:<br>
    <blockquote type="cite"
cite="mid:CAMYG4G=nX0OZ+C6pO7M-FVvVvEP1T2yGH47v5jwtyBrv47fkyg@mail.gmail.com">
      <div dir="ltr">
        <div class="gmail_extra">
          <div class="gmail_quote">On Fri, Apr 27, 2018 at 2:09 AM,
            Danyang Su <span dir="ltr"><<a
                href="mailto:danyang.su@gmail.com" target="_blank"
                moz-do-not-send="true">danyang.su@gmail.com</a>></span>
            wrote:<br>
            <blockquote class="gmail_quote" style="margin:0px 0px 0px
              0.8ex;border-left:1px solid
              rgb(204,204,204);padding-left:1ex">
              <div bgcolor="#FFFFFF">
                <p>Hi Matt,</p>
                <p>Sorry if this is a stupid question. <br>
                </p>
                <p>In the previous code for unstructured grid, I create
                  labels to mark the original node/cell index from VTK
                  file and then distribute it so that each subdomain has
                  a copy of its original node and cell index, as well as
                  the PETSc numbering. Now I am trying to get avoid of
                  using large number of keys in DMSetLabelValue since
                  this costs lot of time for large problem. <br>
                </p>
                <p>I can get the coordinates of subdomain after
                  distribution by using DMGetCoordinatesLocal and
                  DMGetCoordinateDM.</p>
                <p>How can I get the vertex index of each cell after
                  distribution? Would you please give me a hint or
                  functions that I can use.<br>
                </p>
              </div>
            </blockquote>
            <div>You can permute the vectors back to the natural
              ordering using</div>
            <div><br>
            </div>
            <div>  <a
href="http://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/DMPLEX/DMPlexNaturalToGlobalBegin.html"
                moz-do-not-send="true">http://www.mcs.anl.gov/petsc/petsc-master/docs/manualpages/DMPLEX/DMPlexNaturalToGlobalBegin.html</a></div>
            <div><br>
            </div>
            <div>which says you have to call DMPlexSetUseNaturalSF()
              before distributing the mesh. It is tested in</div>
            <div><br>
            </div>
            <div>  src/dm/impls/plex/examples/tests/ex15.c</div>
            <div><br>
            </div>
            <div>so you can see how its intended to work. It is very new
              and has not been tested by many people.</div>
            <div><br>
            </div>
            <div>I can see how you might want this for small tests. Why
              would you want it for production models?</div>
          </div>
        </div>
      </div>
    </blockquote>
    Hi Matt,<br>
    <br>
    This is indeed what I need. As some of years old cases import
    initial conditions from external files, which are in natural
    ordering as the original mesh. Just want to make the code compatible
    to the old input files.<br>
    <br>
    Thanks,<br>
    <br>
    Danyang<br>
    <blockquote type="cite"
cite="mid:CAMYG4G=nX0OZ+C6pO7M-FVvVvEP1T2yGH47v5jwtyBrv47fkyg@mail.gmail.com">
      <div dir="ltr">
        <div class="gmail_extra">
          <div class="gmail_quote">
            <div><br>
            </div>
            <div>  Thanks,</div>
            <div><br>
            </div>
            <div>    Matt <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 bgcolor="#FFFFFF">
                <p> </p>
                Thanks,<br>
                <p>Danyang<br>
                </p>
                <br>
                <div class="gmail-m_8572319983782213150moz-cite-prefix">On
                  18-04-25 02:12 PM, Danyang Su wrote:<br>
                </div>
                <blockquote type="cite"> On 2018-04-25 09:47 AM, Matthew
                  Knepley wrote:<br>
                  <blockquote type="cite">
                    <div dir="ltr">
                      <div class="gmail_extra">
                        <div class="gmail_quote">On Wed, Apr 25, 2018 at
                          12:40 PM, Danyang Su <span dir="ltr"><<a
                              href="mailto:danyang.su@gmail.com"
                              target="_blank" moz-do-not-send="true">danyang.su@gmail.com</a>></span>
                          wrote:<br>
                          <blockquote class="gmail_quote"
                            style="margin:0px 0px 0px
                            0.8ex;border-left:1px solid
                            rgb(204,204,204);padding-left:1ex">
                            <div bgcolor="#FFFFFF">
                              <p>Hi Matthew,</p>
                              <p>In the worst case, every node/cell may
                                have different label. <br>
                              </p>
                            </div>
                          </blockquote>
                          <div>Do not use Label for this. Its not an
                            appropriate thing. If every cell is
                            different, just use the cell number.</div>
                          <div>Labels are for mapping a relatively small
                            number of keys (like material IDs) to sets
                            of points (cells, vertices, etc.) </div>
                          <div>Its not a great data structure for a
                            permutation.</div>
                        </div>
                      </div>
                    </div>
                  </blockquote>
                  Yes. If there is small number of keys, it runs very
                  fast, even for more than one million DMSetLabelValue
                  calls. The performance just deteriorates as the number
                  of keys increases. <br>
                  <br>
                  I cannot get avoid of DMSetLabelValue as node/cell
                  index of original mesh is needed for the previous
                  input file that uses some of global node/cell index to
                  set value. But if I can get the natural order of
                  nodes/cells from DMPlex, I can discard the use of
                  DMSetLabelValue. Is there any function can do this
                  job? <br>
                  <br>
                  Thanks,<br>
                  <br>
                  Danyang<br>
                  <blockquote type="cite">
                    <div dir="ltr">
                      <div class="gmail_extra">
                        <div class="gmail_quote">
                          <div><br>
                          </div>
                          <div>However, I still do not believe these
                            numbers. The old code does a string
                            comparison every time. I will setup a test.</div>
                          <div><br>
                          </div>
                          <div>   Matt</div>
                          <blockquote class="gmail_quote"
                            style="margin:0px 0px 0px
                            0.8ex;border-left:1px solid
                            rgb(204,204,204);padding-left:1ex">
                            <div bgcolor="#FFFFFF">
                              <p>Below is one of the worst scenario with
                                102299 nodes and 102299 different labels
                                for test. I found the time cost increase
                                during the loop. The first 9300 loop
                                takes least time (<0.5) while the
                                last 9300 loops takes much more time
                                (>7.7), as shown below. If I use
                                larger mesh with >1 million nodes, it
                                runs very very slowly in this part. The
                                PETSc is configured with optimization
                                on.</p>
                              <p>Configure options --with-cc=gcc
                                --with-cxx=g++ --with-fc=gfortran
                                --download-mpich --download-scalapack
                                --download-parmetis --download-metis
                                --download-ptscotch
                                --download-fblaslapack --download-hypre
                                --download-superlu_dist
                                --download-hdf5=yes --download-ctetgen
                                --with-debugging=0 COPTFLAGS="-O3
                                -march=native -mtune=native"
                                CXXOPTFLAGS="-O3 -march=native
                                -mtune=native" FOPTFLAGS="-O3
                                -march=native -mtune=native"</p>
                              <p> </p>
                              <table cellspacing="0" width="816"
                                border="0" height="245">
                                <colgroup span="6" width="85"></colgroup>
                                <tbody>
                                  <tr>
                                    <td align="left" height="32">istart</td>
                                    <td align="left">iend</td>
                                    <td align="left">progress</td>
                                    <td align="left">CPU_Time</td>
                                    <td align="left">time cost - old
                                      (sec) </td>
                                    <td align="left">time cost - new
                                      (sec)</td>
                                  </tr>
                                  <tr>
                                    <td align="right" height="17">0</td>
                                    <td align="right">9299</td>
                                    <td align="right">0</td>
                                    <td align="right">1524670045.51166</td>
                                    <td align="left"><br>
                                    </td>
                                    <td align="left"><br>
                                    </td>
                                  </tr>
                                  <tr>
                                    <td align="right" height="17">9300</td>
                                    <td align="right">18599</td>
                                    <td align="right">0.100010753</td>
                                    <td align="right">1524670045.99605</td>
                                    <td align="right">0.4843890667</td>
                                    <td align="right">0.497246027</td>
                                  </tr>
                                  <tr>
                                    <td align="right" height="17">18600</td>
                                    <td align="right">27899</td>
                                    <td align="right">0.200010747</td>
                                    <td align="right">1524670047.32635</td>
                                    <td align="right">1.330302</td>
                                    <td align="right">1.3820912838</td>
                                  </tr>
                                  <tr>
                                    <td align="right" height="17">27900</td>
                                    <td align="right">37199</td>
                                    <td align="right">0.300010741</td>
                                    <td align="right">1524670049.3066</td>
                                    <td align="right">1.9802515507</td>
                                    <td align="right">2.2439446449</td>
                                  </tr>
                                  <tr>
                                    <td align="right" height="17">37200</td>
                                    <td align="right">46499</td>
                                    <td align="right">0.400010765</td>
                                    <td align="right">1524670052.1594</td>
                                    <td align="right">2.852804184</td>
                                    <td align="right">3.0739262104</td>
                                  </tr>
                                  <tr>
                                    <td align="right" height="17">46500</td>
                                    <td align="right">55799</td>
                                    <td align="right">0.500010729</td>
                                    <td align="right">1524670055.90961</td>
                                    <td align="right">3.7502081394</td>
                                    <td align="right">3.9270553589</td>
                                  </tr>
                                  <tr>
                                    <td align="right" height="17">55800</td>
                                    <td align="right">65099</td>
                                    <td align="right">0.600010753</td>
                                    <td align="right">1524670060.47654</td>
                                    <td align="right">4.5669286251</td>
                                    <td align="right">4.7571902275</td>
                                  </tr>
                                  <tr>
                                    <td align="right" height="17">65100</td>
                                    <td align="right">74399</td>
                                    <td align="right">0.700010777</td>
                                    <td align="right">1524670066.0941</td>
                                    <td align="right">5.6175630093</td>
                                    <td align="right">5.7428796291</td>
                                  </tr>
                                  <tr>
                                    <td align="right" height="17">74400</td>
                                    <td align="right">83699</td>
                                    <td align="right">0.800010741</td>
                                    <td align="right">1524670072.53886</td>
                                    <td align="right">6.44475317</td>
                                    <td align="right">6.5761549473</td>
                                  </tr>
                                  <tr>
                                    <td align="right" height="17">83700</td>
                                    <td align="right">92998</td>
                                    <td align="right">0.900010765</td>
                                    <td align="right">1524670079.99072</td>
                                    <td align="right">7.4518604279</td>
                                    <td align="right">7.4606924057</td>
                                  </tr>
                                  <tr>
                                    <td align="right" height="17">92999</td>
                                    <td align="right">102298</td>
                                    <td align="right">1</td>
                                    <td align="right">1524670087.71066</td>
                                    <td align="right">7.7199423313</td>
                                    <td align="right">8.2424075603</td>
                                  </tr>
                                </tbody>
                              </table>
                              <br>
                              <br>
                              old code<br>
                              <br>
                                      do ipoint = 0, istart-1<br>
                                        !c output time cost, use 1
                              processor to test<br>
                                        if (b_enable_output .and. rank
                              == 0) then<br>
                                          if (mod(ipoint,iprogress) == 0
                              .or. ipoint == istart-1) then<br>
                                           
                              !write(*,'(f3.1,1x)',advance="<wbr>no")
                              (ipoint+1.0)/istart<br>
                                            write(*,*) ipoint,
                              (ipoint+1.0)/istart,"time",MPI<wbr>_Wtime()<br>
                                          end if<br>
                                        end if<br>
                              <br>
                                        call
                              DMSetLabelValue(dmda_flow%da,"<wbr>cid_lg2g",ipoint,        
                              &<br>
                                                            
                              ipoint+1,ierr)<br>
                                        CHKERRQ(ierr)<br>
                                      end do<br>
                              <br>
                              <br>
                              new code<br>
                              <br>
                                      call
                              DMCreateLabel(dmda_flow%da,'ci<wbr>d_lg2g',ierr)<br>
                                      CHKERRQ(ierr)<br>
                              <br>
                                      call
                              DMGetLabel(dmda_flow%da,'cid_l<wbr>g2g',label,
                              ierr)<br>
                                      CHKERRQ(ierr)<br>
                              <br>
                                      do ipoint = 0, istart-1<br>
                                        !c output time cost, use 1
                              processor to test<br>
                                        if (b_enable_output .and. rank
                              == 0) then<br>
                                          if (mod(ipoint,iprogress) == 0
                              .or. ipoint == istart-1) then<br>
                                           
                              !write(*,'(f3.1,1x)',advance="<wbr>no")
                              (ipoint+1.0)/istart<br>
                                            write(*,*) ipoint,
                              (ipoint+1.0)/istart,"time",MPI<wbr>_Wtime()<br>
                                          end if<br>
                                        end if<br>
                              <br>
                                        call
                              DMLabelSetValue(label,ipoint,i<wbr>point+1,ierr)<br>
                                        CHKERRQ(ierr)<br>
                                      end do<br>
                              <br>
                              Thanks,<br>
                              <br>
                              Danyang<br>
                              <br>
                              <div
                                class="gmail-m_8572319983782213150m_-7982693774986852993moz-cite-prefix">On
                                2018-04-25 03:16 AM, Matthew Knepley
                                wrote:<br>
                              </div>
                              <blockquote type="cite">
                                <div dir="ltr">
                                  <div class="gmail_extra">
                                    <div class="gmail_quote">On Tue, Apr
                                      24, 2018 at 11:57 PM, Danyang Su <span
                                        dir="ltr"><<a
                                          href="mailto:danyang.su@gmail.com"
                                          target="_blank"
                                          moz-do-not-send="true">danyang.su@gmail.com</a>></span>
                                      wrote:<br>
                                      <blockquote class="gmail_quote"
                                        style="margin:0px 0px 0px
                                        0.8ex;border-left:1px solid
                                        rgb(204,204,204);padding-left:1ex">Hi
                                        All,<br>
                                        <br>
                                        I use DMPlex in unstructured
                                        grid code and recently found
                                        DMSetLabelValue takes a lot of
                                        time for large problem, e.g.,
                                        num. of cells > 1 million. In
                                        my code, I use<br>
                                      </blockquote>
                                      <div><br>
                                      </div>
                                      <div>I read your code wrong. For
                                        large loop, you should not use
                                        the convenience function. You
                                        should use</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">
                                        DMPlexCreateFromCellList ()</blockquote>
                                      <div><br>
                                      </div>
                                      <div>DMGetLabel(dm, name,
                                        &label)</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"> </blockquote>
                                      <blockquote class="gmail_quote"
                                        style="margin:0px 0px 0px
                                        0.8ex;border-left:1px solid
                                        rgb(204,204,204);padding-left:1ex">
                                        <br>
                                        Loop over all cells/nodes{<br>
                                        <br>
                                        DMSetLabelValue<br>
                                      </blockquote>
                                      <div><br>
                                      </div>
                                      <div>Replace this by
                                        DMLabelSetValue(label, point,
                                        val)</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">
                                        }<br>
                                        <br>
                                        DMPlexDistribute<br>
                                        <br>
                                        The code works fine except
                                        DMSetLabelValue takes a lot of
                                        time for large problem. I use
                                        DMSetLabelValue to set material
                                        id for all the nodes or cells so
                                        that each subdomain has a copy
                                        of material id. Is there any
                                        other functions that can be used
                                        more efficient, e.g. set labels
                                        by array, not 1 by 1?<br>
                                      </blockquote>
                                      <div><br>
                                      </div>
                                      <div>That should take much less
                                        time.</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">
                                        Thanks,<br>
                                        <br>
                                        Danyang<br>
                                        <br>
                                      </blockquote>
                                    </div>
                                    <br>
                                    <br clear="all">
                                    <span class="gmail-HOEnZb"><font
                                        color="#888888"> <span
                                          class="gmail-m_8572319983782213150HOEnZb"><font
                                            color="#888888">
                                            <div><br>
                                            </div>
                                            -- <br>
                                            <div
                                              class="gmail-m_8572319983782213150m_-7982693774986852993gmail_signature">
                                              <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.caam.rice.edu/%7Emk51/"
                                                        target="_blank"
moz-do-not-send="true">https://www.cse.buffalo.edu/~k<wbr>nepley/</a><br>
                                                    </div>
                                                  </div>
                                                </div>
                                              </div>
                                            </div>
                                          </font></span></font></span></div>
                                  <span class="gmail-HOEnZb"><font
                                      color="#888888"> </font></span></div>
                                <span class="gmail-HOEnZb"><font
                                    color="#888888"> </font></span></blockquote>
                              <span class="gmail-HOEnZb"><font
                                  color="#888888"> <br>
                                </font></span></div>
                            <span class="gmail-HOEnZb"><font
                                color="#888888"> </font></span></blockquote>
                          <span class="gmail-HOEnZb"><font
                              color="#888888"> </font></span></div>
                        <span class="gmail-HOEnZb"><font color="#888888">
                            <br>
                            <br clear="all">
                            <div><br>
                            </div>
                            -- <br>
                            <div
                              class="gmail-m_8572319983782213150gmail_signature">
                              <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.caam.rice.edu/%7Emk51/"
                                        target="_blank"
                                        moz-do-not-send="true">https://www.cse.buffalo.edu/~<wbr>knepley/</a><br>
                                    </div>
                                  </div>
                                </div>
                              </div>
                            </div>
                          </font></span></div>
                      <span class="gmail-HOEnZb"><font color="#888888">
                        </font></span></div>
                    <span class="gmail-HOEnZb"><font color="#888888"> </font></span></blockquote>
                  <span class="gmail-HOEnZb"><font color="#888888"> <br>
                    </font></span></blockquote>
                <br>
              </div>
            </blockquote>
          </div>
          <br>
          <br clear="all">
          <div><br>
          </div>
          -- <br>
          <div class="gmail_signature">
            <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.caam.rice.edu/%7Emk51/"
                      target="_blank" moz-do-not-send="true">https://www.cse.buffalo.edu/~knepley/</a><br>
                  </div>
                </div>
              </div>
            </div>
          </div>
        </div>
      </div>
    </blockquote>
    <br>
  </body>
</html>