[petsc-users] unsorted local columns in 3.8?

Mark Adams mfadams at lbl.gov
Thu Nov 2 08:51:25 CDT 2017


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,nr
>>>>> eal,ierr)
>>>>>      call MatSeqAIJSetPreallocation(AA,PETSC_NULL_INTEGER,d_nnz,ierr)
>>>>>      call MatMPIAIJSetPreallocation(AA,PETSC_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/353e499d/attachment.html>


More information about the petsc-users mailing list