<html>
  <head>
    <meta http-equiv="Content-Type" content="text/html; charset=UTF-8">
  </head>
  <body text="#000000" bgcolor="#FFFFFF">
    <p><br>
    </p>
    <div class="moz-cite-prefix">On 2019-10-10 5:44 p.m., Matthew
      Knepley wrote:<br>
    </div>
    <blockquote type="cite"
cite="mid:CAMYG4GkonhQpuJop6rTt664rs5TDMTmnPsh6ydq-zbgYb5R6CA@mail.gmail.com">
      <meta http-equiv="content-type" content="text/html; charset=UTF-8">
      <div dir="ltr">
        <div dir="ltr">On Thu, Oct 10, 2019 at 7:53 PM Danyang Su <<a
            href="mailto:danyang.su@gmail.com" moz-do-not-send="true">danyang.su@gmail.com</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 bgcolor="#FFFFFF">
              <p>On 2019-10-10 4:28 p.m., Matthew Knepley wrote:<br>
              </p>
              <blockquote type="cite">
                <div dir="ltr">
                  <div dir="ltr">On Thu, Oct 10, 2019 at 4:26 PM Danyang
                    Su <<a href="mailto:danyang.su@gmail.com"
                      target="_blank" moz-do-not-send="true">danyang.su@gmail.com</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 bgcolor="#FFFFFF">
                        <p>Hi All,</p>
                        <p>Your guess is right. The memory problem
                          occurs after DMPlexCreateFromCellList and
                          DMPlexDistribute. The mesh related memory in
                          the master processor is not released after
                          that. <br>
                        </p>
                        <p>The pseudo code I use is <br>
                        </p>
                        <p>if (rank == 0) then         !only the master
                          processor read the mesh file and create cell
                          list<br>
                        </p>
                        <p>        call
                          DMPlexCreateFromCellList(Petsc_Comm_World,ndim,num_cells,
                          &<br>
                                                               
                          num_nodes,num_nodes_per_cell,    &<br>
                                                               
                          Petsc_False,dmplex_cells,ndim,   &  !use
                          Petsc_True to create intermediate mesh
                          entities (faces, edges),<br>
                                                               
                          dmplex_verts,dmda_flow%da,ierr)     !not work
                          for prism for the current 3.8 version.<br>
                                  CHKERRQ(ierr)</p>
                        <p>else                                 !slave
                          processors pass zero cells<br>
                        </p>
                        <p>        call
                          DMPlexCreateFromCellList(Petsc_Comm_World,ndim,0,0,      
                          &<br>
                                                               
                          num_nodes_per_cell,              &<br>
                                                               
                          Petsc_False,dmplex_cells,ndim,   &  !use
                          Petsc_True to create intermediate mesh
                          entities (faces, edges),<br>
                                                               
                          dmplex_verts,dmda_flow%da,ierr)     !not work
                          for prism for the current 3.8 version.<br>
                                  CHKERRQ(ierr)<br>
                        </p>
                        <p>end if<br>
                        </p>
                        <p>call DMPlexDistribute</p>
                        <p>call DMDestroy(dmda_flow%da,ierr)<br>
                          CHKERRQ(ierr)</p>
                        <p>!c set the global mesh as distributed mesh<br>
                          dmda_flow%da = distributedMesh</p>
                        <p><br>
                        </p>
                        <p>After calling the above functions, the memory
                          usage for the test case (no. points 953,433,
                          nprocs 160) is shown below:<br>
                          rank 0 PETSc memory current MB    1610.39
                          PETSc memory maximum MB    1690.42<br>
                          rank 151 PETSc memory current MB     105.00
                          PETSc memory maximum MB     104.94<br>
                          rank 98 PETSc memory current MB     106.02
                          PETSc memory maximum MB     105.95<br>
                          rank 18 PETSc memory current MB     106.17
                          PETSc memory maximum MB     106.17</p>
                        <p>Is there any function available in the master
                          version that can release this memory?<br>
                        </p>
                      </div>
                    </blockquote>
                    <div>DMDestroy() releases this memory, UNLESS you
                      are holding other objects that refer to it, like a
                      vector from that DM.</div>
                  </div>
                </div>
              </blockquote>
              <p>Well, I have some labels set before distribution. After
                distribution, the labels values are collected but not
                destroyed. I will try this to see if it makes big
                difference.</p>
            </div>
          </blockquote>
          <div>Labels should be destroyed with the DM. Just make a small
            code that does nothing but distribute the mesh and end. If
            you</div>
          <div>run with -malloc_test you should see if everythign is
            destroyed properly.</div>
          <div><br>
          </div>
          <div>  Thanks,</div>
          <div><br>
          </div>
          <div>    Matt <br>
          </div>
        </div>
      </div>
    </blockquote>
    <p>Attached is the output run with -malloc_test using 2 processor.
      It's a big file. How can I quick check if something is not
      properly destroyed?</p>
    <p>Thanks,</p>
    <p>Danyang<br>
    </p>
    <blockquote type="cite"
