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

Mark Adams mfadams at lbl.gov
Thu Jul 7 06:00:15 CDT 2022


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/
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20220707/cb30c914/attachment-0001.html>


More information about the petsc-users mailing list