[petsc-users] KSP_DIVERGED_INDEFINITE_PC message even with PCFactorSetShiftType

Matthew Knepley knepley at gmail.com
Sun Sep 28 18:54:51 CDT 2014


On Sun, Sep 28, 2014 at 6:24 PM, Evan Um <evanum at gmail.com> wrote:

> Dear PETSC Users,
>
> I try to solve a diffusion problem in the time domain. Its system matrix
> is theoretically SPD. I use KSPCG solver. Direct MUMPS solver generates a
> preconditioner. Thus, within 1-2 iterations, a solution is expected to
> converge. Such convergence is observed in most test problems. However, in
> some test problems, I get non-convergence with a message:
> KSP_DIVERGED_INDEFINITE_PC. This problem is alleviated when I use
> PCFactorSetShiftType(pc_fetd_dt,MAT_SHIFT_POSITIVE_DEFINITE). However, the
> indefiniteness error message occurs later during time-stepping. I am aware
> of PETSC function PCFactorSetShiftAmount. However, before I try this, I
> couldn't find the information about a default shift amount when I use
> PCFactorSetShiftType(pc_fetd_dt,MAT_SHIFT_POSITIVE_DEFINITE). How can I
> find the default shift amount?
>
-ksp_view should tell you the shift amount


> More generally, how can I stably solve SPD problems that are close to
> indefiniteness using KSP options? To get the reference solution, I was able
> to solve the same problem above using SuiteSparse (serial) but its solution
> time was not practical.
>
This sounds like a discretization problem. However, you can always change
solvers to GMRES. If you are only doing 1-2
iterates, there is no advantage to CG (other than it minimizing a different
norm).

  Thanks,

    Matt


> In advance, thanks for your advice.
>
> Regards,
> Evan
>
> Sample Code:
>
> KSPCreate(PETSC_COMM_WORLD, &ksp_fetd_dt);
>
> KSPSetOperators(ksp_fetd_dt, A_dt, A_dt);
>
> KSPSetType (ksp_fetd_dt, KSPPREONLY);
>
> KSPGetPC(ksp_fetd_dt, &pc_fetd_dt);
>
> MatSetOption(A_dt, MAT_SPD, PETSC_TRUE);
>
> PCSetType(pc_fetd_dt, PCCHOLESKY);
>
> PCFactorSetMatSolverPackage(pc_fetd_dt, MATSOLVERMUMPS);
>
> PCFactorSetUpMatSolverPackage(pc_fetd_dt);
>
> PCFactorGetMatrix(pc_fetd_dt, &F_dt);
>
> KSPSetType(ksp_fetd_dt, KSPCG);
>
> PCFactorSetShiftType(pc_fetd_dt,MAT_SHIFT_POSITIVE_DEFINITE);
>
>
>
>
>



-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140928/519542d6/attachment.html>


More information about the petsc-users mailing list