cite="mid:CAMYG4GkonhQpuJop6rTt664rs5TDMTmnPsh6ydq-zbgYb5R6CA@mail.gmail.com">
      <div dir="ltr">
        <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 bgcolor="#FFFFFF">
              <p>Thanks,</p>
              <p>danyang<br>
              </p>
              <blockquote type="cite">
                <div dir="ltr">
                  <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>
                        <p>Thanks,</p>
                        <p>Danyang<br>
                        </p>
                        <div>On 2019-10-10 11:09 a.m., Mark Adams via
                          petsc-users wrote:<br>
                        </div>
                        <blockquote type="cite">
                          <div dir="ltr">Now that I think about it, the
                            partitioning and distribution can be done
                            with existing API, I would assume, like is
                            done with matrices. 
                            <div><br>
                            </div>
                            <div>I'm still wondering what the H5 format
                              is. I assume that it is not built for a
                              hardwired number of processes to read in
                              parallel and that the parallel read is
                              somewhat scalable.</div>
                          </div>
                          <br>
                          <div class="gmail_quote">
                            <div dir="ltr" class="gmail_attr">On Thu,
                              Oct 10, 2019 at 12:13 PM Mark Adams <<a
                                href="mailto:mfadams@lbl.gov"
                                target="_blank" moz-do-not-send="true">mfadams@lbl.gov</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">A related question, what is
                                the state of having something like a
                                distributed  DMPlexCreateFromCellList
                                method, but maybe your H5 efforts would
                                work. My bone modeling code is old and a
                                pain, but the apps specialized serial
                                mesh generator could write an H5 file
                                instead of the current FEAP file. Then
                                you reader, SNES and a large deformation
                                plasticity element in PetscFE could
                                replace my code, in the future.
                                <div><br>
                                </div>
                                <div>How does your H5 thing work? Is it
                                  basically a flat file (not
                                  partitioned) that is read in in
                                  parallel by slicing the cell lists,
                                  etc, using file seek or something
                                  equivalent, then reconstructing a
                                  local graph on each processor to give
                                  to say Parmetis, then completes the
                                  distribution with this
                                  reasonable partitioning? (this is what
                                  our current code does)</div>
                                <div><br>
                                </div>
                                <div>Thanks,</div>
                                <div>Mark</div>
                              </div>
                              <br>
                              <div class="gmail_quote">
                                <div dir="ltr" class="gmail_attr">On
                                  Thu, Oct 10, 2019 at 9:30 AM Dave May
                                  via petsc-users <<a
                                    href="mailto:petsc-users@mcs.anl.gov"
                                    target="_blank"
                                    moz-do-not-send="true">petsc-users@mcs.anl.gov</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><br>
                                  </div>
                                  <div><br>
                                    <div class="gmail_quote">
                                      <div dir="ltr">On Thu 10. Oct 2019
                                        at 15:15, Matthew Knepley <<a
href="mailto:knepley@gmail.com" target="_blank" moz-do-not-send="true">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 Thu, Oct 10,
                                            2019 at 9:10 AM Dave May
                                            <<a
                                              href="mailto:dave.mayhem23@gmail.com"
                                              target="_blank"
                                              moz-do-not-send="true">dave.mayhem23@gmail.com</a>>
                                            wrote:<br>
                                          </div>
                                        </div>
                                        <div dir="ltr">
                                          <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>On Thu 10. Oct 2019
                                                at 15:04, Matthew
                                                Knepley <<a
                                                  href="mailto:knepley@gmail.com"
                                                  target="_blank"
                                                  moz-do-not-send="true">knepley@gmail.com</a>>
                                                wrote:<br>
                                              </div>
                                              <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">
                                                      <div dir="ltr">On
                                                        Thu, Oct 10,
                                                        2019 at 8:41 AM
                                                        Dave May <<a
