ddemidov / amgcl

C++ library for solving large sparse linear systems with algebraic multigrid method
http://amgcl.readthedocs.org/
MIT License
728 stars 111 forks source link

Preconditioner failure #106

Closed JabirSS closed 3 years ago

JabirSS commented 5 years ago

Dear @ddemidov, Thank you for creating amgcl! I've been using it to solve the electric potential Poisson equation in a particle simulation and the performance has been quite impressive. However, recently I've been running into issues with the preconditioner as the simulation progresses as it fails to create any levels and gmres takes over 1500 iterations to yield a solution with acceptable tolerance (1.0e-4). The matrix has 68921 rows/cols, around 6 million non-zero elements and does not appear to be singular. I have uploaded the matrix in mtx format to the following link and would greatly appreciate any suggestions you may have regarding optimal parameters and combination of options to help solve it faster.

link: https://drive.google.com/file/d/14LrTqWf8Vkpo0iphj9c3Bgz5Xfe1Sq4g/view?usp=sharing

Thanks!

ddemidov commented 5 years ago

it fails to create any levels

I can not reproduce this, I get normal looking hierarchy with your matrix (but I can see that it does not converge with the default settings):

$ ./solver -A problems/issue-106/Matrix_at_preconditioner_failure.mtx 
Solver
======
Type:             BiCGStab
Unknowns:         68921
Memory footprint: 3.68 M

Preconditioner
==============
Number of levels:    3
Operator complexity: 1.21
Grid complexity:     1.26
Memory footprint:    122.16 M

level     unknowns       nonzeros      memory
---------------------------------------------
    0        68921        5966673     97.77 M (82.63%)
    1        16676        1196948     21.01 M (16.58%)
    2          910          57596      3.38 M ( 0.80%)

Iterations: 100
Error:      0.01114

[Profile:      24.111 s] (100.00%)
[  reading:    16.998 s] ( 70.50%)
[  setup:       0.265 s] (  1.10%)
[  solve:       6.848 s] ( 28.40%)

Switching to ILU0 for relaxation seems to work though:

$ ./solver -A problems/issue-106/Matrix_at_preconditioner_failure.mtx -p precond.relax.type=ilu0
Solver
======
Type:             BiCGStab
Unknowns:         68921
Memory footprint: 3.68 M

Preconditioner
==============
Number of levels:    3
Operator complexity: 1.21
Grid complexity:     1.26
Memory footprint:    233.10 M

level     unknowns       nonzeros      memory
---------------------------------------------
    0        68921        5966673    190.08 M (82.63%)
    1        16676        1196948     39.63 M (16.58%)
    2          910          57596      3.38 M ( 0.80%)

Iterations: 41
Error:      6.08542e-09

[Profile:      21.187 s] (100.00%)
[  reading:    15.323 s] ( 72.32%)
[  setup:       1.253 s] (  5.92%)
[  solve:       4.610 s] ( 21.76%)
ddemidov commented 5 years ago

Adjusting precond.coarsening.aggr.eps_strong also helps, even with the default spai0 relaxation:

$ ./solver -B -A problems/issue-106/Matrix_at_preconditioner_failure.bin -p precond.coarsening.aggr.eps_strong=0.05
Solver
======
Type:             BiCGStab
Unknowns:         68921
Memory footprint: 3.68 M

Preconditioner
==============
Number of levels:    3
Operator complexity: 1.14
Grid complexity:     1.14
Memory footprint:    114.81 M

level     unknowns       nonzeros      memory
---------------------------------------------
    0        68921        5966673     99.75 M (87.55%)
    1         9512         819936     14.23 M (12.03%)
    2          388          28728    846.79 K ( 0.42%)

Iterations: 18
Error:      1.81686e-09

[Profile:       1.464 s] (100.00%)
[  reading:     0.088 s] (  6.00%)
[  setup:       0.223 s] ( 15.20%)
[  solve:       1.153 s] ( 78.74%)

But it seems bicgstab is the only solver that is able to converge with this system. Can you also provide your right-hand-side (I am using a vector of 1 in these examples, this could affect the convergence)?

ddemidov commented 5 years ago

