[petsc-users] [PCGAMG + AGG + GMRES] Non-Exact Dirichlet Boundary Conditions

Mark Adams mfadams at lbl.gov
Thu Jul 7 08:54:38 CDT 2022


GANG has failed here. Compare the real data.
Run with -ksp_converged_reason to start getting an idea what is going on.
GAMG/AGG failures are often from bad eigen estimates in the Chebyshev
smoother.
Is your problem symmetric?

On Thu, Jul 7, 2022 at 9:19 AM Karabelas, Elias (elias.karabelas at uni-graz.at)
<elias.karabelas at uni-graz.at> wrote:

> Am 07.07.22 um 15:05 schrieb Mark Adams:
>
>
>
> On Thu, Jul 7, 2022 at 7:03 AM Karabelas, Elias (
> elias.karabelas at uni-graz.at) <elias.karabelas at uni-graz.at> wrote:
>
>> OK got it, well the "classical" option of GAMG removes this issue, and
>> also HYPRE does that out-of-the box.
>>
>
> Humm, so you are using hypre in PETSc with a Krylov KSP or using hypre
> directly.
>
> GAMG uses a polynomial smoother by default and hypre does not, I believe.
> That could be the difference, but I would think the outer Krylov in PETSc
> would still pollute this.
>
> So the outer Krylov doesn't seem to pollute either for hypre or using gamg
> classical
>
> this is an example output from those two (vtu file) (gmres as krylov and
> either hypre or gamg classical as pctype)
>
> <DataArray type="Float32" format="ascii" NumberOfComponents="3"
> Name="displacement">
>           0 0 0
>           0.182029 0 0
>           0.223618 0 0
>           0.250511 0 0
>           0.209325 0 0
>           0 0.182029 0
>           0.1305 0.1305 0
>
>
> and this for gamg with smoothed aggregation (gmres as krylov and gamg agg)
>
> <DataArray type="Float32" format="ascii" NumberOfComponents="3"
> Name="displacement">
>           0 0 0
>           0.0905508 -5.2019e-05 6.4198e-05
>           0.173765 -4.34019e-05 5.86626e-05
>           0.256981 -4.14403e-05 2.78863e-05
>           0.233807 -3.94787e-05 -2.88998e-06
>           -4.78464e-05 0.0905198 6.21817e-05
>           0.0819423 0.0820083 4.89063e-05
>
> So the difference is noticeable.
>
>
>
>
>
> If you are using hypre directly, I could imagine that they take the
> residual that they compute to check convergence and do a quick update of
> the solution, with the diagonal that they have access to, with u = u + D^-1
> r
> (this is the sort of clever thing they would do)
>
>
>> I would slightly disagree with non-exact DirichletBCs, especially in
>> mechanics where I really assume a clamped node is clamped.
>>
>>
> I'm not saying exactness is not important (to you), I'm just saying many
> people seem to live with it.
>
>
>
>> Best
>> Elias
>>
>> Am 07.07.22 um 13:00 schrieb Mark Adams:
>>
>> I think PCCOMPOSITE is overkill here.
>>
>> First, I would only bother with this if this error is a problem. People
>> use your method all the time and accept error at the scale of the solver
>> tolerance.
>>
>> But if you want the exact solution, well, you could just clobber the
>> solution values after the solve if you have access to the Diri BC data
>> on hand.
>>
>> I was suggesting just creating a "Richardson/jacobi" solver, hardwire one
>> iteration, use an initial guess and solve it.
>> Not great, because this would have to grab the diagonal. If you had the
>> diagonal (eg, from the AMG smoother that does exist) then this would be
>> pretty cheap (a residual calculation mostly).
>>
>> Mark
>>
>>
>> On Thu, Jul 7, 2022 at 3:14 AM Karabelas, Elias (
>> elias.karabelas at uni-graz.at) <elias.karabelas at uni-graz.at> wrote:
>>
>>> Hi Mark,
>>> thanks for the answer, but I'm struggling to translate your suggestion
>>> into solver options.
>>> From scrolling through the user manual I think this points towards
>>> PCCOMPOSITE.
>>> However, the User Manual is not very precise with composite PCs, so how
>>> would I achieve this on top of my existing solver options?
>>> Best regards
>>> Elias
>>>
>>>
>>>
>>> Am 06.07.22 um 16:41 schrieb Mark Adams:
>>>
>>> And one iteration of undamped Jacobi after the solve should fix this.
>>>
>>> On Wed, Jul 6, 2022 at 8:42 AM Karabelas, Elias (
>>> elias.karabelas at uni-graz.at) <elias.karabelas at uni-graz.at> wrote:
>>>
>>>> Dear Matt,
>>>>
>>>> thanks for the fast response. That makes perfect sense to me.
>>>>
>>>> Best regards
>>>> Elias
>>>>
>>>> Am 06.07.22 um 14:35 schrieb Matthew Knepley:
>>>>
>>>> On Wed, Jul 6, 2022 at 7:46 AM Karabelas, Elias (
>>>> elias.karabelas at uni-graz.at) <elias.karabelas at uni-graz.at> wrote:
>>>>
>>>>> Dear all,
>>>>>
>>>>> I don't know if this is a bug, but I observed that when using GMRES
>>>>> with AGG-PCGAMG as preconditioner Dirichlet boundary conditions don't seem
>>>>> to be exactly fulfilled.
>>>>>
>>>>> My Matrix has zero rows and cols with 1 on the diagonal where I have
>>>>> dirichlet-bcs in my FE-mesh and I would expect that the eqs in this rows
>>>>> can be exactly fulfilled (as u_i = g_i) there.
>>>>>
>>>> I would not expect aggregation to be exact here, but only within the
>>>> iteration tolerance. If instead you eliminate those variables, you can
>>>> maintain algebraic exactness.
>>>> This is what we do in examples, like SNES ex56.
>>>>
>>>>   Thanks,
>>>>
>>>>      Matt
>>>>
>>>>> However, when I solve A*x = b with the above solver I only get u_i =
>>>>> g_i + error in that part of the vector. Switching from pc_gamg_type agg to
>>>>> pc_gamg_type classical cures this problem, but the classical is not
>>>>> advertised in the user manual.
>>>>>
>>>>> These are the options I'm currently using:
>>>>>
>>>>> -ksp_type gmres
>>>>> -ksp_pc_side right
>>>>> -pc_type gamg
>>>>> -pc_gamg_type agg [or classical]
>>>>> -pc_gamg_sym_graph 1
>>>>> -pc_gamg_square_graph 1
>>>>> -pc_gamg_agg_nsmooths 1
>>>>> -pc_gamg_threshold 0.01
>>>>> -pc_mg_cycles v
>>>>>
>>>>> Iteration counts are basically the same.
>>>>>
>>>>> Best regards
>>>>>
>>>>> Elias
>>>>>
>>>>> --
>>>>> Dr. Elias Karabelas
>>>>> Research Associate
>>>>> University of Graz
>>>>> Institute of Mathematics and Scientific Computing
>>>>> Heinrichstraße 36
>>>>> A-8010 Graz
>>>>> Austria
>>>>>
>>>>> Phone: +43 316 380 8546
>>>>> Email: elias.karabelas at uni-graz.at
>>>>> Web:  https://ccl.medunigraz.at/
>>>>>
>>>>>
>>>>
>>>> --
>>>> 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.cse.buffalo.edu/~knepley/>
>>>>
>>>>
>>>> --
>>>> Dr. Elias Karabelas
>>>> Research Associate
>>>> University of Graz
>>>> Institute of Mathematics and Scientific Computing
>>>> Heinrichstraße 36
>>>> A-8010 Graz
>>>> Austria
>>>>
>>>> Phone: +43 316 380 8546
>>>> Email: elias.karabelas at uni-graz.at
>>>> Web:  https://ccl.medunigraz.at/
>>>>
>>>>
>>> --
>>> Dr. Elias Karabelas
>>> Research Associate
>>> University of Graz
>>> Institute of Mathematics and Scientific Computing
>>> Heinrichstraße 36
>>> A-8010 Graz
>>> Austria
>>>
>>> Phone: +43 316 380 8546
>>> Email: elias.karabelas at uni-graz.at
>>> Web:  https://ccl.medunigraz.at/
>>>
>>>
>> --
>> Dr. Elias Karabelas
>> Research Associate
>> University of Graz
>> Institute of Mathematics and Scientific Computing
>> Heinrichstraße 36
>> A-8010 Graz
>> Austria
>>
>> Phone: +43 316 380 8546
>> Email: elias.karabelas at uni-graz.at
>> Web:  https://ccl.medunigraz.at/
>>
>>
> --
> Dr. Elias Karabelas
> Research Associate
> University of Graz
> Institute of Mathematics and Scientific Computing
> Heinrichstraße 36
> A-8010 Graz
> Austria
>
> Phone: +43 316 380 8546
> Email: elias.karabelas at uni-graz.at
> Web:  https://ccl.medunigraz.at/
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20220707/632e5ff5/attachment.html>


More information about the petsc-users mailing list