href="mailto:dave.mayhem23@gmail.com" target="_blank"
                                                          moz-do-not-send="true">dave.mayhem23@gmail.com</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>
                                                          <div>On Thu
                                                          10. Oct 2019
                                                          at 14:34,
                                                          Matthew
                                                          Knepley <<a
href="mailto:knepley@gmail.com" target="_blank" moz-do-not-send="true">knepley@gmail.com</a>>
                                                          wrote:<br>
                                                          </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 dir="ltr">On
                                                          Thu, Oct 10,
                                                          2019 at 8:31
                                                          AM Dave May
                                                          <<a
                                                          href="mailto:dave.mayhem23@gmail.com"
target="_blank" moz-do-not-send="true">dave.mayhem23@gmail.com</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">
                                                          <div dir="ltr">On
                                                          Thu, 10 Oct
                                                          2019 at 13:21,
                                                          Matthew
                                                          Knepley via
                                                          petsc-users
                                                          <<a
                                                          href="mailto:petsc-users@mcs.anl.gov"
target="_blank" moz-do-not-send="true">petsc-users@mcs.anl.gov</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">
                                                          <div dir="ltr">On
                                                          Wed, Oct 9,
                                                          2019 at 5:10
                                                          PM Danyang Su
                                                          via
                                                          petsc-users
                                                          <<a
                                                          href="mailto:petsc-users@mcs.anl.gov"
target="_blank" moz-do-not-send="true">petsc-users@mcs.anl.gov</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
                                                          bgcolor="#FFFFFF">
                                                          <p>Dear All,</p>
                                                          <p>I have a
                                                          question
                                                          regarding the
                                                          maximum memory
                                                          usage for the
                                                          scaling test.
                                                          My code is
                                                          written in
                                                          Fortran with
                                                          support for
                                                          both
                                                          structured
                                                          grid (DM) and
                                                          unstructured
                                                          grid (DMPlex).
                                                          It looks like
                                                          memory
                                                          consumption is
                                                          much larger
                                                          when DMPlex is
                                                          used and
                                                          finally causew
                                                          out_of_memory
                                                          problem. <br>
                                                          </p>
                                                          <p>Below are
                                                          some test
                                                          using both
                                                          structured
                                                          grid and
                                                          unstructured
                                                          grid. The
                                                          memory
                                                          consumption by
                                                          the code is
                                                          estimated
                                                          based on all
                                                          allocated
                                                          arrays and
                                                          PETSc memory
                                                          consumption is
                                                          estimated
                                                          based on
                                                          PetscMemoryGetMaximumUsage.
                                                          <br>
                                                          </p>
                                                          <p>I just
                                                          wonder why the
                                                          PETSc memory
                                                          consumption
                                                          does not
                                                          decrease when
                                                          number of
                                                          processors
                                                          increases. For
                                                          structured
                                                          grid (scenario
                                                          7-9), the
                                                          memory
                                                          consumption
                                                          decreases as
                                                          number of
                                                          processors
                                                          increases.
                                                          However, for
                                                          unstructured
                                                          grid case
                                                          (scenario
                                                          14-16), the
                                                          memory for
                                                          PETSc part
                                                          remains
                                                          unchanged.
                                                          When I run a
                                                          larger case,
                                                          the code
                                                          crashes
                                                          because memory
                                                          is ran out.
                                                          The same case
                                                          works on
                                                          another
                                                          cluster with
                                                          480GB memory
                                                          per node. Does
                                                          this make
                                                          sense?<br>
                                                          </p>
                                                          </div>
                                                          </blockquote>
                                                          <div>We would
                                                          need a finer
                                                          breakdown of
                                                          where memory
                                                          is being used.
                                                          I did this for
                                                          a paper:</div>
                                                          <div><br>
                                                          </div>
                                                          <div>  <a
