[petsc-dev] (PETSc-3.1-p6 & SUPERLU_dist) + ("-ksp_view" CLI option & KSPSetOptionsPrefix()/PCSetOptionsPrefix() ) = bug?
Satish Balay
balay at mcs.anl.gov
Wed Dec 15 13:46:29 CST 2010
sorry missed the attachment.. its attached now.
satish
On Wed, 15 Dec 2010, Satish Balay wrote:
> Can you try the attached patch - and see if this goes away?
>
> patch -Np1 < superlu_dist.patch
>
> Satish
>
> On Wed, 15 Dec 2010, Filippo Spiga wrote:
>
> > Dear all,
> > Sorry for the subject of this email, I had no better ideas to introduce
> > the topic (-: I spent a couple of hours yesterday night trying to figure out
> > about a strange behavior of PETSc-3.1-p6 using SUPERLU_dist. This problem
> > arises if I combine the usage of KSPSetOptionsPrefix()/PCSetOptionsPrefix()
> > and "-ksp_view" option. This is the error that I obtain:
> >
> > KSP Object:(lop_)
> > type: gmres
> > GMRES: restart=30, using Classical (unmodified) Gram-Schmidt
> > Orthogonalization with no iterative refinement
> > GMRES: happy breakdown tolerance 1e-30
> > maximum iterations=1000, initial guess is zero
> > tolerances: relative=1e-14, absolute=1e-14, divergence=10000
> > left preconditioning
> > using PRECONDITIONED norm type for convergence test
> > PC Object:(lop_)
> > type: lu
> > LU: out-of-place factorization
> > tolerance for zero pivot 1e-12
> > matrix ordering: nd
> > factor fill ratio given 0, needed 0
> > Factored matrix follows:
> > Matrix Object:
> > type=seqaij, rows=242, cols=242
> > package used to perform factorization: superlu_dist
> > total: nonzeros=0, allocated nonzeros=242
> > SuperLU_DIST run parameters:
> > Process grid nprow 1 x npcol 1
> > Equilibrate matrix TRUE
> > Matrix input mode 0
> > Replace tiny pivots TRUE
> > Use iterative refinement FALSE
> > Processors in row 1 col partition 1
> > Row permutation LargeDiag
> > [0|09:06:54]: unknown: MatFactorInfo_SuperLU_DIST() line 705 in
> > src/mat/impls/aij/mpi/superlu_dist/superlu_dist.c: Unknown column
> > permutation
> >
> >
> > Inside the code, I have these lines:
> >
> > ierr = KSPCreate(PETSC_COMM_WORLD,&ksp); CHKERRQ(ierr);
> > ierr = KSPSetOptionsPrefix(ksp,"lop_"); CHKERRQ(ierr);
> > ierr = KSPGetPC(ksp,&precond); CHKERRQ(ierr);
> > ierr = PCSetOptionsPrefix(precond,"lop_"); CHKERRQ(ierr);
> > ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
> > ierr = KSPSolve(ksp,b,*u); CHKERRQ(ierr);
> >
> > so all the options that I pass to PETSc (using "PetscOptionsSetValue()" or
> > the command line) must have the "lop_" prefix. These are my options:
> >
> > -lop_ksp_rtol 1.e-10 -lop_ksp_atol 1.e-10 -lop_ksp_monitor
> > -lop_ksp_converged_reason -lop_ksp_type gmres -lop_ksp_monitor_true_residual
> > -lop_pc_type lu -lop_pc_factor_mat_solver_package superlu_dist -lop_ksp_view
> >
> >
> > If I remove "-lop_ksp_view", the error that I reported disappears (and also
> > disappear the final summary
> > If I remove -lop_pc_type lu -lop_pc_factor_mat_solver_package superlu_dist"
> > and I use HYPRE ("-lop_pc_type hypre") and I keep the option "-lop_ksp_view"
> > as valid option, the problem does not appear.
> > If I remove "-lop_ksp_view" from command line but after "KSPSolve()" I add
> > "KSPView()", the problem persists.
> > If I suppress the prefix ("-ksp_rtol 1.e-10 -ksp_atol ... -ksp_view"), the
> > problem does not appear.
> >
> > PETSc-3.1-p3, PETSc-3.1-p4 and PETSc-3.1-p5 are not affected of this problem
> > (other machines that I use have these versions installed, the code runs
> > without problems). In the meanwhile I am recompiling PETSc-3.1-p6 on another
> > machine.
> >
> > Do you have any suggestion to skip/solve this problem? (if it is a real
> > problem, maybe I did something wrong...)
> >
> > Many thanks in advance,
> > Cheers
> >
> > --
> > Filippo SPIGA, MSc Computer Science
> > ~ homepage: http://tinyurl.com/fspiga ~
> >
> > «Nobody will drive us out of Cantor's paradise.»
> > -- David Hilbert
> >
> > *****
> > Disclaimer: "Please note this message and any attachment are CONFIDENTIAL an
> > may be privileged or otherwise protected from disclosure. The contents are
> > not to be disclosed to anyone other than the addressee. Unauthorized
> > recipients are requested to preserve this confidentiality and to advise the
> > sender immediately of any error in transmission."
> >
>
-------------- next part --------------
diff -r aabb5f301af6 src/mat/impls/aij/mpi/superlu_dist/superlu_dist.c
--- a/src/mat/impls/aij/mpi/superlu_dist/superlu_dist.c Sat Nov 20 16:12:27 2010 -0600
+++ b/src/mat/impls/aij/mpi/superlu_dist/superlu_dist.c Wed Dec 15 13:36:17 2010 -0600
@@ -701,8 +701,6 @@
ierr = PetscViewerASCIIPrintf(viewer," Column permutation MMD_ATA\n");CHKERRQ(ierr);
} else if (options.ColPerm == PARMETIS) {
ierr = PetscViewerASCIIPrintf(viewer," Column permutation PARMETIS\n");CHKERRQ(ierr);
- } else {
- SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown column permutation");
}
ierr = PetscViewerASCIIPrintf(viewer," Parallel symbolic factorization %s \n",PetscTruths[options.ParSymbFact != NO]);CHKERRQ(ierr);
More information about the petsc-dev
mailing list