<div dir="ltr">Thanks Mark, both your suggestions of using pc hypre or turning off repartitioning does indeed make the error go away.</div><div class="gmail_extra"><br><div class="gmail_quote">On Thu, Nov 2, 2017 at 8:51 AM, Mark Adams <span dir="ltr"><<a href="mailto:mfadams@lbl.gov" target="_blank">mfadams@lbl.gov</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><br><div class="gmail_extra"><br><div class="gmail_quote"><span class="">On Wed, Nov 1, 2017 at 9:36 PM, Randy Michael Churchill <span dir="ltr"><<a href="mailto:rchurchi@pppl.gov" target="_blank">rchurchi@pppl.gov</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">Doing some additional testing, the issue goes away when removing the gamg preconditioner line from the petsc.rc:<div>-pc_type gamg</div></div></blockquote><div><br></div></span><div>Yea, this is GAMG setup. </div><div> </div><div>This is the code.  findices is create with ISCreateStride, so it is sorted ...</div><div><br></div><div>Michael is repartitioning the coarse grids. Maybe we don't have a regression test with this...</div><div><br></div><div>I will try to reproduce this.</div><div><br></div><div>Michael: you can use hypre for now, or turn repartitioning off (eg, -fsa_fieldsplit_lambda_upper_<wbr>pc_gamg_repartition false), but I'm not sure this will fix this.<br></div><div><br></div><div>You don't have hypre parameters for all of your all of your solvers. I think 'boomeramg' is the default pc_hypre_type. That should be good enough for you.</div><div><br></div><div><br></div><div><div>    {</div><div>      IS       findices;</div><div>      PetscInt Istart,Iend;</div><div>      Mat      Pnew;</div><div><br></div><div>      ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);</div><div>#if defined PETSC_GAMG_USE_LOG</div><div>      ierr = PetscLogEventBegin(petsc_gamg_<wbr>setup_events[SET15],0,0,0,0);<wbr>CHKERRQ(ierr);</div><div>#endif</div><div>      ierr = ISCreateStride(comm,Iend-<wbr>Istart,Istart,1,&findices);<wbr>CHKERRQ(ierr);</div><div>      ierr = ISSetBlockSize(findices,f_bs);<wbr>CHKERRQ(ierr);</div><div>      ierr = MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);</div><div>      ierr = ISDestroy(&findices);CHKERRQ(<wbr>ierr);</div><div><br></div><div>#if defined PETSC_GAMG_USE_LOG</div><div>      ierr = PetscLogEventEnd(petsc_gamg_<wbr>setup_events[SET15],0,0,0,0);<wbr>CHKERRQ(ierr);</div><div>#endif</div><div>      ierr = MatDestroy(a_P_inout);CHKERRQ(<wbr>ierr);</div><div><br></div><div>      /* output - repartitioned */</div><div>      *a_P_inout = Pnew;</div><div>    }</div></div><div><div class="h5"><div><br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div><br></div></div><div class="gmail_extra"><div><div class="m_-8266480047035692828gmail-h5"><br><div class="gmail_quote">On Wed, Nov 1, 2017 at 8:23 PM, Hong <span dir="ltr"><<a href="mailto:hzhang@mcs.anl.gov" target="_blank">hzhang@mcs.anl.gov</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">Randy:</div><div class="gmail_quote">Thanks, I'll check it tomorrow.</div><div class="gmail_quote"><span class="m_-8266480047035692828gmail-m_6217542122699884941HOEnZb"><font color="#888888">Hong</font></span><div><div class="m_-8266480047035692828gmail-m_6217542122699884941h5"><br><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>OK, this might not be completely satisfactory, because it doesn't show the partitioning or how the matrix is created, but this reproduces the problem. I wrote out my matrix, Amat, from the larger simulation, and load it in this script. This must be run with MPI rank greater than 1. This may be some combination of my petsc.rc, because when I use the PetscInitialize with it, it throws the error, but when using default (PETSC_NULL_CHARACTER) it runs fine.</div><div><br></div><div class="gmail_extra"><div><div class="m_-8266480047035692828gmail-m_6217542122699884941m_3858210591366423390h5"><br><div class="gmail_quote">On Tue, Oct 31, 2017 at 9:58 AM, Hong <span dir="ltr"><<a href="mailto:hzhang@mcs.anl.gov" target="_blank">hzhang@mcs.anl.gov</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">Randy:</div><div class="gmail_quote">It could be a bug or a missing feature in our new MatCreateSubMatrix_MPIAIJ_Same<wbr>RowDist(). </div><div class="gmail_quote">It would be helpful if you can provide us a simple example that produces this example. </div><span class="m_-8266480047035692828gmail-m_6217542122699884941m_3858210591366423390m_3719836770456546225m_1134399805551081211gmail-HOEnZb"><font color="#888888"><div class="gmail_quote">Hong</div></font></span><div><div class="m_-8266480047035692828gmail-m_6217542122699884941m_3858210591366423390m_3719836770456546225m_1134399805551081211gmail-h5"><div class="gmail_quote"><br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">I'm running a Fortran code that was just changed over to using petsc 3.8 (previously petsc 3.7.6). An error was thrown during a KSPSetUp() call. The error is "unsorted iscol_local is not implemented yet" (see full error below). I tried to trace down the difference in the source files, but where the error occurs (MatCreateSubMatrix_MPIAIJ_Sam<wbr>eRowDist()) doesn't seem to have existed in v3.7.6, so I'm unsure how to compare. It seems the error is that the order of the columns locally are unsorted, though I don't think I specify a column order in the creation of the matrix:<div>     call MatCreate(this%comm,AA,ierr)</div><div>     call MatSetSizes(AA,npetscloc,npets<wbr>cloc,nreal,nreal,ierr)</div><div>     call MatSetType(AA,MATAIJ,ierr)</div><div>     call MatSetup(AA,ierr)</div><div>     call MatGetOwnershipRange(AA,low,hi<wbr>gh,ierr)</div><div>     allocate(d_nnz(npetscloc),o_n<wbr>nz(npetscloc))</div><div>     call getNNZ(grid,npetscloc,low,high<wbr>,d_nnz,o_nnz,this%xgc_petsc,nr<wbr>eal,ierr)</div><div>     call MatSeqAIJSetPreallocation(AA,P<wbr>ETSC_NULL_INTEGER,d_nnz,ierr)</div><div>     call MatMPIAIJSetPreallocation(AA,P<wbr>ETSC_NULL_INTEGER,d_nnz,PETSC_<wbr>NULL_INTEGER,o_nnz,ierr)</div><div>     deallocate(d_nnz,o_nnz)</div><div>     call MatSetOption(AA,MAT_IGNORE_OFF<wbr>_PROC_ENTRIES,PETSC_TRUE,ierr)</div><div>     call MatSetOption(AA,MAT_KEEP_NONZE<wbr>RO_PATTERN,PETSC_TRUE,ierr)</div><div>     call MatSetup(AA,ierr) <div>

