[petsc-users] KSP solver increases the solution time

Matthew Knepley knepley at gmail.com
Tue Apr 26 10:26:54 CDT 2011


2011/4/26 Alejandro Marcos Aragón <alejandro.aragon at gmail.com>

> Thank you Barry, your suggestion really helped speed up the program. The
> maximum number of iterations is 48. I still don't know what the asm
> pre-conditioner is but I guess I just need to read the manual. I'm trying to
> add code to do what you suggested automatically, and I found that I can add:
>
>             PC pc;
>             ierr = KSPGetPC(ksp_,&pc); CHKERR(ierr);
>             ierr = PCSetType(pc, PCASM); CHKERR(ierr);
>
> However, I cannot find the function to replace the -sub_pc_type. Can you
> point me where to look? My system may not be symmetric so I can't use the cg
> option.
>

1) Hard coding it in your program does not make sense. You gain nothing, and
lose a lot of flexibility.

2) You can do this using


http://www.mcs.anl.gov/petsc/petsc-as/snapshots/petsc-current/docs/manualpages/PC/PCASMGetSubKSP.html

    Matt


> Thanks again for your response.
>
>>
>
> On Apr 26, 2011, at 2:23 PM, Barry Smith wrote:
>
>
> System solved in 3664 iterations... 42.683 s
>
>
>   The preconditioner is simply not up to the task it has been assigned.
>  This number of iterations is problematic.
>
>    Have you tried -pc_type asm -sub_pc_type lu     If that works well you
> can try -pc_type asm -sub_pc_type ilu and see if that still works.
>
>   If the matrix is indeed symmetric positive definite you will want to use
> -ksp_type cg
>
>
>
>    Barry
>
>
> On Apr 26, 2011, at 1:37 AM, Alejandro Marcos Aragón wrote:
>
> Hi all,
>
>
> I'm using the standard configuration of the KSP solver, but the time it
> takes to solve a large system of equations is increasing (because of the
> increasing number of iterations?). These are my timing lines and the log
> from the KSP solver in two consecutive solves:
>
>
> [1]   Solving system... 0.154203 s
>
>
> KSP Object:
>
> type: gmres
>
>   GMRES: restart=30, using Classical (unmodified) Gram-Schmidt
> Orthogonalization with no iterative refinement
>
>   GMRES: happy breakdown tolerance 1e-30
>
> maximum iterations=10000
>
> tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
>
> left preconditioning
>
> using nonzero initial guess
>
> using PRECONDITIONED norm type for convergence test
>
> PC Object:
>
> type: bjacobi
>
>   block Jacobi: number of blocks = 3
>
>   Local solve is same for all blocks, in the following KSP and PC objects:
>
> KSP Object:(sub_)
>
>   type: preonly
>
>   maximum iterations=10000, initial guess is zero
>
>   tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
>
>   left preconditioning
>
>   using PRECONDITIONED norm type for convergence test
>
> PC Object:(sub_)
>
>   type: ilu
>
>     ILU: out-of-place factorization
>
>     0 levels of fill
>
>     tolerance for zero pivot 1e-12
>
>     using diagonal shift to prevent zero pivot
>
>     matrix ordering: natural
>
>     factor fill ratio given 1, needed 1
>
>       Factored matrix follows:
>
>         Matrix Object:
>
>           type=seqaij, rows=2020, cols=2020
>
>           package used to perform factorization: petsc
>
>           total: nonzeros=119396, allocated nonzeros=163620
>
>             using I-node routines: found 676 nodes, limit used is 5
>
>   linear system matrix = precond matrix:
>
>   Matrix Object:
>
>     type=seqaij, rows=2020, cols=2020
>
>     total: nonzeros=119396, allocated nonzeros=163620
>
>       using I-node routines: found 676 nodes, limit used is 5
>
> linear system matrix = precond matrix:
>
> Matrix Object:
>
>   type=mpiaij, rows=6058, cols=6058
>
>   total: nonzeros=365026, allocated nonzeros=509941
>
>     using I-node (on process 0) routines: found 676 nodes, limit used is 5
>
>
>
> [1]   System solved in 51 iterations... 0.543215 s
>
> ...
>
> ...
>
> ...
>
>
> [1]   Solving system... 0.302414 s
>
>
> KSP Object:
>
> type: gmres
>
>   GMRES: restart=30, using Classical (unmodified) Gram-Schmidt
> Orthogonalization with no iterative refinement
>
>   GMRES: happy breakdown tolerance 1e-30
>
> maximum iterations=10000
>
> tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
>
> left preconditioning
>
> using nonzero initial guess
>
> using PRECONDITIONED norm type for convergence test
>
> PC Object:
>
> type: bjacobi
>
>   block Jacobi: number of blocks = 3
>
>   Local solve is same for all blocks, in the following KSP and PC objects:
>
> KSP Object:(sub_)
>
>   type: preonly
>
>   maximum iterations=10000, initial guess is zero
>
>   tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
>
>   left preconditioning
>
>   using PRECONDITIONED norm type for convergence test
>
> PC Object:(sub_)
>
>   type: ilu
>
>     ILU: out-of-place factorization
>
>     0 levels of fill
>
>     tolerance for zero pivot 1e-12
>
>     using diagonal shift to prevent zero pivot
>
>     matrix ordering: natural
>
>     factor fill ratio given 1, needed 1
>
>       Factored matrix follows:
>
>         Matrix Object:
>
>           type=seqaij, rows=2020, cols=2020
>
>           package used to perform factorization: petsc
>
>           total: nonzeros=119396, allocated nonzeros=163620
>
>             using I-node routines: found 676 nodes, limit used is 5
>
>   linear system matrix = precond matrix:
>
>   Matrix Object:
>
>     type=seqaij, rows=2020, cols=2020
>
>     total: nonzeros=119396, allocated nonzeros=163620
>
>       using I-node routines: found 676 nodes, limit used is 5
>
> linear system matrix = precond matrix:
>
> Matrix Object:
>
>   type=mpiaij, rows=6058, cols=6058
>
>   total: nonzeros=365026, allocated nonzeros=509941
>
>     using I-node (on process 0) routines: found 676 nodes, limit used is 5
>
>
> [1]   System solved in 3664 iterations... 42.683 s
>
>
>
>
> As you can see, the second iteration takes more than 40 seconds to solve.
> Could some explain why this is happening and why he number of iterations is
> increasing dramatically between solves? Thank you all,
>
>
> Alejandro M. Aragón, Ph.D.
>
>
>
>


-- 
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/20110426/b48f8601/attachment-0001.htm>


More information about the petsc-users mailing list