[petsc-users] KSP solver increases the solution time

Alejandro Marcos Aragón alejandro.aragon at gmail.com
Wed Apr 27 02:17:27 CDT 2011


I understand that, but I'm trying to provide default behavior for the solver because the default one (no parameters) works very bad in my case.
However, I'm stuck because I can't set the same parameters that I obtain with command line arguments "-pc_type asm -sub_pc_type lu".

Can someone point me where is the error with the following code?

...
...
            PetscInitialize(&argc, &argv,NULL,NULL);
            PetscErrorCode ierr = MatCreate(PETSC_COMM_WORLD,&A_);CHKERR(ierr);
            
            // create linear solver context
            ierr = KSPCreate(PETSC_COMM_WORLD,&ksp_);CHKERR(ierr);
            
            // initial nonzero guess
            ierr = KSPSetInitialGuessNonzero(ksp_,PETSC_TRUE); CHKERR(ierr);
            
            // set runtime options
            ierr = KSPSetFromOptions(ksp_);CHKERR(ierr);

            // set the default preconditioner for this program to be ASM
            PC pc;
            ierr = KSPGetPC(ksp_,&pc); CHKERR(ierr);
            ierr = PCSetType(pc, PCASM); CHKERR(ierr);

            KSP        *subksp;       /* array of KSP contexts for local subblocks */
            PetscInt   nlocal,first;  /* number of local subblocks, first local subblock */
            PC         subpc;          /* PC context for subblock */
                                    
            /* 
             Call KSPSetUp() to set the block Jacobi data structures (including
             creation of an internal KSP context for each block).
             
             Note: KSPSetUp() MUST be called before PCASMGetSubKSP().
             */
            ierr = KSPSetUp(ksp_);CHKERR(ierr);
            
            /*
             Extract the array of KSP contexts for the local blocks
             */
            ierr = PCASMGetSubKSP(pc,&nlocal,&first,&subksp);CHKERR(ierr);
            
            /*
             Loop over the local blocks, setting various KSP options
             for each block.  
             */
            for (int i=0; i<nlocal; i++) {
                ierr = KSPGetPC(subksp[i],&subpc);CHKERR(ierr);
                ierr = PCSetType(subpc,PCLU);CHKERR(ierr);
            }

This is the error I get:

User explicitly sets subdomain solvers.
[0]PETSC ERROR: --------------------- Error Message ------------------------------------
[0]PETSC ERROR: Null argument, when expecting valid pointer!
[0]PETSC ERROR: Null Object: Parameter # 1!
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: Petsc Release Version 3.1.0, Patch 7, Mon Dec 20 14:26:37 CST 2010
[0]PETSC ERROR: See docs/changes/index.html for recent updates.
[0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
[0]PETSC ERROR: See docs/index.html for manual pages.
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: yafeq/a.out on a darwin10. named lsmspc26.epfl.ch by aaragon Wed Apr 27 09:14:25 2011
[0]PETSC ERROR: Libraries linked from /Users/aaragon/Local/lib
[0]PETSC ERROR: Configure run at Thu Apr  7 17:01:26 2011
[0]PETSC ERROR: Configure options --prefix=/Users/aaragon/Local --with-mpi-include=/Users/aaragon/Local/include --with-mpi-lib=/Users/aaragon/Local/lib/libmpich.a --with-superlu=1 --with-superlu-include=/Users/aaragon/Local/include/superlu --with-superlu-lib=/Users/aaragon/Local/lib/libsuperlu.a --with-superlu_dist=1 --with-superlu_dist-include=/Users/aaragon/Local/include/superlu_dist --with-superlu_dist-lib=/Users/aaragon/Local/lib/libsuperlu_dist.a --with-parmetis=1 --download-parmetis=ifneeded
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: MatGetVecs() line 7265 in src/mat/interface/matrix.c
[0]PETSC ERROR: KSPGetVecs() line 806 in src/ksp/ksp/interface/iterativ.c
[0]PETSC ERROR: KSPSetUp_GMRES() line 94 in src/ksp/ksp/impls/gmres/gmres.c
[0]PETSC ERROR: KSPSetUp() line 199 in src/ksp/ksp/interface/itfunc.c
[0]PETSC ERROR: User provided function() line 397 in "unknowndirectory/"/Users/aaragon/Local/include/cpputils/solver.hpp
[0]PETSC ERROR: --------------------- Error Message ------------------------------------
[0]PETSC ERROR: Operation done in wrong order!
[0]PETSC ERROR: Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here!
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: Petsc Release Version 3.1.0, Patch 7, Mon Dec 20 14:26:37 CST 2010
[0]PETSC ERROR: See docs/changes/index.html for recent updates.
[0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
[0]PETSC ERROR: See docs/index.html for manual pages.
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: yafeq/a.out on a darwin10. named lsmspc26.epfl.ch by aaragon Wed Apr 27 09:14:25 2011
[0]PETSC ERROR: Libraries linked from /Users/aaragon/Local/lib
[0]PETSC ERROR: Configure run at Thu Apr  7 17:01:26 2011
[0]PETSC ERROR: Configure options --prefix=/Users/aaragon/Local --with-mpi-include=/Users/aaragon/Local/include --with-mpi-lib=/Users/aaragon/Local/lib/libmpich.a --with-superlu=1 --with-superlu-include=/Users/aaragon/Local/include/superlu --with-superlu-lib=/Users/aaragon/Local/lib/libsuperlu.a --with-superlu_dist=1 --with-superlu_dist-include=/Users/aaragon/Local/include/superlu_dist --with-superlu_dist-lib=/Users/aaragon/Local/lib/libsuperlu_dist.a --with-parmetis=1 --download-parmetis=ifneeded
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: PCASMGetSubKSP_ASM() line 644 in src/ksp/pc/impls/asm/asm.c
[0]PETSC ERROR: PCASMGetSubKSP() line 926 in src/ksp/pc/impls/asm/asm.c
[0]PETSC ERROR: User provided function() line 402 in "unknowndirectory/"/Users/aaragon/Local/include/cpputils/solver.hpp
[0]PETSC ERROR: --------------------- Error Message ------------------------------------
[0]PETSC ERROR: Invalid argument!
[0]PETSC ERROR: Wrong type of object: Parameter # 1!
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: Petsc Release Version 3.1.0, Patch 7, Mon Dec 20 14:26:37 CST 2010

and the error continues...


On Apr 26, 2011, at 5:26 PM, Matthew Knepley wrote:

> 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/20110427/a143dc5f/attachment-0001.htm>


More information about the petsc-users mailing list