<div><br></div><div><br></div><div>[62]PETSC ERROR: --------------------- Error Message ------------------------------<wbr>------------------------------<wbr>--</div><div>[62]PETSC ERROR: No support for this operation for this object type</div><div>[62]PETSC ERROR: unsorted iscol_local is not implemented yet</div><div>[62]PETSC ERROR: See <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank">http://www.mcs.anl.gov/petsc/d<wbr>ocumentation/faq.html</a> for trouble shooting.</div><div>[62]PETSC ERROR: Petsc Release Version 3.8.0, unknown[62]PETSC ERROR: #1 MatCreateSubMatrix_MPIAIJ_Same<wbr>RowDist() line 3418 in /global/u1/r/rchurchi/petsc/3.<wbr>8.0/src/mat/impls/aij/mpi/mpia<wbr>ij.c
</div><div>[62]PETSC ERROR: #2 MatCreateSubMatrix_MPIAIJ() line 3247 in /global/u1/r/rchurchi/petsc/3.<wbr>8.0/src/mat/impls/aij/mpi/mpia<wbr>ij.c
</div><div>[62]PETSC ERROR: #3 MatCreateSubMatrix() line 7872 in /global/u1/r/rchurchi/petsc/3.<wbr>8.0/src/mat/interface/matrix.c
</div><div>[62]PETSC ERROR: #4 PCGAMGCreateLevel_GAMG() line 383 in /global/u1/r/rchurchi/petsc/3.<wbr>8.0/src/ksp/pc/impls/gamg/gamg<wbr>.c
</div><div>[62]PETSC ERROR: #5 PCSetUp_GAMG() line 561 in /global/u1/r/rchurchi/petsc/3.<wbr>8.0/src/ksp/pc/impls/gamg/gamg<wbr>.c
</div><div>[62]PETSC ERROR: #6 PCSetUp() line 924 in /global/u1/r/rchurchi/petsc/3.<wbr>8.0/src/ksp/pc/interface/preco<wbr>n.c
</div><div>[62]PETSC ERROR: #7 KSPSetUp() line 378 in /global/u1/r/rchurchi/petsc/3.<wbr>8.0/src/ksp/ksp/interface/itfu<wbr>nc.c </div><span class="m_-8266480047035692828gmail-m_6217542122699884941m_3858210591366423390m_3719836770456546225m_1134399805551081211gmail-m_2087216136149895940gmail-HOEnZb"><font color="#888888"><div><br></div>-- <br><div class="m_-8266480047035692828gmail-m_6217542122699884941m_3858210591366423390m_3719836770456546225m_1134399805551081211gmail-m_2087216136149895940gmail-m_-4997150650000030250gmail_signature"><div dir="ltr">R. Michael Churchill</div></div>
</font></span></div></div></div>
</blockquote></div><br></div></div></div></div>
</blockquote></div><br><br clear="all"><div><br></div></div></div><span class="m_-8266480047035692828gmail-m_6217542122699884941m_3858210591366423390HOEnZb"><font color="#888888">-- <br><div class="m_-8266480047035692828gmail-m_6217542122699884941m_3858210591366423390m_3719836770456546225m_1134399805551081211gmail_signature"><div dir="ltr">R. Michael Churchill</div></div>
</font></span></div></div>
</blockquote></div></div></div><br></div></div>
</blockquote></div><br><br clear="all"><div><br></div></div></div><span class="m_-8266480047035692828gmail-HOEnZb"><font color="#888888">-- <br><div class="m_-8266480047035692828gmail-m_6217542122699884941gmail_signature"><div dir="ltr">R. Michael Churchill</div></div>
</font></span></div>
</blockquote></div></div></div><br></div></div>
</blockquote></div><br><br clear="all"><div><br></div>-- <br><div class="gmail_signature" data-smartmail="gmail_signature"><div dir="ltr">R. Michael Churchill</div></div>
</div>