href="https://agupubs.onlinelibrary.wiley.com/doi/full/10.1002/jgrb.50217"
target="_blank" moz-do-not-send="true">https://agupubs.onlinelibrary.wiley.com/doi/full/10.1002/jgrb.50217</a></div>
                                                          <div><br>
                                                          </div>
                                                          <div>If the
                                                          subdomains,
                                                          the halo sizes
                                                          can overwhelm
                                                          the basic
                                                          storage. It
                                                          looks like the
                                                          subdomains are
                                                          big here,</div>
                                                          <div>but
                                                          things are not
                                                          totally clear
                                                          to me. It
                                                          would be
                                                          helpful to
                                                          send the
                                                          output of
                                                          -log_view for
                                                          each case
                                                          since</div>
                                                          <div>PETSc
                                                          tries to keep
                                                          track of
                                                          allocated
                                                          memory. </div>
                                                          </div>
                                                          </div>
                                                          </blockquote>
                                                          <div><br>
                                                          </div>
                                                          <div>Matt -
                                                          I'd guess that
                                                          there is a
                                                          sequential
                                                          (non-partitioned)
                                                          mesh hanging
                                                          around in
                                                          memory.</div>
                                                          <div>Is it
                                                          possible that
                                                          he's created
                                                          the PLEX
                                                          object which
                                                          is loaded
                                                          sequentially
                                                          (stored and
                                                          retained in
                                                          memory and
                                                          never
                                                          released), and
                                                          then
                                                          afterwards
                                                          distributed?</div>
                                                          <div>This can
                                                          never happen
                                                          with the DMDA
                                                          and the table
                                                          verifies this.</div>
                                                          <div>If his
                                                          code using the
                                                          DMDA and
                                                          DMPLEX are as
                                                          identical as
                                                          possible
                                                          (albeit the DM
                                                          used), then a
                                                          sequential
                                                          mesh held in
                                                          memory seems
                                                          the likely
                                                          cause.</div>
                                                          </div>
                                                          </div>
                                                          </blockquote>
                                                          <div><br>
                                                          </div>
                                                          <div>Dang it,
                                                          Dave is always
                                                          right.</div>
                                                          <div><br>
                                                          </div>
                                                          <div>How to
                                                          prevent this?
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </blockquote>
                                                          <div
                                                          dir="auto"><br>
                                                          </div>
                                                          </div>
                                                          <div>
                                                          <div
                                                          dir="auto">I
                                                          thought
                                                          you/Lawrence/Vaclav/others...
                                                          had developed
                                                          and provided
                                                          support  for a
                                                          parallel
                                                          DMPLEX load
                                                          via a suitably
                                                          defined plex
                                                          specific H5
                                                          mesh file.</div>
                                                          </div>
                                                        </blockquote>
                                                        <div><br>
                                                        </div>
                                                        <div>We have,
                                                          but these
                                                          tests looked
                                                          like generated
                                                          meshes.</div>
                                                      </div>
                                                    </div>
                                                  </blockquote>
                                                  <div dir="auto"><br>
                                                  </div>
                                                  <div dir="auto">Great. </div>
                                                  <div dir="auto"><br>
                                                  </div>
                                                  <div dir="auto">So
                                                    would a solution to
                                                    the problem be to
                                                    have the user modify
                                                    their code in the
                                                    follow way:</div>
                                                  <div dir="auto">* they
                                                    move the mesh gen
                                                    stage into a
                                                    seperate exec which
                                                    they call offline
                                                    (on a fat node with
                                                    lots of memory), and
                                                    dump the appropriate
                                                    file</div>
                                                  <div dir="auto">* they
                                                    change their
                                                    existing application
                                                    to simply load that
                                                    file in parallel.</div>
                                                </div>
                                              </div>
                                            </blockquote>
                                            <div><br>
                                            </div>
                                          </div>
                                        </div>
                                        <div dir="ltr">
                                          <div class="gmail_quote">
                                            <div>Yes.</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>
                                                <div class="gmail_quote">
                                                  <div dir="auto">If
                                                    there were examples
                                                    illustrating how to
                                                    create the file
                                                    which can be loaded
                                                    in parallel I think
                                                    it would be very
                                                    helpful for the user
                                                    (and many others)</div>
                                                </div>
                                              </div>
                                            </blockquote>
                                            <div><br>
                                            </div>
                                            <div>I think Vaclav is going
                                              to add his examples as
                                              soon as we fix this
                                              parallel interpolation
                                              bug. I am praying for time
                                              in the latter</div>
                                            <div>part of October to do
                                              this.</div>
                                          </div>
                                        </div>
                                      </blockquote>
                                      <div dir="auto"><br>
                                      </div>
                                      <div dir="auto"><br>
                                      </div>
                                      <div dir="auto">Excellent news -
                                        thanks for the update and info.</div>
                                      <div dir="auto"><br>
                                      </div>
                                      <div dir="auto">Cheers</div>
                                      <div dir="auto">Dave</div>
                                      <div dir="auto"><br>
                                      </div>
                                      <div dir="auto"><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 class="gmail_quote">
                                            <div><br>
                                            </div>
                                            <div>  Thanks,</div>
                                            <div><br>
                                            </div>
                                            <div>    Matt</div>
                                          </div>
                                        </div>
                                        <div dir="ltr">
                                          <div class="gmail_quote">
                                            <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>
                                                <div class="gmail_quote">
                                                  <div dir="auto">Cheers</div>
                                                  <div dir="auto">Dave</div>
                                                  <div dir="auto"><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
                                                        class="gmail_quote">
                                                        <div><br>
                                                        </div>
                                                        <div>  Thanks,</div>
                                                        <div><br>
                                                        </div>
                                                        <div>    Matt</div>
                                                      </div>
                                                    </div>
                                                    <div dir="ltr">
                                                      <div
                                                        class="gmail_quote">
                                                        <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>
                                                          <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">
                                                          <div
                                                          class="gmail_quote">
                                                          <div>Since it
                                                          looks like you
                                                          are okay with
                                                          fairly regular
                                                          meshes, I
                                                          would
                                                          construct the</div>
                                                          <div>coarsest
                                                          mesh you can,
                                                          and then use</div>
                                                          <div><br>
                                                          </div>
                                                          <div> 
                                                          -dm_refine
                                                          <k></div>
                                                          <div><br>
                                                          </div>
                                                          <div>which is
                                                          activated by
                                                          DMSetFromOptions().
                                                          Make sure to
                                                          call it after
