[petsc-users] unsorted local columns in 3.8?

Mark Adams mfadams at lbl.gov
Thu Nov 2 09:28:47 CDT 2017


I am able to reproduce this with snes ex56 with 2 processors and adding
-pc_gamg_repartition true

I'm not sure how to fix it.

10:26 1 knepley/feature-plex-boxmesh-create *=
~/Codes/petsc/src/snes/examples/tutorials$ make
PETSC_DIR=/Users/markadams/Codes/petsc PETSC_ARCH=arch-macosx-gnu-g runex

/Users/markadams/Codes/petsc/arch-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_esteig
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
[0] 27 global equations, 9 vertices
[0] 27 equations in vector, 9 vertices
  0 SNES Function norm 122.396
    0 KSP Residual norm 122.396
    1 KSP Residual norm 20.4696
    2 KSP Residual norm 3.95009
    3 KSP Residual norm 0.176181
    4 KSP Residual norm 0.0208781
    5 KSP Residual norm 0.00278873
    6 KSP Residual norm 0.000482741
    7 KSP Residual norm 4.68085e-05
    8 KSP Residual norm 5.42381e-06
    9 KSP Residual norm 5.12785e-07
   10 KSP Residual norm 2.60389e-08
   11 KSP Residual norm 4.96201e-09
   12 KSP Residual norm 1.989e-10
  Linear solve converged due to CONVERGED_RTOL iterations 12
  1 SNES Function norm 1.990e-10
Nonlinear solve converged due to CONVERGED_FNORM_RELATIVE iterations 1
DM Object: Mesh (ex56_) 2 MPI processes
  type: plex
Mesh in 3 dimensions:
  0-cells: 12 12
  1-cells: 20 20
  2-cells: 11 11
  3-cells: 2 2
Labels:
  boundary: 1 strata with value/size (1 (39))
  Face Sets: 5 strata with value/size (1 (2), 2 (2), 3 (2), 5 (1), 6 (1))
  marker: 1 strata with value/size (1 (27))
  depth: 4 strata with value/size (0 (12), 1 (20), 2 (11), 3 (2))
[0] 441 global equations, 147 vertices
[0] 441 equations in vector, 147 vertices
  0 SNES Function norm 49.7106
    0 KSP Residual norm 49.7106
    1 KSP Residual norm 12.9252
    2 KSP Residual norm 2.38019
    3 KSP Residual norm 0.426307
    4 KSP Residual norm 0.0692155
    5 KSP Residual norm 0.0123092
    6 KSP Residual norm 0.00184874
    7 KSP Residual norm 0.000320761
    8 KSP Residual norm 5.48957e-05
    9 KSP Residual norm 9.90089e-06
   10 KSP Residual norm 1.5127e-06
   11 KSP Residual norm 2.82192e-07
   12 KSP Residual norm 4.62364e-08
   13 KSP Residual norm 7.99573e-09
   14 KSP Residual norm 1.3028e-09
   15 KSP Residual norm 2.174e-10
  Linear solve converged due to CONVERGED_RTOL iterations 15
  1 SNES Function norm 2.174e-10
Nonlinear solve converged due to CONVERGED_FNORM_RELATIVE iterations 1
DM Object: Mesh (ex56_) 2 MPI processes
  type: plex
Mesh in 3 dimensions:
  0-cells: 45 45
  1-cells: 96 96
  2-cells: 68 68
  3-cells: 16 16
Labels:
  marker: 1 strata with value/size (1 (129))
  Face Sets: 5 strata with value/size (1 (18), 2 (18), 3 (18), 5 (9), 6 (9))
  boundary: 1 strata with value/size (1 (141))
  depth: 4 strata with value/size (0 (45), 1 (96), 2 (68), 3 (16))
[0] 4725 global equations, 1575 vertices
[0] 4725 equations in vector, 1575 vertices
  0 SNES Function norm 17.9091
[0]PETSC ERROR: --------------------- Error Message
--------------------------------------------------------------
[0]PETSC ERROR: No support for this operation for this object type
[0]PETSC ERROR: unsorted iscol_local is not implemented yet
[1]PETSC ERROR: --------------------- Error Message
--------------------------------------------------------------
[1]PETSC ERROR: No support for this operation for this object type
[1]PETSC ERROR: unsorted iscol_local is not implemented yet


On Thu, Nov 2, 2017 at 9:51 AM, Mark Adams <mfadams at lbl.gov> wrote:

>
>
> On Wed, Nov 1, 2017 at 9:36 PM, Randy Michael Churchill <rchurchi at pppl.gov
> > wrote:
>
>> Doing some additional testing, the issue goes away when removing the gamg
>> preconditioner line from the petsc.rc:
>> -pc_type gamg
>>
>
> Yea, this is GAMG setup.
>
> This is the code.  findices is create with ISCreateStride, so it is sorted
> ...
>
> Michael is repartitioning the coarse grids. Maybe we don't have a
> regression test with this...
>
> I will try to reproduce this.
>
> Michael: you can use hypre for now, or turn repartitioning off (eg,
> -fsa_fieldsplit_lambda_upper_pc_gamg_repartition false), but I'm not sure
> this will fix this.
>
> 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.
>
>
>     {
>       IS       findices;
>       PetscInt Istart,Iend;
>       Mat      Pnew;
>
>       ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr);
> #if defined PETSC_GAMG_USE_LOG
>       ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);
> CHKERRQ(ierr);
> #endif
>       ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);
> CHKERRQ(ierr);
>       ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr);
>       ierr = MatCreateSubMatrix(Pold, findices, new_eq_indices,
> MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr);
>       ierr = ISDestroy(&findices);CHKERRQ(ierr);
>
> #if defined PETSC_GAMG_USE_LOG
>       ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);
> CHKERRQ(ierr);
> #endif
>       ierr = MatDestroy(a_P_inout);CHKERRQ(ierr);
>
>       /* output - repartitioned */
>       *a_P_inout = Pnew;
>     }
>
>
>>
>> On Wed, Nov 1, 2017 at 8:23 PM, Hong <hzhang at mcs.anl.gov> wrote:
>>
>>> Randy:
>>> Thanks, I'll check it tomorrow.
>>> Hong
>>>
>>> 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.
>>>>
>>>>
>>>> On Tue, Oct 31, 2017 at 9:58 AM, Hong <hzhang at mcs.anl.gov> wrote:
>>>>
>>>>> Randy:
>>>>> It could be a bug or a missing feature in our new
>>>>> MatCreateSubMatrix_MPIAIJ_SameRowDist().
>>>>> It would be helpful if you can provide us a simple example that
>>>>> produces this example.
>>>>> Hong
>>>>>
>>>>> 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_SameRowDist()) 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:
>>>>>>      call MatCreate(this%comm,AA,ierr)
>>>>>>      call MatSetSizes(AA,npetscloc,npetscloc,nreal,nreal,ierr)
>>>>>>      call MatSetType(AA,MATAIJ,ierr)
>>>>>>      call MatSetup(AA,ierr)
>>>>>>      call MatGetOwnershipRange(AA,low,high,ierr)
>>>>>>      allocate(d_nnz(npetscloc),o_nnz(npetscloc))
>>>>>>      call getNNZ(grid,npetscloc,low,high
>>>>>> ,d_nnz,o_nnz,this%xgc_petsc,nreal,ierr)
>>>>>>      call MatSeqAIJSetPreallocation(AA,PETSC_NULL_INTEGER,d_nnz,ierr)
>>>>>>      call MatMPIAIJSetPreallocation(AA,P
>>>>>> ETSC_NULL_INTEGER,d_nnz,PETSC_NULL_INTEGER,o_nnz,ierr)
>>>>>>      deallocate(d_nnz,o_nnz)
>>>>>>      call MatSetOption(AA,MAT_IGNORE_OFF
>>>>>> _PROC_ENTRIES,PETSC_TRUE,ierr)
>>>>>>      call MatSetOption(AA,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE,ierr)
>>>>>>      call MatSetup(AA,ierr)
>>>>>>
>>>>>>
>>>>>> [62]PETSC ERROR: --------------------- Error Message
>>>>>> --------------------------------------------------------------
>>>>>> [62]PETSC ERROR: No support for this operation for this object type
>>>>>> [62]PETSC ERROR: unsorted iscol_local is not implemented yet
>>>>>> [62]PETSC ERROR: See http://www.mcs.anl.gov/petsc/d
>>>>>> ocumentation/faq.html for trouble shooting.
>>>>>> [62]PETSC ERROR: Petsc Release Version 3.8.0, unknown[62]PETSC ERROR:
>>>>>> #1 MatCreateSubMatrix_MPIAIJ_SameRowDist() line 3418 in
>>>>>> /global/u1/r/rchurchi/petsc/3.8.0/src/mat/impls/aij/mpi/mpiaij.c
>>>>>> [62]PETSC ERROR: #2 MatCreateSubMatrix_MPIAIJ() line 3247 in
>>>>>> /global/u1/r/rchurchi/petsc/3.8.0/src/mat/impls/aij/mpi/mpiaij.c
>>>>>> [62]PETSC ERROR: #3 MatCreateSubMatrix() line 7872 in
>>>>>> /global/u1/r/rchurchi/petsc/3.8.0/src/mat/interface/matrix.c
>>>>>> [62]PETSC ERROR: #4 PCGAMGCreateLevel_GAMG() line 383 in
>>>>>> /global/u1/r/rchurchi/petsc/3.8.0/src/ksp/pc/impls/gamg/gamg.c
>>>>>> [62]PETSC ERROR: #5 PCSetUp_GAMG() line 561 in
>>>>>> /global/u1/r/rchurchi/petsc/3.8.0/src/ksp/pc/impls/gamg/gamg.c
>>>>>> [62]PETSC ERROR: #6 PCSetUp() line 924 in
>>>>>> /global/u1/r/rchurchi/petsc/3.8.0/src/ksp/pc/interface/precon.c
>>>>>> [62]PETSC ERROR: #7 KSPSetUp() line 378 in
>>>>>> /global/u1/r/rchurchi/petsc/3.8.0/src/ksp/ksp/interface/itfunc.c
>>>>>>
>>>>>> --
>>>>>> R. Michael Churchill
>>>>>>
>>>>>
>>>>>
>>>>
>>>>
>>>> --
>>>> R. Michael Churchill
>>>>
>>>
>>>
>>
>>
>> --
>> R. Michael Churchill
>>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20171102/a399ee2b/attachment.html>


More information about the petsc-users mailing list