[petsc-users] Solving Linear Systems with Scalar Real and Complex

Matthew Knepley knepley at gmail.com
Fri Aug 3 09:11:18 CDT 2018


On Fri, Aug 3, 2018 at 9:42 AM Pierpaolo Minelli <pierpaolo.minelli at cnr.it>
wrote:

> I tried to use these options with GAMG:
>
> -mg_levels_ksp_type richardson
>
> -mg_levels_pc_type sor
>
>
> but also in this case i don’t obtain a convergence.
>

Can you send the convergence and solver view so we can check?

Also, if you send the matrix we can try running it.


> I used also direct solvers (MUMPS ad SuperLU) and they works fine, but
> slower (solve time) then iterative method with ML and hypre.
>
> If i check with -ksp_monitor_true_residual, is it advisable to continue
> using the iterative method?
>

Yes.

   Matt


> Pierpaolo
>
>
> Il giorno 03 ago 2018, alle ore 15:16, Mark Adams <mfadams at lbl.gov> ha
> scritto:
>
> So this is a complex valued indefinite Helmholtz operator (very hard to
> solve scalably) with axisymmetric coordinates. ML, hypre and GAMG all
> performed about the same, with a big jump in residual initially and
> essentially not solving it. You scaled it and this fixed ML and hypre but
> not GAMG.
>
> From this output I can see that the eigenvalue estimates are strange. Your
> equations look fine so I have to assume that the complex values are the
> problem. If this is symmetric the CG is a much better solver and eigen
> estimator. But this is not a big deal, especially since you have two
> options that work. I would suggest not using cheby smoother, it uses these
> bad eigen estimates, and it is basically not smoothing on some levels. You
> can use this instead:
>
> -mg_levels_ksp_type richardson
> -mg_levels_pc_type sor
>
> Note, if you have a large shift these equations are very hard to solve
> iteratively and you should just use a direct solver. Direct solvers in 2D
> are not bad,
>
> Mark
>
>
> On Fri, Aug 3, 2018 at 3:02 AM Pierpaolo Minelli <pierpaolo.minelli at cnr.it>
> wrote:
>
>> In this simulation I'm solving two equations in a two-dimensional domain
>> (z,r) at each time step. The first is an equation derived from the Maxwell
>> equation. Taking advantage of the fact that the domain is axialsymmetric
>> and applying a temporal harmonic approximation, the equation I am solving
>> is the following:
>>
>>
>> The second equation is a Poisson’s  equation in cylindrical coordinates
>> (z,r) in the field of real numbers.
>>
>> This is the output obtained using these options (Note that at this moment
>> of development I am only using a processor):
>>
>> *-pc_type gamg -pc_gamg_agg_nsmooths 1 -pc_gamg_reuse_interpolation true
>> -pc_gamg_square_graph 1 -pc_gamg_threshold 0. -ksp_rtol 1.e-7
>> -ksp_max_it 30 -ksp_monitor_true_residual -info | grep GAMG*
>>
>> [0] PCSetUp_GAMG(): level 0) N=321201, n data rows=1, n data cols=1,
>> nnz/row (ave)=5, np=1
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 4.97011 nnz ave. (N=321201)
>> [0] PCGAMGCoarsen_AGG(): Square Graph on level 1 of 1 to square
>> [0] PCGAMGProlongator_AGG(): New grid 45754 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=1.972526e+00
>> min=3.461411e-03 PC=jacobi
>> [0] PCSetUp_GAMG(): 1) N=45754, n data cols=1, nnz/row (ave)=10, 1 active
>> pes
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 10.7695 nnz ave. (N=45754)
>> [0] PCGAMGProlongator_AGG(): New grid 7893 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=5.686837e+00
>> min=5.062501e-01 PC=jacobi
>> [0] PCSetUp_GAMG(): 2) N=7893, n data cols=1, nnz/row (ave)=23, 1 active
>> pes
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 23.2179 nnz ave. (N=7893)
>> [0] PCGAMGProlongator_AGG(): New grid 752 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=1.504451e+01
>> min=2.124898e-02 PC=jacobi
>> [0] PCSetUp_GAMG(): 3) N=752, n data cols=1, nnz/row (ave)=30, 1 active
>> pes
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 30.7367 nnz ave. (N=752)
>> [0] PCGAMGProlongator_AGG(): New grid 56 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=7.781296e+00
>> min=2.212257e-02 PC=jacobi
>> [0] PCSetUp_GAMG(): 4) N=56, n data cols=1, nnz/row (ave)=22, 1 active pes
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 22.9643 nnz ave. (N=56)
>> [0] PCGAMGProlongator_AGG(): New grid 6 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=1.525086e+00
>> min=1.375043e-01 PC=jacobi
>> [0] PCSetUp_GAMG(): 5) N=6, n data cols=1, nnz/row (ave)=6, 1 active pes
>> [0] PCSetUp_GAMG(): 6 levels, grid complexity = 1.43876
>> [0] PCSetUp_GAMG(): level 0) N=321201, n data rows=1, n data cols=1,
>> nnz/row (ave)=5, np=1
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 4.97011 nnz ave. (N=321201)
>> [0] PCGAMGCoarsen_AGG(): Square Graph on level 1 of 1 to square
>> [0] PCGAMGProlongator_AGG(): New grid 45754 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=1.972526e+00
>> min=3.461411e-03 PC=jacobi
>> [0] PCSetUp_GAMG(): 1) N=45754, n data cols=1, nnz/row (ave)=10, 1 active
>> pes
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 10.7695 nnz ave. (N=45754)
>> [0] PCGAMGProlongator_AGG(): New grid 7893 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=5.686837e+00
>> min=5.062501e-01 PC=jacobi
>> [0] PCSetUp_GAMG(): 2) N=7893, n data cols=1, nnz/row (ave)=23, 1 active
>> pes
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 23.2179 nnz ave. (N=7893)
>> [0] PCGAMGProlongator_AGG(): New grid 752 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=1.504451e+01
>> min=2.124898e-02 PC=jacobi
>> [0] PCSetUp_GAMG(): 3) N=752, n data cols=1, nnz/row (ave)=30, 1 active
>> pes
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 30.7367 nnz ave. (N=752)
>> [0] PCGAMGProlongator_AGG(): New grid 56 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=7.781296e+00
>> min=2.212257e-02 PC=jacobi
>> [0] PCSetUp_GAMG(): 4) N=56, n data cols=1, nnz/row (ave)=22, 1 active pes
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 22.9643 nnz ave. (N=56)
>> [0] PCGAMGProlongator_AGG(): New grid 6 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=1.525086e+00
>> min=1.375043e-01 PC=jacobi
>> [0] PCSetUp_GAMG(): 5) N=6, n data cols=1, nnz/row (ave)=6, 1 active pes
>> [0] PCSetUp_GAMG(): 6 levels, grid complexity = 1.43876
>> [0] PCSetUp_GAMG(): level 0) N=271201, n data rows=1, n data cols=1,
>> nnz/row (ave)=5, np=1
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 4.97455 nnz ave. (N=271201)
>> [0] PCGAMGCoarsen_AGG(): Square Graph on level 1 of 1 to square
>> [0] PCGAMGProlongator_AGG(): New grid 38501 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=1.933798e+00
>> min=4.684075e-02 PC=jacobi
>> [0] PCSetUp_GAMG(): 1) N=38501, n data cols=1, nnz/row (ave)=10, 1 active
>> pes
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 10.7732 nnz ave. (N=38501)
>> [0] PCGAMGProlongator_AGG(): New grid 6664 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=1.623029e+00
>> min=1.250957e-02 PC=jacobi
>> [0] PCSetUp_GAMG(): 2) N=6664, n data cols=1, nnz/row (ave)=23, 1 active
>> pes
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 23.2098 nnz ave. (N=6664)
>> [0] PCGAMGProlongator_AGG(): New grid 620 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=1.763329e+00
>> min=1.611776e-02 PC=jacobi
>> [0] PCSetUp_GAMG(): 3) N=620, n data cols=1, nnz/row (ave)=29, 1 active
>> pes
>> [0] PCGAMGFilterGraph():   100.% nnz after filtering, with threshold 0.,
>> 29.6129 nnz ave. (N=620)
>> [0] PCGAMGProlongator_AGG(): New grid 46 nodes
>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=1.497611e+00
>> min=2.630403e-02 PC=jacobi
>> [0] PCSetUp_GAMG(): 4) N=46, n data cols=1, nnz/row (ave)=20, 1 active pes
>> [0] PCSetUp_GAMG(): 5 levels, grid complexity = 1.43639
>>
>>
>>
>>
>> Il giorno 02 ago 2018, alle ore 17:39, Mark Adams <mfadams at lbl.gov> ha
>> scritto:
>>
>> It looks like ML and hypre are working well now. If you want to debug
>> GAMG you can run with -info, which is very noisy, and grep on GAMG or send
>> the whole output.
>>
>> BTW, what equations are you solving?
>>
>> On Thu, Aug 2, 2018 at 5:12 AM Pierpaolo Minelli <
>> pierpaolo.minelli at cnr.it> wrote:
>>
>>> Thank you very much for the correction.
>>> By rebalancing the matrix coefficients, making them dimensionless, in my
>>> complex problem i managed to obtain a better result with both ML and HYPRE.
>>> GAMG instead seems to be unable to converge and in fact I get unexpected
>>> results. I report the outputs of the three simulations again (i tried to
>>> use also -pc_gamg_square_graph 20 without any improvement).
>>>
>>> Pierpaolo
>>>
>>> *-pc_type hypre -ksp_rtol 1.e-7 -ksp_monitor_true_residual*
>>>
>>>   0 KSP preconditioned resid norm 1.984853668904e-02 true resid norm
>>> 2.979865703850e-03 ||r(i)||/||b|| 1.000000000000e+00
>>>   1 KSP preconditioned resid norm 1.924446712661e-04 true resid norm
>>> 1.204260204811e-04 ||r(i)||/||b|| 4.041323752460e-02
>>>   2 KSP preconditioned resid norm 5.161509100765e-06 true resid norm
>>> 2.810809726926e-06 ||r(i)||/||b|| 9.432672496933e-04
>>>   3 KSP preconditioned resid norm 9.297326931238e-08 true resid norm
>>> 4.474617977876e-08 ||r(i)||/||b|| 1.501617328625e-05
>>>   4 KSP preconditioned resid norm 1.910271882670e-09 true resid norm
>>> 9.637470658283e-10 ||r(i)||/||b|| 3.234196308186e-07
>>>   0 KSP preconditioned resid norm 2.157687745805e+04 true resid norm
>>> 3.182001523188e+03 ||r(i)||/||b|| 1.000000000000e+00
>>>   1 KSP preconditioned resid norm 1.949268476386e+02 true resid norm
>>> 1.243419788627e+02 ||r(i)||/||b|| 3.907665598415e-02
>>>   2 KSP preconditioned resid norm 5.078054475792e+00 true resid norm
>>> 2.745355604400e+00 ||r(i)||/||b|| 8.627763325675e-04
>>>   3 KSP preconditioned resid norm 8.663802743529e-02 true resid norm
>>> 4.254290979292e-02 ||r(i)||/||b|| 1.336985839979e-05
>>>   4 KSP preconditioned resid norm 1.795605563039e-03 true resid norm
>>> 9.040507428245e-04 ||r(i)||/||b|| 2.841138623714e-07
>>>   0 KSP preconditioned resid norm 6.728304961395e+02 true resid norm
>>> 1.879478105170e+02 ||r(i)||/||b|| 1.000000000000e+00
>>>   1 KSP preconditioned resid norm 2.190497539532e+01 true resid norm
>>> 4.630095820203e+02 ||r(i)||/||b|| 2.463500802413e+00
>>>   2 KSP preconditioned resid norm 8.425561564252e-01 true resid norm
>>> 7.012565302251e+01 ||r(i)||/||b|| 3.731123700223e-01
>>>   3 KSP preconditioned resid norm 3.029848345705e-02 true resid norm
>>> 4.379018464663e+00 ||r(i)||/||b|| 2.329911932795e-02
>>>   4 KSP preconditioned resid norm 7.374025528575e-04 true resid norm
>>> 1.337183702137e-01 ||r(i)||/||b|| 7.114654320570e-04
>>>   5 KSP preconditioned resid norm 3.009400175162e-05 true resid norm
>>> 7.731135032616e-03 ||r(i)||/||b|| 4.113447776459e-05
>>>
>>> *-pc_type ml -ksp_rtol 1.e-7 -ksp_monitor_true_residual*
>>>
>>>   0 KSP preconditioned resid norm 1.825767020538e-02 true resid norm
>>> 2.979865703850e-03 ||r(i)||/||b|| 1.000000000000e+00
>>>   1 KSP preconditioned resid norm 6.495628259383e-04 true resid norm
>>> 3.739440526742e-04 ||r(i)||/||b|| 1.254902367550e-01
>>>   2 KSP preconditioned resid norm 4.971875712015e-05 true resid norm
>>> 2.118856024328e-05 ||r(i)||/||b|| 7.110575559127e-03
>>>   3 KSP preconditioned resid norm 3.726806462912e-06 true resid norm
>>> 1.370355844514e-06 ||r(i)||/||b|| 4.598716790303e-04
>>>   4 KSP preconditioned resid norm 2.496898447120e-07 true resid norm
>>> 9.494701893753
>>>
>> <clip_image002.png><clip_image002.png>
>
>
>

-- 
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

https://www.cse.buffalo.edu/~knepley/ <http://www.caam.rice.edu/~mk51/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20180803/cecc5e8b/attachment.html>


More information about the petsc-users mailing list