<div dir="ltr"><div dir="ltr">On Wed, Mar 13, 2019 at 12:04 PM Manuel Colera Rico <<a href="mailto:m.colera@upm.es">m.colera@upm.es</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>After adding that line the problem gets fixed.</p>
    <p></p></div></blockquote><div>Junchao, I submitted a PR against maint for this as you.</div><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>Regards,</p>
    <p>Manuel</p>
    <p>---<br>
    </p>
    <div class="gmail-m_-5289257232827223512moz-cite-prefix">On 3/13/19 3:13 PM, Zhang, Junchao
      wrote:<br>
    </div>
    <blockquote type="cite">
      
      <div dir="ltr">
        <div dir="ltr">Manuel, </div>
        <div dir="ltr">  Could you try to add this line<br>
               sbaij->free_imax_ilen = PETSC_TRUE;<br>
           after line 2431 in
          /opt/PETSc_library/petsc-3.10.4/src/mat/impls/sbaij/seq/sbaij.c<br>
          <br>
        </div>
        <div dir="ltr"> PS: Matt, this bug looks unrelated to my
          VecRestoreArrayRead_Nest fix.
          <div><font color="#500050"><br clear="all">
            </font>
            <div>
              <div dir="ltr" class="gmail-m_-5289257232827223512gmail_signature">
                <div dir="ltr">--Junchao Zhang</div>
              </div>
            </div>
            <br>
          </div>
        </div>
      </div>
      <br>
      <div class="gmail_quote">
        <div dir="ltr" class="gmail_attr">On Wed, Mar 13, 2019 at 9:05
          AM Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:<br>
        </div>
        <blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">
          <div dir="ltr">
            <div dir="ltr">
              <div dir="ltr">On Wed, Mar 13, 2019 at 9:44 AM Manuel
                Colera Rico via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">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">
                  Yes:<br>
                  <br>
                  [ 0]8416 bytes MatCreateSeqSBAIJWithArrays() line 2431
                  in <br>
/opt/PETSc_library/petsc-3.10.4/src/mat/impls/sbaij/seq/sbaij.c<br>
                  [ 0]8416 bytes MatCreateSeqSBAIJWithArrays() line 2431
                  in <br>
/opt/PETSc_library/petsc-3.10.4/src/mat/impls/sbaij/seq/sbaij.c<br>
                  [ 0]4544 bytes MatCreateSeqSBAIJWithArrays() line 2431
                  in <br>
/opt/PETSc_library/petsc-3.10.4/src/mat/impls/sbaij/seq/sbaij.c<br>
                  [ 0]4544 bytes MatCreateSeqSBAIJWithArrays() line 2431
                  in <br>
/opt/PETSc_library/petsc-3.10.4/src/mat/impls/sbaij/seq/sbaij.c<br>
                </blockquote>
                <div><br>
                </div>
                <div>Junchao, do imax and ilen get missed in the Destroy
                  with the user provides arrays?</div>
                <div><br>
                </div>
                <div>  <a href="https://bitbucket.org/petsc/petsc/src/06a3e802b3873ffbfd04b71a0821522327dd9b04/src/mat/impls/sbaij/seq/sbaij.c#lines-2431" target="_blank">https://bitbucket.org/petsc/petsc/src/06a3e802b3873ffbfd04b71a0821522327dd9b04/src/mat/impls/sbaij/seq/sbaij.c#lines-2431</a></div>
                <div><br>
                </div>
                <div>    Matt</div>
                <div> </div>
                <blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">
                  I have checked that I have destroyed all the MatNest
                  matrices and all <br>
                  the submatrices individually.<br>
                  <br>
                  Manuel<br>
                  <br>
                  ---<br>
                  <br>
                  On 3/13/19 2:28 PM, Jed Brown wrote:<br>
                  > Is there any output if you run with -malloc_dump?<br>
                  ><br>
                  > Manuel Colera Rico via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>>
                  writes:<br>
                  ><br>
                  >> Hi, Junchao,<br>
                  >><br>
                  >> I have installed the newest version of PETSc
                  and it works fine. I just<br>
                  >> get the following memory leak warning:<br>
                  >><br>
                  >> Direct leak of 28608 byte(s) in 12 object(s)
                  allocated from:<br>
                  >>       #0 0x7f1ddd5caa38 in
                  __interceptor_memalign<br>
                  >>
                  ../../../../gcc-8.1.0/libsanitizer/asan/asan_malloc_linux.cc:111<br>
                  >>       #1 0x7f1ddbef1213 in PetscMallocAlign<br>
                  >>
