[petsc-users] PCASM in transposed format

Barry Smith bsmith at mcs.anl.gov
Mon Sep 22 22:47:43 CDT 2014


   Antoine,

    That is one nasty matrix!  You are actually getting essentially garbage during the solution process with and without the transpose. There is no reason to think that the additive Schwarz method, or any standard iterative method will work much at all on this matrix. 

    You should run with -ksp_monitor_true_residual to see that, though in the “preconditioned” “norm” it looks like reasonable convergence, it is actually not giving reasonable convergence in the two norm of the true residual regardless of transpose.

    I computed A’*A - A*A’ and the matrix is extremely far from being normal, non-normal matrices are hard to solve with iterative methods. 

     The diagonal entries of the matrix are microscopic compared to some non-diagonal entries, this is also a problem for iterative methods

    I also ran with -ksp_monitor_singular_value and -ksp_plot_eigenvalues -ksp_gmres_restart   the condition number of the matrix is not very large, maybe 10^5 but the distribution of the eigenvalues (er singular values) is nasty, lots and lots lined up on the imaginary axis.

    I ran with the PETSc built in sparse LU solver and got satisfactory results

./ex10 -f0 ~/Datafiles/Matrices/A_and_rhs.bin  -pc_type lu -ksp_monitor_true 
  0 KSP preconditioned resid norm 4.620896262080e-02 true resid norm 1.367267281011e-05 ||r(i)||/||b|| 1.000000000000e+00
  1 KSP preconditioned resid norm 6.143025615547e-13 true resid norm 3.207400714524e-16 ||r(i)||/||b|| 2.345847632770e-11

pretty quickly. 

    My recommendation is either to reformulate the problem to not get such a nasty matrix or use a direct solver. Iterative methods are going to hopeless for yo.

   Barry


On Sep 22, 2014, at 1:28 PM, Antoine De Blois <antoine.deblois at aero.bombardier.com> wrote:

> Dear all,
>  
> I am using the ASM preconditioner to solve a transposed system through MatSolveTranspose. Strangely, the results I obtain differ in each call. Is there a non-deterministic operations within ASM? If I use it in a non-transposed way, I get correct results… If I use GASM, then the results are always the same, both in transposed and non-transposed formats.
>  
> Below is a log of my calls (note that some of them actually diverged).
> I can give you my matrix and rhs if you want.
>  
> Regards,
> Antoine
>  
> $ mpirun -np 4 /gpfs/fs1/aero/SOFTWARE/TOOLS/PROGRAMMING/petsc/src/ksp/ksp/examples/tutorials/ex10 -f0 A_and_rhs.bin  -trans -pc_type asm              
> Number of iterations = 690
> Residual norm 0.000617544
>  
> $ mpirun -np 4 /gpfs/fs1/aero/SOFTWARE/TOOLS/PROGRAMMING/petsc/src/ksp/ksp/examples/tutorials/ex10 -f0 A_and_rhs.bin  -trans -pc_type asm
> Number of iterations = 475
> Residual norm 0.000727253
>  
> $ mpirun -np 4 /gpfs/fs1/aero/SOFTWARE/TOOLS/PROGRAMMING/petsc/src/ksp/ksp/examples/tutorials/ex10 -f0 A_and_rhs.bin  -trans -pc_type asm
> Number of iterations = 10000
> Residual norm 1.3866
>  
> $ mpirun -np 4 /gpfs/fs1/aero/SOFTWARE/TOOLS/PROGRAMMING/petsc/src/ksp/ksp/examples/tutorials/ex10 -f0 A_and_rhs.bin  -trans -pc_type asm
> Number of iterations = 568
> Residual norm 0.000684401
>  
> $ mpirun -np 4 /gpfs/fs1/aero/SOFTWARE/TOOLS/PROGRAMMING/petsc/src/ksp/ksp/examples/tutorials/ex10 -f0 A_and_rhs.bin  -trans -pc_type asm
> Number of iterations = 540
> Residual norm 0.000555548
>  
> $ mpirun -np 4 /gpfs/fs1/aero/SOFTWARE/TOOLS/PROGRAMMING/petsc/src/ksp/ksp/examples/tutorials/ex10 -f0 A_and_rhs.bin  -trans -pc_type asm
> Number of iterations = 10000
> Residual norm 1.30198
>  
> $ mpirun -np 4 /gpfs/fs1/aero/SOFTWARE/TOOLS/PROGRAMMING/petsc/src/ksp/ksp/examples/tutorials/ex10 -f0 A_and_rhs.bin  -trans -pc_type asm
> Number of iterations = 207
> Residual norm 0.000555849
>  
>  
> ---------------
>  
> $ mpirun -np 4 /gpfs/fs1/aero/SOFTWARE/TOOLS/PROGRAMMING/petsc/src/ksp/ksp/examples/tutorials/ex10 -f0 A_and_rhs.bin  -trans -pc_type gasm
> Number of iterations = 297
> Residual norm 0.000600143
>  
> $ mpirun -np 4 /gpfs/fs1/aero/SOFTWARE/TOOLS/PROGRAMMING/petsc/src/ksp/ksp/examples/tutorials/ex10 -f0 A_and_rhs.bin  -trans -pc_type gasm
> Number of iterations = 297
> Residual norm 0.000600143
>  
>  
>  
> Antoine DeBlois
> Specialiste ingenierie, MDO lead / Engineering Specialist, MDO lead
> Aéronautique / Aerospace
> 514-855-5001, x 50862
> antoine.deblois at aero.bombardier.com
> 
> 2351 Blvd Alfred-Nobel
> Montreal, Qc
> H4S 1A9
> 
> <image001.jpg>
> CONFIDENTIALITY NOTICE - This communication may contain privileged or confidential information.
> If you are not the intended recipient or received this communication by error, please notify the sender
> and delete the message without copying



More information about the petsc-users mailing list