[petsc-users] PetscLogFlops Error

Barry Smith bsmith at petsc.dev
Sun Jul 19 00:26:22 CDT 2020



> On Jul 18, 2020, at 11:27 PM, Karl Yang <y.juntao at hotmail.com> wrote:
> 
> 
> Hello,
> 
> I was using FGMRES solver with fieldsplit preconditioner for solving a fluid equation.
> I have the code ran for coarser meshes and it works fine. But when I got it to run for very fine meshes I encountered the PetscLogFlops error.
> 
> The following is the output of error and the output from kspview.
> 
> [7]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> [7]PETSC ERROR: Argument out of range
> [7]PETSC ERROR: Cannot log negative flops
> [7]PETSC ERROR: See https://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> [7]PETSC ERROR: Petsc Release Version 3.12.5, Mar, 29, 2020
> [7]PETSC ERROR: ./testSolve on a arch-linux2-c-debug named hsw222 by yjuntao Sat Jul 18 20:30:05 2020
> [7]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-mpich --download-fblaslapack --with-cuda
> [7]PETSC ERROR: #1 PetscLogFlops() line 232 in /usr/local/petsc/petsc-3.12.5/include/petsclog.h
> [7]PETSC ERROR: #2 MatSolve_SeqAIJ() line 3200 in /usr/local/petsc/petsc-3.12.5/src/mat/impls/aij/seq/aijfact.c
> [7]PETSC ERROR: #3 MatSolve() line 3290 in /usr/local/petsc/petsc-3.12.5/src/mat/interface/matrix.c
> [7]PETSC ERROR: #4 PCApply_LU() line 177 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/impls/factor/lu/lu.c
> [7]PETSC ERROR: #5 PCApply() line 444 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/interface/precon.c
> [7]PETSC ERROR: #6 KSP_PCApply() line 281 in /usr/local/petsc/petsc-3.12.5/include/petsc/private/kspimpl.h
> [7]PETSC ERROR: #7 KSPInitialResidual() line 65 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itres.c
> [7]PETSC ERROR: #8 KSPSolve_GMRES() line 236 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/impls/gmres/gmres.c
> [7]PETSC ERROR: #9 KSPSolve() line 760 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itfunc.c
> [7]PETSC ERROR: #10 PCApply_FieldSplit_Schur() line 1180 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/impls/fieldsplit/fieldsplit.c
> [7]PETSC ERROR: #11 PCApply() line 444 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/interface/precon.c
> [7]PETSC ERROR: #12 KSP_PCApply() line 281 in /usr/local/petsc/petsc-3.12.5/include/petsc/private/kspimpl.h
> [7]PETSC ERROR: #13 KSPFGMRESCycle() line 166 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/impls/gmres/fgmres/fgmres.c
> [7]PETSC ERROR: #14 KSPSolve_FGMRES() line 291 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/impls/gmres/fgmres/fgmres.c
> [7]PETSC ERROR: #15 KSPSolve() line 760 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itfunc.c
> 

PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
{
  Mat_SeqAIJ        *a    = (Mat_SeqAIJ*)A->data;
  IS                iscol = a->col,isrow = a->row;
  PetscErrorCode    ierr;
  PetscInt          i,n=A->rmap->n,*vi,*ai=a->i,*aj=a->j,*adiag = a->diag,nz;
  const PetscInt    *rout,*cout,*r,*c;
  PetscScalar       *x,*tmp,sum;
  const PetscScalar *b;
  const MatScalar   *aa = a->a,*v;

  PetscFunctionBegin;
  if (!n) PetscFunctionReturn(0);

......

  ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
  PetscFunctionReturn(0);

 Could maybe be integer overflow where a->nz fits but 2*a->nz doesn't fit so it results in a negative integer that is then converted to a negative floating point number.

 It should have a 2.0*a->nz but we must have missed this place in the code. We can fix this.

> There is some issues during kspview if PCLSC is used as shown below. But the code could be ran if I do not show ksp details.
> In order to show the whole ksp setup, I also attached another kspview output without PCLSC.
> ======================kspview output ==============================
> KSP Object: 1 MPI processes
>   type: fgmres
>     restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
>     happy breakdown tolerance 1e-30
>   maximum iterations=10000, initial guess is zero
>   tolerances:  relative=1e-08, absolute=1e-50, divergence=10000.
>   right preconditioning
>   using UNPRECONDITIONED norm type for convergence test
> PC Object: 1 MPI processes
>   type: fieldsplit
>     FieldSplit with Schur preconditioner, blocksize = 1, factorization FULL
>     Preconditioner for the Schur complement formed from S itself
>     Split info:
>     Split number 0 Defined by IS
>     Split number 1 Defined by IS
>     KSP solver for A00 block
>       KSP Object: (fieldsplit_0_) 1 MPI processes
>         type: gmres
>           restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
>           happy breakdown tolerance 1e-30
>         maximum iterations=10000, initial guess is zero
>         tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
>         left preconditioning
>         using DEFAULT norm type for convergence test
>       PC Object: (fieldsplit_0_) 1 MPI processes
>         type: lu
>         PC has not been set up so information may be incomplete
>           out-of-place factorization
>           tolerance for zero pivot 2.22045e-14
>           matrix ordering: nd
>         linear system matrix = precond matrix:
>         Mat Object: (fieldsplit_0_) 1 MPI processes
>           type: seqaijcusparse
>           rows=8396802, cols=8396802
>           total: nonzeros=75440146, allocated nonzeros=75440146
>           total number of mallocs used during MatSetValues calls=0
>             not using I-node routines
>     KSP solver for S = A11 - A10 inv(A00) A01
>       KSP Object: (fieldsplit_1_) 1 MPI processes
>         type: gmres
>           restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
>           happy breakdown tolerance 1e-30
>         maximum iterations=10000, initial guess is zero
>         tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
>         left preconditioning
>         using DEFAULT norm type for convergence test
>       PC Object: (fieldsplit_1_) 1 MPI processes
>         type: lsc
>         PC has not been set up so information may be incomplete
> [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> [0]PETSC ERROR: Null argument, when expecting valid pointer
> [0]PETSC ERROR: Null Object: Parameter # 1
> [0]PETSC ERROR: See https://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> [0]PETSC ERROR: Petsc Release Version 3.12.5, Mar, 29, 2020
> [0]PETSC ERROR: ./testSolve on a arch-linux2-c-debug named hsw221 by yjuntao Sat Jul 18 21:08:42 2020
> [0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-mpich --download-fblaslapack --with-cuda
> [0]PETSC ERROR: #1 KSPView() line 111 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itcreate.c
> [0]PETSC ERROR: #2 PCView_LSC() line 142 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/impls/lsc/lsc.c
> [0]PETSC ERROR: #3 PCView() line 1588 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/interface/precon.c
> [0]PETSC ERROR: #4 KSPView() line 217 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itcreate.c
> [0]PETSC ERROR: #5 PCView_FieldSplit_Schur() line 238 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/impls/fieldsplit/fieldsplit.c
> [0]PETSC ERROR: #6 PCView() line 1588 in /usr/local/petsc/petsc-3.12.5/src/ksp/pc/interface/precon.c
> [0]PETSC ERROR: #7 KSPView() line 217 in /usr/local/petsc/petsc-3.12.5/src/ksp/ksp/interface/itcreate.c

   This is something we can fix. Since the LSC has not been setup yet it should skip the KSPView on the KSP that has not yet been created.

static PetscErrorCode PCView_LSC(PC pc,PetscViewer viewer)
{
  PC_LSC         *jac = (PC_LSC*)pc->data;
  PetscErrorCode ierr;
  PetscBool      iascii;

  PetscFunctionBegin;
  ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
  if (iascii) {
    ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
    ierr = KSPView(jac->kspL,viewer);CHKERRQ(ierr);
    ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
  }
  PetscFunctionReturn(0);
}

> 
> =====================KSP without PCLSC======================================
> KSP Object: 1 MPI processes
>   type: fgmres
>     restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
>     happy breakdown tolerance 1e-30
>   maximum iterations=10000, initial guess is zero
>   tolerances:  relative=1e-08, absolute=1e-50, divergence=10000.
>   right preconditioning
>   using UNPRECONDITIONED norm type for convergence test
> PC Object: 1 MPI processes
>   type: fieldsplit
>     FieldSplit with Schur preconditioner, blocksize = 1, factorization FULL
>     Preconditioner for the Schur complement formed from S itself
>     Split info:
>     Split number 0 Defined by IS
>     Split number 1 Defined by IS
>     KSP solver for A00 block
>       KSP Object: (fieldsplit_0_) 1 MPI processes
>         type: gmres
>           restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
>           happy breakdown tolerance 1e-30
>         maximum iterations=10000, initial guess is zero
>         tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
>         left preconditioning
>         using DEFAULT norm type for convergence test
>       PC Object: (fieldsplit_0_) 1 MPI processes
>         type: lu
>         PC has not been set up so information may be incomplete
>           out-of-place factorization
>           tolerance for zero pivot 2.22045e-14
>           matrix ordering: nd
>         linear system matrix = precond matrix:
>         Mat Object: (fieldsplit_0_) 1 MPI processes
>           type: seqaijcusparse
>           rows=8396802, cols=8396802
>           total: nonzeros=75440146, allocated nonzeros=75440146
>           total number of mallocs used during MatSetValues calls=0
>             not using I-node routines
>     KSP solver for S = A11 - A10 inv(A00) A01
>       KSP Object: (fieldsplit_1_) 1 MPI processes
>         type: gmres
>           restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
>           happy breakdown tolerance 1e-30
>         maximum iterations=10000, initial guess is zero
>         tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
>         left preconditioning
>         using DEFAULT norm type for convergence test
>       PC Object: (fieldsplit_1_) 1 MPI processes
>         type: none
>         PC has not been set up so information may be incomplete
>         linear system matrix = precond matrix:
>         Mat Object: (fieldsplit_1_) 1 MPI processes
>           type: schurcomplement
>           rows=1050625, cols=1050625
>             Schur complement A11 - A10 inv(A00) A01
>             A11
>               Mat Object: (fieldsplit_1_) 1 MPI processes
>                 type: seqaijcusparse
>                 rows=1050625, cols=1050625
>                 total: nonzeros=0, allocated nonzeros=0
>                 total number of mallocs used during MatSetValues calls=0
>                   not using I-node routines
>             A10
>               Mat Object: 1 MPI processes
>                 type: seqaijcusparse
>                 rows=1050625, cols=8396802
>                 total: nonzeros=44060674, allocated nonzeros=44060674
>                 total number of mallocs used during MatSetValues calls=0
>                   not using I-node routines
>             KSP of A00
>               KSP Object: (fieldsplit_0_) 1 MPI processes
>                 type: gmres
>                   restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
>                   happy breakdown tolerance 1e-30
>                 maximum iterations=10000, initial guess is zero
>                 tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
>                 left preconditioning
>                 using DEFAULT norm type for convergence test
>               PC Object: (fieldsplit_0_) 1 MPI processes
>                 type: lu
>                 PC has not been set up so information may be incomplete
>                   out-of-place factorization
>                   tolerance for zero pivot 2.22045e-14
>                   matrix ordering: nd
>                 linear system matrix = precond matrix:
>                 Mat Object: (fieldsplit_0_) 1 MPI processes
>                   type: seqaijcusparse
>                   rows=8396802, cols=8396802
>                   total: nonzeros=75440146, allocated nonzeros=75440146
>                   total number of mallocs used during MatSetValues calls=0
>                     not using I-node routines
>             A01
>               Mat Object: 1 MPI processes
>                 type: seqaijcusparse
>                 rows=8396802, cols=1050625
>                 total: nonzeros=43995146, allocated nonzeros=43995146
>                 total number of mallocs used during MatSetValues calls=0
>                   not using I-node routines
>   linear system matrix = precond matrix:
>   Mat Object: 1 MPI processes
>     type: seqaijcusparse
>     rows=9447427, cols=9447427
>     total: nonzeros=163495966, allocated nonzeros=163643398
>     total number of mallocs used during MatSetValues calls=0
>       not using I-node routines
> 
> Regards
> JT



More information about the petsc-users mailing list