[petsc-users] Command line option to set number iterations for Chebyshev eigenvalue estimation?

Dave May dave.mayhem23 at gmail.com
Thu Jan 23 09:02:00 CST 2014


It should be just "est".
However, I misspelled "cheby". Those options should be

-mg_levels_ksp_type chebyshev
-mg_levels_ksp_chebyshev_estimate_eigenvalues 0,0.2,0,1.1






On 23 January 2014 15:56, Matthew Knepley <knepley at gmail.com> wrote:

> On Thu, Jan 23, 2014 at 8:54 AM, Dave May <dave.mayhem23 at gmail.com> wrote:
>
>> If you don't have any prefix set on the outer KSP, you should be able to
>> simply do
>>
>> -mg_levels_ksp_type chebychev
>> -mg_levels_ksp_chebychev_estimate_eigenvalues 0,0.2,0,1.1
>> -mg_levels_est_ksp_type cg
>> -mg_levels_est_ksp_max_it X
>>
>
> Are you sure the est_ is not specest_?
>
>   Thanks,
>
>      Matt
>
>
>> Cheers,
>>   Dave
>>
>>
>>
>>
>>
>> On 23 January 2014 15:49, Garth N. Wells <gnw20 at cam.ac.uk> wrote:
>>
>>> On 2014-01-23 14:32, Matthew Knepley wrote:
>>>
>>>> On Thu, Jan 23, 2014 at 5:35 AM, Garth N. Wells <gnw20 at cam.ac.uk>
>>>> wrote:
>>>>
>>>>  Is there are command line option to set the number of iterations
>>>>> used for the eigenvalue estimation inside the Chebyshev
>>>>> preconditioner?  "-gamg_est_ksp_max_it" does the trick with GAMG,
>>>>> but I'd also like to change the number of iterations when using
>>>>> Chebyshev smoothing with ML.
>>>>>
>>>>
>>>> In GAMG, we use PETSc itself to do the estimate, and give it the
>>>> prefix gamg_est_. In ML, as far as I can see from the code, we use the
>>>> internal
>>>> estimator, and I don't think they give us this knob.
>>>>
>>>>
>>> I'm not so interested in the what happens inside ML, but the Chebyshev
>>> smoother that PETSc applies. With the default 10 GMRES iterations, the
>>> eigenvalue estimate is not good enough for my problem, which leads to the
>>> outer CG breaking down with an indefinite operator error.
>>>
>>> I can get it working by changing
>>>
>>>   -mg_levels_ksp_chebyshev_estimate_eigenvalues
>>>
>>> but I'd rather do a few more iterations to get a better eigenvalue
>>> estimate. It would also be nice to be able to switch from GMRES to CG. I
>>> can do this with GAMG via "-gamg_est_ksp_type"
>>>
>>> Garth
>>>
>>>
>>>     Matt
>>>>
>>>>
>>>>  Garth
>>>>>
>>>>
>>>> --
>>>> 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
>>>>
>>>
>>
>
>
> --
> 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/20140123/f3cc1ea5/attachment-0001.html>


More information about the petsc-users mailing list