(/opt/PETSc_library/petsc-3.10.4/mcr_20190313/lib/libpetsc.so.3.10+0x150213)<br>
                  >><br>
                  >> Thank you,<br>
                  >><br>
                  >> Manuel<br>
                  >><br>
                  >> ---<br>
                  >><br>
                  >> On 3/12/19 7:08 PM, Zhang, Junchao wrote:<br>
                  >>> Hi, Manuel,<br>
                  >>>    I recently fixed a problem in
                  VecRestoreArrayRead. Basically, I<br>
                  >>> added VecRestoreArrayRead_Nest. Could you
                  try the master branch of<br>
                  >>> PETSc to see if it fixes your problem?<br>
                  >>>    Thanks.<br>
                  >>><br>
                  >>> --Junchao Zhang<br>
                  >>><br>
                  >>><br>
                  >>> On Mon, Mar 11, 2019 at 6:56 AM Manuel
                  Colera Rico via petsc-users<br>
                  >>> <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>
                  <mailto:<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>>>
                  wrote:<br>
                  >>><br>
                  >>>      Hello,<br>
                  >>><br>
                  >>>      I need to solve a 2*2 block linear
                  system. The matrices A_00, A_01,<br>
                  >>>      A_10, A_11 are constructed
                  separately via<br>
                  >>>      MatCreateSeqAIJWithArrays and<br>
                  >>>      MatCreateSeqSBAIJWithArrays. Then, I
                  construct the full system matrix<br>
                  >>>      with MatCreateNest, and use
                  MatNestGetISs and PCFieldSplitSetIS to<br>
                  >>>      set<br>
                  >>>      up the PC, trying to follow the
                  procedure described here:<br>
                  >>>      <a href="https://www.mcs.anl.gov/petsc/petsc-current/src/snes/examples/tutorials/ex70.c.html" rel="noreferrer" target="_blank">
https://www.mcs.anl.gov/petsc/petsc-current/src/snes/examples/tutorials/ex70.c.html</a>.<br>
                  >>><br>
                  >>>      However, when I run the code with
                  Leak Sanitizer, I get the<br>
                  >>>      following error:<br>
                  >>><br>
                  >>>     
                  =================================================================<br>
                  >>>      ==54927==ERROR: AddressSanitizer:
                  attempting free on address which<br>
                  >>>      was<br>
                  >>>      not malloc()-ed: 0x627000051ab8 in
                  thread T0<br>
                  >>>           #0 0x7fbd95c08f30 in
                  __interceptor_free<br>
                  >>>     
                  ../../../../gcc-8.1.0/libsanitizer/asan/asan_malloc_linux.cc:66<br>
                  >>>           #1 0x7fbd92b99dcd in
                  PetscFreeAlign<br>
                  >>>     
(/opt/PETSc_library/petsc/manuel_OpenBLAS_petsc/lib/libpetsc.so.3.8+0x146dcd)<br>
                  >>>           #2 0x7fbd92ce0178 in
                  VecRestoreArray_Nest<br>
                  >>>     
(/opt/PETSc_library/petsc/manuel_OpenBLAS_petsc/lib/libpetsc.so.3.8+0x28d178)<br>
                  >>>           #3 0x7fbd92cd627d in
                  VecRestoreArrayRead<br>
                  >>>     
(/opt/PETSc_library/petsc/manuel_OpenBLAS_petsc/lib/libpetsc.so.3.8+0x28327d)<br>
                  >>>           #4 0x7fbd92d1189e in
                  VecScatterBegin_SSToSS<br>
                  >>>     
(/opt/PETSc_library/petsc/manuel_OpenBLAS_petsc/lib/libpetsc.so.3.8+0x2be89e)<br>
                  >>>           #5 0x7fbd92d1a414 in
                  VecScatterBegin<br>
                  >>>     
(/opt/PETSc_library/petsc/manuel_OpenBLAS_petsc/lib/libpetsc.so.3.8+0x2c7414)<br>
                  >>>           #6 0x7fbd934a999c in
                  PCApply_FieldSplit<br>
                  >>>     
(/opt/PETSc_library/petsc/manuel_OpenBLAS_petsc/lib/libpetsc.so.3.8+0xa5699c)<br>
                  >>>           #7 0x7fbd93369071 in PCApply<br>
                  >>>     
(/opt/PETSc_library/petsc/manuel_OpenBLAS_petsc/lib/libpetsc.so.3.8+0x916071)<br>
                  >>>           #8 0x7fbd934efe77 in
                  KSPInitialResidual<br>
                  >>>     
(/opt/PETSc_library/petsc/manuel_OpenBLAS_petsc/lib/libpetsc.so.3.8+0xa9ce77)<br>
                  >>>           #9 0x7fbd9350272c in
                  KSPSolve_GMRES<br>
                  >>>     
(/opt/PETSc_library/petsc/manuel_OpenBLAS_petsc/lib/libpetsc.so.3.8+0xaaf72c)<br>
                  >>>           #10 0x7fbd934e3c01 in KSPSolve<br>
                  >>>     
(/opt/PETSc_library/petsc/manuel_OpenBLAS_petsc/lib/libpetsc.so.3.8+0xa90c01)<br>
                  >>><br>
                  >>>      Disabling Leak Sanitizer also
                  outputs an "invalid pointer" error.<br>
                  >>><br>
                  >>>      Did I forget something when writing
                  the code?<br>
                  >>><br>
                  >>>      Thank you,<br>
                  >>><br>
                  >>>      Manuel<br>
                  >>><br>
                  >>>      ---<br>
                  >>><br>
                </blockquote>
              </div>
              <br clear="all">
              <div><br>
              </div>
              -- <br>
              <div dir="ltr" class="gmail-m_-5289257232827223512gmail-m_-161896209400255159gmail_signature">
                <div dir="ltr">
                  <div>
                    <div dir="ltr">
                      <div>
                        <div dir="ltr">
                          <div>What most experimenters take for granted
                            before they begin their experiments is
                            infinitely more interesting than any results
                            to which their experiments lead.<br>
                            -- Norbert Wiener</div>
                          <div><br>
                          </div>
                          <div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br>
                          </div>
                        </div>
                      </div>
                    </div>
                  </div>
                </div>
              </div>
            </div>
          </div>
        </blockquote>
      </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">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>