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

Evan Um evanum at gmail.com
Fri Dec 5 14:11:34 CST 2014


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 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;
}
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20141205/075fac75/attachment.html>


More information about the petsc-users mailing list