DMPlexDistribute(). It will regularly</div>
                                                          <div>refine in
                                                          parallel and
                                                          should show
                                                          good memory
                                                          scaling as
                                                          Dave says.</div>
                                                          <div><br>
                                                          </div>
                                                          <div>  Thanks,</div>
                                                          <div><br>
                                                          </div>
                                                          <div>   
                                                           Matt </div>
                                                          </div>
                                                          </div>
                                                          <div dir="ltr">
                                                          <div
                                                          class="gmail_quote">
                                                          <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
                                                          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">
                                                          <div
                                                          class="gmail_quote">
                                                          <div><br>
                                                          </div>
                                                          <div>  Thanks,</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> </p>
                                                          <p> </p>
                                                          <table
                                                          cellspacing="0"
                                                          border="0">
                                                          <colgroup
                                                          width="85"
                                                          span="5"></colgroup>
                                                          <colgroup
                                                          width="83"></colgroup>
                                                          <colgroup
                                                          width="125"></colgroup>
                                                          <colgroup
                                                          width="211"></colgroup>
                                                          <colgroup
                                                          width="125"></colgroup>
                                                          <colgroup
                                                          width="175"></colgroup>
                                                          <colgroup
                                                          width="119"></colgroup>
                                                          <tbody>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">scenario</td>
                                                          <td
                                                          align="right">no.
                                                          points</td>
                                                          <td
                                                          align="right">cell
                                                          type</td>
                                                          <td
                                                          align="right">DMPLex</td>
                                                          <td
                                                          align="right">nprocs</td>
                                                          <td
                                                          align="right">no.
                                                          nodes</td>
                                                          <td
                                                          align="right">mem
                                                          per node GB</td>
                                                          <td
                                                          align="right">solver</td>
                                                          <td
                                                          align="right">Rank
                                                          0 memory MB</td>
                                                          <td
                                                          align="right">Rank
                                                          0 petsc memory
                                                          MB</td>
                                                          <td
                                                          align="right">Runtime
                                                          (sec)</td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">1</td>
                                                          <td
                                                          align="right">2121</td>
                                                          <td
                                                          align="right">rectangle</td>
                                                          <td
                                                          align="right">no</td>
                                                          <td
                                                          align="right">40</td>
                                                          <td
                                                          align="right">1</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right">0.21</td>
                                                          <td
                                                          align="right">41.6</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">2</td>
                                                          <td
                                                          align="right">8241</td>
                                                          <td
                                                          align="right">rectangle</td>
                                                          <td
                                                          align="right">no</td>
                                                          <td
                                                          align="right">40</td>
                                                          <td
                                                          align="right">1</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right">0.59</td>
                                                          <td
                                                          align="right">51.84</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">3</td>
                                                          <td
                                                          align="right">32481</td>
                                                          <td
                                                          align="right">rectangle</td>
                                                          <td
                                                          align="right">no</td>
                                                          <td
                                                          align="right">40</td>
                                                          <td
                                                          align="right">1</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right">1.95</td>
                                                          <td
                                                          align="right">59.1</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">4</td>
                                                          <td
                                                          align="right">128961</td>
                                                          <td
                                                          align="right">rectangle</td>
                                                          <td
                                                          align="right">no</td>
                                                          <td
                                                          align="right">40</td>
                                                          <td
                                                          align="right">1</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right">7.05</td>
                                                          <td
                                                          align="right">89.71</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">5</td>
                                                          <td
                                                          align="right">513921</td>
                                                          <td
                                                          align="right">rectangle</td>
                                                          <td
                                                          align="right">no</td>
                                                          <td
                                                          align="right">40</td>
                                                          <td
                                                          align="right">1</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right">26.76</td>
                                                          <td
                                                          align="right">110.58</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">6</td>
                                                          <td
                                                          align="right">2051841</td>
                                                          <td
                                                          align="right">rectangle</td>
                                                          <td
                                                          align="right">no</td>
                                                          <td
                                                          align="right">40</td>
                                                          <td
                                                          align="right">1</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right">104.21</td>
                                                          <td
                                                          align="right">232.05</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right"><b>7</b></td>
                                                          <td
                                                          align="right"><b>8199681</b></td>
                                                          <td
                                                          align="right"><b>rectangle</b></td>
                                                          <td
                                                          align="right"><b>no</b></td>
                                                          <td
                                                          align="right"><b>40</b></td>
                                                          <td
                                                          align="right"><b>1</b></td>
                                                          <td
                                                          align="right"><b>200</b></td>
                                                          <td
                                                          align="right"><b>GMRES,Hypre
                                                          preconditioner</b></td>
                                                          <td
                                                          align="right"><b>411.26</b></td>
                                                          <td
                                                          align="right"><b>703.27</b></td>
                                                          <td
                                                          align="right"><b>140.29</b></td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right"><b>8</b></td>
                                                          <td
                                                          align="right"><b>8199681</b></td>
                                                          <td
                                                          align="right"><b>rectangle</b></td>
                                                          <td
                                                          align="right"><b>no</b></td>
                                                          <td
                                                          align="right"><b>80</b></td>
                                                          <td
                                                          align="right"><b>2</b></td>
                                                          <td
                                                          align="right"><b>200</b></td>
                                                          <td
                                                          align="right"><b>GMRES,Hypre
                                                          preconditioner</b></td>
                                                          <td
                                                          align="right"><b>206.6</b></td>
                                                          <td
                                                          align="right"><b>387.25</b></td>
                                                          <td
                                                          align="right"><b>62.04</b></td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right"><b>9</b></td>
                                                          <td
                                                          align="right"><b>8199681</b></td>
                                                          <td
                                                          align="right"><b>rectangle</b></td>
                                                          <td
                                                          align="right"><b>no</b></td>
                                                          <td
                                                          align="right"><b>160</b></td>
                                                          <td
                                                          align="right"><b>4</b></td>
                                                          <td
                                                          align="right"><b>200</b></td>
                                                          <td
                                                          align="right"><b>GMRES,Hypre
                                                          preconditioner</b></td>
                                                          <td
                                                          align="right"><b>104.28</b></td>
                                                          <td
                                                          align="right"><b>245.3</b></td>
                                                          <td
                                                          align="right"><b>32.76</b></td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">10</td>
                                                          <td
                                                          align="right">2121</td>
                                                          <td
                                                          align="right">triangle</td>
                                                          <td
                                                          align="right">yes</td>
                                                          <td
                                                          align="right">40</td>
                                                          <td
                                                          align="right">1</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right">0.49</td>
                                                          <td
                                                          align="right">61.78</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">11</td>
                                                          <td
                                                          align="right">15090</td>
                                                          <td
                                                          align="right">triangle</td>
                                                          <td
                                                          align="right">yes</td>
                                                          <td
                                                          align="right">40</td>
                                                          <td
                                                          align="right">1</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right">2.32</td>
                                                          <td
                                                          align="right">96.61</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">12</td>
                                                          <td
                                                          align="right">59847</td>
                                                          <td
                                                          align="right">triangle</td>
                                                          <td
                                                          align="right">yes</td>
                                                          <td
                                                          align="right">40</td>
                                                          <td
                                                          align="right">1</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right">8.28</td>
                                                          <td
                                                          align="right">176.14</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">13</td>
                                                          <td
                                                          align="right">238568</td>
                                                          <td
                                                          align="right">triangle</td>
                                                          <td
                                                          align="right">yes</td>
                                                          <td
                                                          align="right">40</td>
                                                          <td
                                                          align="right">1</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right">31.89</td>
                                                          <td
                                                          align="right">573.73</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right"><b>14</b></td>
                                                          <td
                                                          align="right"><b>953433</b></td>
                                                          <td
                                                          align="right"><b>triangle</b></td>
                                                          <td
                                                          align="right"><b>yes</b></td>
                                                          <td
                                                          align="right"><b>40</b></td>
                                                          <td
                                                          align="right"><b>1</b></td>
                                                          <td
                                                          align="right"><b>200</b></td>
                                                          <td
                                                          align="right"><b>GMRES,Hypre
                                                          preconditioner</b></td>
                                                          <td
                                                          align="right"><b>119.23</b></td>
                                                          <td
                                                          align="right"><b>2102.54</b></td>
                                                          <td
                                                          align="right"><b>44.11</b></td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right"><b>15</b></td>
                                                          <td
                                                          align="right"><b>953433</b></td>
                                                          <td
                                                          align="right"><b>triangle</b></td>
                                                          <td
                                                          align="right"><b>yes</b></td>
                                                          <td
                                                          align="right"><b>80</b></td>
                                                          <td
                                                          align="right"><b>2</b></td>
                                                          <td
                                                          align="right"><b>200</b></td>
                                                          <td
                                                          align="right"><b>GMRES,Hypre
                                                          preconditioner</b></td>
                                                          <td
                                                          align="right"><b>72.99</b></td>
                                                          <td
                                                          align="right"><b>2123.8</b></td>
                                                          <td
                                                          align="right"><b>24.36</b></td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right"><b>16</b></td>
                                                          <td
                                                          align="right"><b>953433</b></td>
                                                          <td
                                                          align="right"><b>triangle</b></td>
                                                          <td
                                                          align="right"><b>yes</b></td>
                                                          <td
                                                          align="right"><b>160</b></td>
                                                          <td
                                                          align="right"><b>4</b></td>
                                                          <td
                                                          align="right"><b>200</b></td>
                                                          <td
                                                          align="right"><b>GMRES,Hypre
                                                          preconditioner</b></td>
                                                          <td
                                                          align="right"><b>48.65</b></td>
                                                          <td
                                                          align="right"><b>2076.25</b></td>
                                                          <td
                                                          align="right"><b>14.87</b></td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">17</td>
                                                          <td
                                                          align="right">55770</td>
                                                          <td
                                                          align="right">prism</td>
                                                          <td
                                                          align="right">yes</td>
                                                          <td
                                                          align="right">40</td>
                                                          <td
                                                          align="right">1</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right">18.46</td>
                                                          <td
                                                          align="right">219.39</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">18</td>
                                                          <td
                                                          align="right">749814</td>
                                                          <td
                                                          align="right">prism</td>
                                                          <td
                                                          align="right">yes</td>
                                                          <td
                                                          align="right">40</td>
                                                          <td
                                                          align="right">1</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right">149.86</td>
                                                          <td
                                                          align="right">2412.39</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right">19</td>
                                                          <td
                                                          align="right">7000050</td>
                                                          <td
                                                          align="right">prism</td>
                                                          <td
                                                          align="right">yes</td>
                                                          <td
                                                          align="right">40
                                                          to 640</td>
                                                          <td
                                                          align="right">1
                                                          to 16</td>
                                                          <td
                                                          align="right">200</td>
                                                          <td
                                                          align="right">GMRES,Hypre
                                                          preconditioner</td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right">out_of_memory</td>
                                                          <td
                                                          align="left"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          <tr>
                                                          <td
                                                          height="17"
                                                          align="right"><b>20</b></td>
                                                          <td
                                                          align="right"><b>7000050</b></td>
                                                          <td
                                                          align="right"><b>prism</b></td>
                                                          <td
                                                          align="right"><b>yes</b></td>
                                                          <td
                                                          align="right"><b>64</b></td>
                                                          <td
                                                          align="right"><b>2</b></td>
                                                          <td
                                                          align="right"><b>480</b></td>
                                                          <td
                                                          align="right"><b>GMRES,Hypre
                                                          preconditioner</b></td>
                                                          <td
                                                          align="right"><b>890.92</b></td>
                                                          <td
                                                          align="right"><b>17214.41</b></td>
                                                          <td
                                                          align="right"><br>
                                                          </td>
                                                          </tr>
                                                          </tbody>
                                                          </table>
                                                          <p> </p>
                                                          <p>The error
                                                          information of
                                                          scenario 19 is
                                                          shown below:<br>
                                                          </p>
                                                          <p>kernel
                                                          messages
                                                          produced
                                                          during job
                                                          executions:<br>
                                                          [Oct 9 10:41]
                                                          mpiexec.hydra
                                                          invoked
                                                          oom-killer:
                                                          gfp_mask=0x200da,
                                                          order=0,
                                                          oom_score_adj=0<br>
                                                          [  +0.010274]
                                                          mpiexec.hydra
                                                          cpuset=/
                                                          mems_allowed=0-1<br>
                                                          [  +0.006680]
                                                          CPU: 2 PID:
                                                          144904 Comm:
                                                          mpiexec.hydra
                                                          Tainted:
                                                          G          
                                                          OE 
                                                          ------------  
