<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">Randy:</div><div class="gmail_quote">I can reproduce it in maint branch (=v3.8) with option '-pc_gamg_mat_partitioning_type current'.</div><div class="gmail_quote">My petsc is built with parmetis, thus it uses '-pc_gamg_mat_partitioning_type parmetis' by default and works well. Replacing it with 'current', I am able to see the error -- will fix it.</div><div class="gmail_quote"><br></div><div class="gmail_quote">Meanwhile, you can use '-pc_gamg_mat_partitioning_type average' to get your work done.</div><div class="gmail_quote">I'll let you know the fix.</div><div class="gmail_quote"><br></div><div class="gmail_quote">Hong</div><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 used:<div>git checkout v3.8</div></div><div class="gmail_extra"><div><div class="gmail-h5"><br><div class="gmail_quote">On Thu, Nov 2, 2017 at 10:07 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">Mark,<div>I can reproduce this in an old branch, but not in current maint and master.</div><div>Which branch are you using to produce this error?</div><div><span class="gmail-m_-6588307739706285246HOEnZb"><font color="#888888">Hong</font></span><div><div class="gmail-m_-6588307739706285246h5"><br><div class="gmail_extra"><br><div class="gmail_quote">On Thu, Nov 2, 2017 at 9:28 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:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">I am able to reproduce this with snes ex56 with 2 processors and adding -pc_gamg_repartition true<div><br></div><div>I'm not sure how to fix it.</div><div><br></div><div><div>10:26 1 knepley/feature-plex-boxmesh-c<wbr>reate *= ~/Codes/petsc/src/snes/example<wbr>s/tutorials$ make PETSC_DIR=/Users/markadams/Cod<wbr>es/petsc PETSC_ARCH=arch-macosx-gnu-g runex       </div><div>/Users/markadams/Codes/petsc/a<wbr>rch-macosx-gnu-g/bin/mpiexec -n 2 ./ex56  -cells 2,2,1 -max_conv_its 3 -petscspace_order 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -ksp_converged_reason -snes_monitor_short -ksp_monitor_short -snes_converged_reason -use_mat_nearnullspace true -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_estei<wbr>g 0,0.05,0,1.05 -mg_levels_pc_type jacobi -petscpartitioner_type simple -mat_block_size 3 -matrap 0 -matptap_scalable -ex56_dm_view -run_type 1 -pc_gamg_repartition true</div><div>[0] 27 global equations, 9 vertices</div><div>[0] 27 equations in vector, 9 vertices</div><div>  0 SNES Function norm 122.396 </div><div>    0 KSP Residual norm 122.396 </div><div>    1 KSP Residual norm 20.4696 </div><div>    2 KSP Residual norm 3.95009 </div><div>    3 KSP Residual norm 0.176181 </div><div>    4 KSP Residual norm 0.0208781 </div><div>    5 KSP Residual norm 0.00278873 </div><div>    6 KSP Residual norm 0.000482741 </div><div>    7 KSP Residual norm 4.68085e-05 </div><div>    8 KSP Residual norm 5.42381e-06 </div><div>    9 KSP Residual norm 5.12785e-07 </div><div>   10 KSP Residual norm 2.60389e-08 </div><div>   11 KSP Residual norm 4.96201e-09 </div><div>   12 KSP Residual norm 1.989e-10 </div><div>  Linear solve converged due to CONVERGED_RTOL iterations 12</div><div>  1 SNES Function norm 1.990e-10 </div><div>Nonlinear solve converged due to CONVERGED_FNORM_RELATIVE iterations 1</div><div>DM Object: Mesh (ex56_) 2 MPI processes</div><div>  type: plex</div><div>Mesh in 3 dimensions:</div><div>  0-cells: 12 12</div><div>  1-cells: 20 20</div><div>  2-cells: 11 11</div><div>  3-cells: 2 2</div><div>Labels:</div><div>  boundary: 1 strata with value/size (1 (39))</div><div>  Face Sets: 5 strata with value/size (1 (2), 2 (2), 3 (2), 5 (1), 6 (1))</div><div>  marker: 1 strata with value/size (1 (27))</div><div>  depth: 4 strata with value/size (0 (12), 1 (20), 2 (11), 3 (2))</div><div>[0] 441 global equations, 147 vertices</div><div>[0] 441 equations in vector, 147 vertices</div><div>  0 SNES Function norm 49.7106 </div><div>    0 KSP Residual norm 49.7106 </div><div>    1 KSP Residual norm 12.9252 </div><div>    2 KSP Residual norm 2.38019 </div><div>    3 KSP Residual norm 0.426307 </div><div>    4 KSP Residual norm 0.0692155 </div><div>    5 KSP Residual norm 0.0123092 </div><div>    6 KSP Residual norm 0.00184874 </div><div>    7 KSP Residual norm 0.000320761 </div><div>    8 KSP Residual norm 5.48957e-05 </div><div>    9 KSP Residual norm 9.90089e-06 </div><div>   10 KSP Residual norm 1.5127e-06 </div><div>   11 KSP Residual norm 2.82192e-07 </div><div>   12 KSP Residual norm 4.62364e-08 </div><div>   13 KSP Residual norm 7.99573e-09 </div><div>   14 KSP Residual norm 1.3028e-09 </div><div>   15 KSP Residual norm 2.174e-10 </div><div>  Linear solve converged due to CONVERGED_RTOL iterations 15</div><div>  1 SNES Function norm 2.174e-10 </div><div>Nonlinear solve converged due to CONVERGED_FNORM_RELATIVE iterations 1</div><div>DM Object: Mesh (ex56_) 2 MPI processes</div><div>  type: plex</div><div>Mesh in 3 dimensions:</div><div>  0-cells: 45 45</div><div>  1-cells: 96 96</div><div>  2-cells: 68 68</div><div>  3-cells: 16 16</div><div>Labels:</div><div>  marker: 1 strata with value/size (1 (129))</div><div>  Face Sets: 5 strata with value/size (1 (18), 2 (18), 3 (18), 5 (9), 6 (9))</div><div>  boundary: 1 strata with value/size (1 (141))</div><div>  depth: 4 strata with value/size (0 (45), 1 (96), 2 (68), 3 (16))</div><div>[0] 4725 global equations, 1575 vertices</div><div>[0] 4725 equations in vector, 1575 vertices</div><div>  0 SNES Function norm 17.9091 </div><div>[0]PETSC ERROR: --------------------- Error Message ------------------------------<wbr>------------------------------<wbr>--</div><div>[0]PETSC ERROR: No support for this operation for this object type</div><div>[0]PETSC ERROR: unsorted iscol_local is not implemented yet</div><div>[1]PETSC ERROR: --------------------- Error Message ------------------------------<wbr>------------------------------<wbr>--</div><div>[1]PETSC ERROR: No support for this operation for this object type</div><div>[1]PETSC ERROR: unsorted iscol_local is not implemented yet</div></div><div><br></div></div><div class="gmail-m_-6588307739706285246m_4896846594077993743HOEnZb"><div class="gmail-m_-6588307739706285246m_4896846594077993743h5"><div class="gmail_extra"><br><div class="gmail_quote">On Thu, Nov 2, 2017 at 9: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:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><br><div class="gmail_extra"><br><div class="gmail_quote"><span>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_p<wbr>c_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);C<wbr>HKERRQ(ierr);</div><div>#endif</div><div>      ierr = ISCreateStride(comm,Iend-Istar<wbr>t,Istart,1,&findices);CHKERRQ(<wbr>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(i<wbr>err);</div><div><br></div><div>#if defined PETSC_GAMG_USE_LOG</div><div>      ierr = PetscLogEventEnd(petsc_gamg_se<wbr>tup_events[SET15],0,0,0,0);CHK<wbr>ERRQ(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="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761h5"><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="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761m_6139592594514371676gmail-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="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761m_6139592594514371676gmail-m_6217542122699884941HOEnZb"><font color="#888888">Hong</font></span><div><div class="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761m_6139592594514371676gmail-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="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761m_6139592594514371676gmail-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="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761m_6139592594514371676gmail-m_6217542122699884941m_3858210591366423390m_3719836770456546225m_1134399805551081211gmail-HOEnZb"><font color="#888888"><div class="gmail_quote">Hong</div></font></span><div><div class="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761m_6139592594514371676gmail-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="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761m_6139592594514371676gmail-m_6217542122699884941m_3858210591366423390m_3719836770456546225m_1134399805551081211gmail-m_2087216136149895940gmail-HOEnZb"><font color="#888888"><div><br></div>-- <br><div class="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761m_6139592594514371676gmail-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="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761m_6139592594514371676gmail-m_6217542122699884941m_3858210591366423390HOEnZb"><font color="#888888">-- <br><div class="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761m_6139592594514371676gmail-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="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761m_6139592594514371676gmail-HOEnZb"><font color="#888888">-- <br><div class="gmail-m_-6588307739706285246m_4896846594077993743m_-130663177837359761m_6139592594514371676gmail-m_6217542122699884941gmail_signature"><div dir="ltr">R. Michael Churchill</div></div>
</font></span></div>
</blockquote></div></div></div><br></div></div>
</blockquote></div><br></div>
</div></div></blockquote></div><br></div></div></div></div></div>
</blockquote></div><br><br clear="all"><div><br></div></div></div><span class="gmail-HOEnZb"><font color="#888888">-- <br><div class="gmail-m_-6588307739706285246gmail_signature"><div dir="ltr">R. Michael Churchill</div></div>
</font></span></div>
</blockquote></div><br></div></div>