[petsc-users] Using MUMPS as a preconditioner for KSPSolve()

Matthew Knepley knepley at gmail.com
Fri Dec 5 14:20:55 CST 2014


On Fri, Dec 5, 2014 at 2:11 PM, Evan Um <evanum at gmail.com> wrote:

> Dear PETSC Users,
>
> I tried to use a Cholesky factor (MUMPS results) as a preconditioner for
> KSPSolve(). An example code is pasted below. When the code runs, the log
> file indicates that Job=3 (i.e. backward/forward substitution) of MUMPS is
> called every time inside the loop. Is there anyway to avoid JOB=3 of MUMPS
> and use the factor as a pure preconditioner for the CG solver inside
> KSPSOLVE()? On my cluster, JOB=3 shows unexpected slow performance (see
> Vol72, Issue 35) and should be avoided. In advance, thanks for
>

What do you mean by the above sentence. When you Cholesky factor a matrix
A, you get

  A = L^T L

When you use this factorization as a preconditioner, you modify your
original problem to

  (L^T L)^{-1} A x = (L^T L)^{-1} b

When you use a Krylov method like CG, this means that at every iterate you
must apply
(L^T L){-1} which entails a forward and backward solve.

If you mean, you only want to apply (L^T L)^{-1} and not run CG, then use
-ksp_type preonly

  Matt


> your help.
>
> Regards,
> Evan
>
>
> 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);
> MatMumpsSetIcntl(F_dt, 4, 1); // Turn on MUMPS's log file.
> KSPSetType(ksp_fetd_dt, KSPCG);
> KSPSetTolerances(ksp_fetd_dt, 1e-9, 1.0e-50, 1.0e10, ksp_iter);
> for (int i=0; i<1000; i++) {
> // Create a new RHS vector B_dt
> KSPSolve(ksp_fetd_dt,B_dt,solution);
> // Output solution time=time2-time1;
> }
>



-- 
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/20141205/848d9ba8/attachment.html>


More information about the petsc-users mailing list