<html><head><meta http-equiv="Content-Type" content="text/html charset=utf-8"></head><body style="word-wrap: break-word; -webkit-nbsp-mode: space; -webkit-line-break: after-white-space;" class="">It should be fixed now in maint and master<div class=""><br class=""></div><div class=""><a href="https://bitbucket.org/petsc/petsc/commits/4c8dd594d1988a0cbe282f8a37d9916f61e0c445" class="">https://bitbucket.org/petsc/petsc/commits/4c8dd594d1988a0cbe282f8a37d9916f61e0c445</a></div><div class=""><br class=""></div><div class="">Thanks for reporting the problem,</div><div class="">Stefano</div><div class=""><br class=""><div><blockquote type="cite" class=""><div class="">On Jun 19, 2017, at 10:46 PM, Stefano Zampini <<a href="mailto:stefano.zampini@gmail.com" class="">stefano.zampini@gmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><div dir="auto" class="">Franck,<div dir="auto" class=""><br class=""></div><div dir="auto" class="">Thanks. I'll​ get back soon with a fix.</div><div dir="auto" class=""><br class=""></div><div dir="auto" class="">Stefano</div><div dir="auto" class=""><br class=""></div></div><div class="gmail_extra"><br class=""><div class="gmail_quote">Il 19 Giu 2017 18:17, "Franck Houssen" <<a href="mailto:franck.houssen@inria.fr" class="">franck.houssen@inria.fr</a>> ha scritto:<br type="attribution" class=""><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div class=""><div style="font-family: 'times new roman', 'new york', times, serif; font-size: 12pt;" class=""><div class="">The problem was difficult to reduce as reducing make things disappear... Luckily, I believe I got it (or at least, it looks "like" the one I "really" have...).<br class=""></div><div class=""><br class=""></div><div class="">Seems that for square matrix, it works fine for csr and dense matrix. But, If I am not mistaken, it does not for dense rectangular matrix (still OK for csr).<br class=""></div><div class=""><br class=""></div><div class="">matISCSRDenseSquare.cpp: 2 procs, global 3x3 matrix, each proc adds a 2x2 local matrix in the global matrix.<br class=""></div><div class="">matISCSRDenseRect.cpp: 2 procs, global <strong class="">2</strong>x3 matrix, each proc adds a <strong class="">1</strong>x2 local <strong class="">vector</strong> in the global matrix.<br class=""></div><div class=""><br class=""></div><div class="">reminder: running debian/testing with gcc-6.3 + petsc-3.7.6</div><div class=""><br class=""></div><div class="">Franck<br class=""></div><div class=""><br class=""></div><div class="">>> mpirun -n 2 ./matISCSRDenseSquare.exe csr; mpirun -n 2 ./matISCSRDenseSquare.exe dense<br class="">csr<br class="">csr<br class="">Mat Object: 2 MPI processes<br class="">  type: is<br class="">  Mat Object:   1 MPI processes<br class="">    type: seqaij<br class="">row 0: (0, 1.)  (1, 0.) <br class="">row 1: (0, 0.)  (1, 1.) <br class="">  Mat Object:   1 MPI processes<br class="">    type: seqaij<br class="">row 0: (0, 1.)  (1, 0.) <br class="">row 1: (0, 0.)  (1, 1.) <br class="">dense<br class="">dense<br class="">Mat Object: 2 MPI processes<br class="">  type: is<br class="">  Mat Object:   1 MPI processes<br class="">    type: seqdense<br class="">1.0000000000000000e+00 0.0000000000000000e+00 <br class="">0.0000000000000000e+00 1.0000000000000000e+00 <br class="">  Mat Object:   1 MPI processes<br class="">    type: seqdense<br class="">1.0000000000000000e+00 0.0000000000000000e+00 <br class="">0.0000000000000000e+00 1.0000000000000000e+00 <br class=""><br class=""></div><div class="">>> mpirun -n 2 ./matISCSRDenseRect.exe csr; mpirun -n 2 ./matISCSRDenseRect.exe dense<br class="">csr<br class="">csr<br class="">Mat Object: 2 MPI processes<br class="">  type: is<br class="">  Mat Object:   1 MPI processes<br class="">    type: seqaij<br class="">row 0: (0, 1.)  (1, 0.) <br class="">  Mat Object:   1 MPI processes<br class="">    type: seqaij<br class="">row 0: (0, 1.)  (1, 0.) <br class="">dense<br class="">dense<br class="">[0]PETSC ERROR: --------------------- Error Message ------------------------------<wbr class="">------------------------------<wbr class="">--<br class="">[0]PETSC ERROR: Argument out of range<br class="">[0]PETSC ERROR: [1]PETSC ERROR: --------------------- Error Message ------------------------------<wbr class="">------------------------------<wbr class="">--<br class="">[1]PETSC ERROR: Argument out of range<br class="">[1]PETSC ERROR: New nonzero at (0,1) caused a malloc<br class="">Use MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_<wbr class="">ERR, PETSC_FALSE) to turn off this check<br class="">[1]PETSC ERROR: See <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank" class="">http://www.mcs.anl.gov/petsc/<wbr class="">documentation/faq.html</a> for trouble shooting.<br class="">[1]PETSC ERROR: Petsc Release Version 3.7.6, Apr, 24, 2017 <br class="">[1]PETSC ERROR: ./matISCSRDenseRect.exe on a arch-linux2-c-debug named yoda by fghoussen Mon Jun 19 18:08:58 2017<br class="">New nonzero at (0,1) caused a malloc<br class="">Use MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_<wbr class="">ERR, PETSC_FALSE) to turn off this check<br class="">[0]PETSC ERROR: See <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank" class="">http://www.mcs.anl.gov/petsc/<wbr class="">documentation/faq.html</a> for trouble shooting.<br class="">[0]PETSC ERROR: Petsc Release Version 3.7.6, Apr, 24, 2017 <br class="">[1]PETSC ERROR: Configure options --prefix=/home/fghoussen/<wbr class="">Documents/INRIA/petsc-3.7.6/<wbr class="">local --with-mpi=1 --with-pthread=1 --download-f2cblaslapack=yes --download-mumps=yes --download-scalapack=yes --download-superlu=yes --download-suitesparse=yes<br class="">[1]PETSC ERROR: #1 MatSetValues_MPIAIJ() line 616 in /home/fghoussen/Documents/<wbr class="">INRIA/petsc-3.7.6/src/mat/<wbr class="">impls/aij/mpi/mpiaij.c<br class="">[1]PETSC ERROR: #2 MatSetValues() line 1190 in /home/fghoussen/Documents/<wbr class="">INRIA/petsc-3.7.6/src/mat/<wbr class="">interface/matrix.c<br class="">[1]PETSC ERROR: #3 MatSetValuesLocal() line 2053 in /home/fghoussen/Documents/<wbr class="">INRIA/petsc-3.7.6/src/mat/<wbr class="">interface/matrix.c<br class="">[1]PETSC ERROR: #4 MatISGetMPIXAIJ_IS() line 365 in /home/fghoussen/Documents/<wbr class="">INRIA/petsc-3.7.6/src/mat/<wbr class="">impls/is/matis.c<br class="">[1]PETSC ERROR: #5 MatISGetMPIXAIJ() line 437 in /home/fghoussen/Documents/<wbr class="">INRIA/petsc-3.7.6/src/mat/<wbr class="">impls/is/matis.c<br class="">[0]PETSC ERROR: ./matISCSRDenseRect.exe on a arch-linux2-c-debug named yoda by fghoussen Mon Jun 19 18:08:58 2017<br class="">[0]PETSC ERROR: Configure options --prefix=/home/fghoussen/<wbr class="">Documents/INRIA/petsc-3.7.6/<wbr class="">local --with-mpi=1 --with-pthread=1 --download-f2cblaslapack=yes --download-mumps=yes --download-scalapack=yes --download-superlu=yes --download-suitesparse=yes<br class="">[0]PETSC ERROR: #1 MatSetValues_MPIAIJ() line 582 in /home/fghoussen/Documents/<wbr class="">INRIA/petsc-3.7.6/src/mat/<wbr class="">impls/aij/mpi/mpiaij.c<br class="">[0]PETSC ERROR: #2 MatSetValues() line 1190 in /home/fghoussen/Documents/<wbr class="">INRIA/petsc-3.7.6/src/mat/<wbr class="">interface/matrix.c<br class="">[0]PETSC ERROR: #3 MatSetValuesLocal() line 2053 in /home/fghoussen/Documents/<wbr class="">INRIA/petsc-3.7.6/src/mat/<wbr class="">interface/matrix.c<br class="">[0]PETSC ERROR: #4 MatISGetMPIXAIJ_IS() line 365 in /home/fghoussen/Documents/<wbr class="">INRIA/petsc-3.7.6/src/mat/<wbr class="">impls/is/matis.c<br class="">[0]PETSC ERROR: #5 MatISGetMPIXAIJ() line 437 in /home/fghoussen/Documents/<wbr class="">INRIA/petsc-3.7.6/src/mat/<wbr class="">impls/is/matis.c<br class="">Mat Object: 2 MPI processes<br class="">  type: is<br class="">  Mat Object:   1 MPI processes<br class="">    type: seqdense<br class="">1.0000000000000000e+00 0.0000000000000000e+00 <br class="">  Mat Object:   1 MPI processes<br class="">    type: seqdense<br class="">1.0000000000000000e+00 0.0000000000000000e+00 <br class=""><br class=""></div><div class="">>> diff matISCSRDenseSquare.cpp matISCSRDenseRect.cpp<br class="">3c3<br class="">< // ~> g++ -o matISCSRDenseSquare.exe matISCSRDenseSquare.cpp -lpetsc -lm; mpirun -n 2 matISCSRDenseSquare.exe<br class="">---<br class="">> // ~> g++ -o matISCSRDenseRect.exe matISCSRDenseRect.cpp -lpetsc -lm; mpirun -n 2 matISCSRDenseRect.exe<br class="">24c24<br class=""><   ISLocalToGlobalMappingCreate(<wbr class="">PETSC_COMM_WORLD, 1, 2, localIdx, PETSC_COPY_VALUES, &rmap);<br class="">---<br class="">>   ISLocalToGlobalMappingCreate(<wbr class="">PETSC_COMM_WORLD, 1, 1, &rank, PETSC_COPY_VALUES, &rmap);<br class="">29c29<br class=""><   MatCreateIS(PETSC_COMM_WORLD, 1, PETSC_DECIDE, PETSC_DECIDE, 3, 3, rmap, cmap, &A);<br class="">---<br class="">>   MatCreateIS(PETSC_COMM_WORLD, 1, PETSC_DECIDE, PETSC_DECIDE, 2, 3, rmap, cmap, &A);<br class="">32,33c32,33<br class=""><   if (matType == "csr") {cout << matType << endl; MatCreateSeqAIJ(PETSC_COMM_<wbr class="">SELF, 2, 2, 2, NULL, &Aloc);}<br class=""><   else {cout << matType << endl; MatCreateSeqDense(PETSC_COMM_<wbr class="">SELF, 2, 2, NULL, &Aloc);}<br class="">---<br class="">>   if (matType == "csr") {cout << matType << endl; MatCreateSeqAIJ(PETSC_COMM_<wbr class="">SELF, 1, 2, 2, NULL, &Aloc);}<br class="">>   else {cout << matType << endl; MatCreateSeqDense(PETSC_COMM_<wbr class="">SELF, 1, 2, NULL, &Aloc);}<br class="">35,36c35,36<br class=""><   PetscScalar localVal[4] = {1., 0., 0., 1.};<br class=""><   MatSetValues(Aloc, 2, localIdx, 2, localIdx, localVal, ADD_VALUES); // Add local 2x2 matrix<br class="">---<br class="">>   PetscScalar localVal[2] = {1., 0.}; PetscInt oneLocalRow = 0;<br class="">>   MatSetValues(Aloc, 1, &oneLocalRow, 2, localIdx, localVal, ADD_VALUES); // Add local row<br class=""><br class=""></div><hr id="m_5842485459567131580zwchr" class=""><blockquote style="border-left-width: 2px; border-left-style: solid; border-left-color: rgb(16, 16, 255); margin-left: 5px; padding-left: 5px; font-weight: normal; font-style: normal; text-decoration: none; font-family: Helvetica, Arial, sans-serif; font-size: 12pt;" class=""><b class="">De: </b>"Stefano Zampini" <<a href="mailto:stefano.zampini@gmail.com" target="_blank" class="">stefano.zampini@gmail.com</a>><br class=""><b class="">À: </b>"Franck Houssen" <<a href="mailto:franck.houssen@inria.fr" target="_blank" class="">franck.houssen@inria.fr</a>><br class=""><b class="">Cc: </b>"PETSc users list" <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank" class="">petsc-users@mcs.anl.gov</a>><br class=""><b class="">Envoyé: </b>Lundi 19 Juin 2017 15:25:35<br class=""><b class="">Objet: </b>Re: [petsc-users] Building MatIS with dense local matrix ?<br class=""><div class=""><br class=""></div><div dir="auto" class="">Can you send a minimal working example so that I can fix the code?<div dir="auto" class=""><br class=""></div><div dir="auto" class="">Thanks</div><div dir="auto" class="">Stefano</div></div><div class="gmail_extra"><br class=""><div class="gmail_quote">Il 19 Giu 2017 15:20, "Franck Houssen" <<a href="mailto:franck.houssen@inria.fr" target="_blank" class="">franck.houssen@inria.fr</a>> ha scritto:<br class=""><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div class=""><div style="font-family: 'times new roman', 'new york', times, serif; font-size: 12pt;" class=""><div class="">Hi,<br class=""></div><div class=""><br class=""></div><div class="">I try to call MatISGetMPIXAIJ on a MatIS (A) that has been feed locally by sequential (Aloc) dense matrix.</div><div class="">Seems this ends up with this error: [0]PETSC ERROR: New nonzero at (0,1) caused a malloc. Is this a known error / limitation ? (not supposed to work with dense matrix ?)<br class=""></div><div class=""><br class=""></div><div class="">This (pseudo code) works fine:<br class=""></div><div class="">MatCreateIS(..., A)</div><div class="">MatCreateSeqAIJ(..., Aloc)<br class=""></div><div class="">MatISSetLocalMat(pcA, pcALoc) <br class=""></div><div class="">MatISGetMPIXAIJ(A, ...) // OK !<br class=""></div><div class=""><br class=""></div><div class="">When I try to replace MatCreateSeqAIJ(..., Aloc) with MatCreateSeqDense(..., Aloc), it does no more work. <br class=""></div><div class=""><br class=""></div><div class="">Franck<br class=""></div><div class=""><br class=""></div><div class="">PS: running debian/testing with gcc-6.3 + petsc-3.7.6<br class=""></div></div></div></blockquote></div></div></blockquote><div class=""><br class=""></div></div></div></blockquote></div></div>
</div></blockquote></div><br class=""></div></body></html>