Using ruge_stuben for coarsening I can reproduce the 'single level' issue:

./solver -B -A problems/issue-106/Matrix_at_preconditioner_failure.bin -p precond.coarsening.type=ruge_stuben
Solver
======
Type:             BiCGStab
Unknowns:         68921
Memory footprint: 3.68 M

Preconditioner
==============
Number of levels:    1
Operator complexity: 1.00
Grid complexity:     1.00
Memory footprint:    93.67 M

level     unknowns       nonzeros      memory
---------------------------------------------
    0        68921        5966673     93.67 M (100.00%)

Iterations: 100
Error:      4.12735e+12

[Profile:       4.566 s] (100.00%)
[  reading:     0.080 s] (  1.75%)
[  setup:       0.067 s] (  1.47%)
[  solve:       4.419 s] ( 96.77%)

This is solved by flipping the sign of the matrix values (ruge_stuben assumes matrix diagonal to be positive):

./solver -A problems/issue-106/flip_sign.mtx -p precond.coarsening.type=ruge_stuben
Solver
======
Type:             BiCGStab
Unknowns:         68921
Memory footprint: 3.68 M

Preconditioner
==============
Number of levels:    3
Operator complexity: 1.21
Grid complexity:     1.21
Memory footprint:    132.05 M

level     unknowns       nonzeros      memory
---------------------------------------------
    0        68921        5966673     99.98 M (82.98%)
    1        12223        1039989     17.46 M (14.46%)
    2         2083         183639     14.61 M ( 2.55%)

Iterations: 12
Error:      1.64684e-09

[Profile:      10.067 s] (100.00%)
[  reading:     8.318 s] ( 82.63%)
[  setup:       0.866 s] (  8.60%)
[  solve:       0.882 s] (  8.76%)
JabirSS commented 5 years ago

Thank you for your swift response! You are abolutely right. I have tried several settings and found ruge_stuben to be the fastest and most consistent until this issue appeared. I'll give it a shot and let you know if any further issues appeared as the simulation progressed

On Wed, Feb 27, 2019, 18:40 Denis Demidov notifications@github.com wrote:

Using ruge_stuben for coarsening I can reproduce the 'single level' issue:

./solver -B -A problems/issue-106/Matrix_at_preconditioner_failure.bin -p precond.coarsening.type=ruge_stuben Solver

Type: BiCGStab Unknowns: 68921 Memory footprint: 3.68 M

Preconditioner

Number of levels: 1 Operator complexity: 1.00 Grid complexity: 1.00 Memory footprint: 93.67 M

level unknowns nonzeros memory

0        68921        5966673     93.67 M (100.00%)

Iterations: 100 Error: 4.12735e+12

[Profile: 4.566 s] (100.00%) [ reading: 0.080 s] ( 1.75%) [ setup: 0.067 s] ( 1.47%) [ solve: 4.419 s] ( 96.77%)

This is solved by flipping the sign of the matrix values (ruge_stuben assumes matrix diagonal to be positive):

./solver -A problems/issue-106/flip_sign.mtx -p precond.coarsening.type=ruge_stuben Solver

Type: BiCGStab Unknowns: 68921 Memory footprint: 3.68 M

Preconditioner

Number of levels: 3 Operator complexity: 1.21 Grid complexity: 1.21 Memory footprint: 132.05 M

level unknowns nonzeros memory

0        68921        5966673     99.98 M (82.98%)
1        12223        1039989     17.46 M (14.46%)
2         2083         183639     14.61 M ( 2.55%)

Iterations: 12 Error: 1.64684e-09

[Profile: 10.067 s] (100.00%) [ reading: 8.318 s] ( 82.63%) [ setup: 0.866 s] ( 8.60%) [ solve: 0.882 s] ( 8.76%)

— You are receiving this because you authored the thread. Reply to this email directly, view it on GitHub https://github.com/ddemidov/amgcl/issues/106#issuecomment-467793223, or mute the thread https://github.com/notifications/unsubscribe-auth/Aku2sIajBWFoat7K7xOMCsQoLa0wtjpmks5vRlKAgaJpZM4bUCF7 .