[petsc-users] Petsc has generated inconsistent data!

Barry Smith bsmith at mcs.anl.gov
Fri Aug 26 16:17:56 CDT 2011


  Are you sure that is the entire error message. It should print the routine and the line number where this happens. 

   Likely it is at

  do {
    ierr = VecDot(R,RP,&rho);CHKERRQ(ierr);       /*   rho <- (r,rp)      */
    beta = (rho/rhoold) * (alpha/omegaold);
    ierr = VecAXPBYPCZ(P,1.0,-omegaold*beta,beta,R,V);CHKERRQ(ierr);  /* p <- r - omega * beta* v + beta * p */
    ierr = KSP_PCApplyBAorAB(ksp,P,V,T);CHKERRQ(ierr);  /*   v <- K p           */
    ierr = VecDot(V,RP,&d1);CHKERRQ(ierr);
    if (d1 == 0.0) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"Divide by zero");
    alpha = rho / d1;                 /*   a <- rho / (v,rp)  */

  Which means bi-cg-stab has broken down. You'll need to consult references to Bi-CG-stab to see why this might happen (it can happen while GMRES is happy). It may be KSPBCGSL can proceed past this point with a problem values for 
   Options Database Keys:
+  -ksp_bcgsl_ell <ell> Number of Krylov search directions
-  -ksp_bcgsl_cxpol Use a convex function of the MR and OR polynomials after the BiCG step
-  -ksp_bcgsl_xres <res> Threshold used to decide when to refresh computed residuals

but most likely the preconditioner or matrix is bogus in some way since I think Bi-CG-stab rarely breaks down in practice.


  Barry



On Aug 26, 2011, at 4:05 PM, Dominik Szczerba wrote:

> When solving my linear system with -ksp_type bcgs I get:
> 
>  0 KSP preconditioned resid norm 1.166190378969e+01 true resid norm
> 1.166190378969e+01 ||Ae||/||Ax|| 1.000000000000e+00
>  1 KSP preconditioned resid norm 5.658835826231e-01 true resid norm
> 5.658835826231e-01 ||Ae||/||Ax|| 4.852411688762e-02
> [1]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [1]PETSC ERROR: Petsc has generated inconsistent data!
> [1]PETSC ERROR: Divide by zero!
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: Petsc Release Version 3.1.0, Patch 8, Thu Mar 17
> 13:37:48 CDT 2011
> [1]PETSC ERROR: See docs/changes/index.html for recent updates.
> [1]PETSC ERROR: [3]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [3]PETSC ERROR: Petsc has generated inconsistent data!
> [3]PETSC ERROR: Divide by zero!
> 
> PETSc Option Table entries:
> -ksp_converged_reason
> -ksp_monitor
> -ksp_monitor_true_residual
> -ksp_norm_type unpreconditioned
> -ksp_right_pc
> -ksp_rtol 1e-3
> -ksp_type bcgs
> -ksp_view
> -log_summary
> -pc_type jacobi
> 
> When solving the same system with GMRES all works fine. This is a
> simple test diffusion problem. How can I find out what the problem is?
> 
> Dominik



More information about the petsc-users mailing list