3.10.0-862.14.4.el7.x86_64 #1<br>
                                                          [  +0.013365]
                                                          Hardware name:
                                                          Lenovo
                                                          ThinkSystem
                                                          SD530
                                                          -[7X21CTO1WW]-/-[7X21CTO1WW]-,
                                                          BIOS
                                                          -[TEE124N-1.40]-
                                                          06/12/2018<br>
                                                          [  +0.012866]
                                                          Call Trace:<br>
                                                          [  +0.003945] 
[<ffffffffb3313754>] dump_stack+0x19/0x1b<br>
                                                          [  +0.006995] 
[<ffffffffb330e91f>] dump_header+0x90/0x229<br>
                                                          [  +0.007121] 
[<ffffffffb2cfa982>] ? ktime_get_ts64+0x52/0xf0<br>
                                                          [  +0.007451] 
[<ffffffffb2d5141f>] ? delayacct_end+0x8f/0xb0<br>
                                                          [  +0.007393] 
[<ffffffffb2d9ac94>] oom_kill_process+0x254/0x3d0<br>
                                                          [  +0.007592] 
[<ffffffffb2d9a73d>] ? oom_unkillable_task+0xcd/0x120<br>
                                                          [  +0.007978] 
[<ffffffffb2d9a7e6>] ? find_lock_task_mm+0x56/0xc0<br>
                                                          [  +0.007729] 
[<ffffffffb2d9b4d6>] <b>out_of_memory+0x4b6/0x4f0</b><br>
                                                          [  +0.007358] 
[<ffffffffb330f423>] __alloc_pages_slowpath+0x5d6/0x724<br>
                                                          [  +0.008190] 
[<ffffffffb2da18b5>] __alloc_pages_nodemask+0x405/0x420</p>
                                                          <p>Thanks,</p>
                                                          <p>Danyang<br>
                                                          </p>
                                                          </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" moz-do-not-send="true">https://www.cse.buffalo.edu/~knepley/</a><br>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </blockquote>
                                                          </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" moz-do-not-send="true">https://www.cse.buffalo.edu/~knepley/</a><br>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </blockquote>
                                                          </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" moz-do-not-send="true">https://www.cse.buffalo.edu/~knepley/</a><br>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </div>
                                                          </div>
                                                        </div>
                                                      </div>
                                                    </div>
                                                  </blockquote>
                                                </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" moz-do-not-send="true">https://www.cse.buffalo.edu/~knepley/</a><br>
                                                      </div>
                                                    </div>
                                                  </div>
                                                </div>
                                              </div>
                                            </div>
                                          </div>
                                        </div>
                                      </blockquote>
                                    </div>
                                  </div>
                                </blockquote>
                              </div>
                            </blockquote>
                          </div>
                        </blockquote>
                      </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" moz-do-not-send="true">https://www.cse.buffalo.edu/~knepley/</a><br>
                              </div>
                            </div>
                          </div>
                        </div>
                      </div>
                    </div>
                  </div>
                </div>
              </blockquote>
            </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" moz-do-not-send="true">https://www.cse.buffalo.edu/~knepley/</a><br>
                    </div>
                  </div>
                </div>
              </div>
            </div>
          </div>
        </div>
      </div>
    </blockquote>
  </body>
</html>