[petsc-users] Block Preconditioning
Matthew Knepley
knepley at gmail.com
Fri May 29 06:18:54 CDT 2015
On Fri, May 29, 2015 at 5:09 AM, Elias Karabelas <karabelaselias at gmail.com>
wrote:
> Dear Matt,
>
> I tried out what you said. I used the following example I found with google
>
>
> https://raw.githubusercontent.com/petsc/petsc/5edff71f342f05166073de0ae93226f0997d7fe9/src/snes/examples/tutorials/ex70.c
>
Please do not start with this one. It exhibits the manual approach, which
is necessary sometimes, but
definitely less elegant than the command line.
> So now I have some questions concerning that
>
>
> a) When I use pc_fieldsplit_precondition_schur selfp I get slower
> convergence than with the user-provided Schur-Complement approximation in
> this example code. However from what I saw this should coincide in this
> case?
>
For the Stokes Equation selfp gives a Schur preconditioner that is not
spectrally equivalent:
B^T diag(A)^{-1} B \approx \Delta whereas S \approx I
so we generally use something that looks like that mass matrix. I think
that is what is in the code.
> b) I get a little confused with the settings of the inner and outer
> solvers in fieldsplit. So if I understood it correctly
> "-pc_fieldsplit_schur_fact_type upper" uses
>
> K -B^T
> 0 S
>
Yes.
> as preconditioner. Where K and S are inverted with the help of a KSP. So
> what I don't get is how to tweak the inner solvers correctly. For Example,
> In my own Implementation of preconditioned GMRES I also use the same
> structure as above as preconditioner and I take K^-1 to be approximated by
> the application of an AMG (in my case BoomerAMG) . Further I take S = B
> diag(K)^-1 B^T and again use an AMG to invert S. This results in my case in
> a good preconditioner.
>
> I tried now the same with this example and used
>
> -fieldsplit_0_ksp_type preonly -fieldsplit_0_pc_type hypre
> -fieldsplit_1_ksp_type preonly -fieldsplit_1_pc_type hypre
>
1) Always start with 'lu' so that you can see what is responsible for
deterioration in convergence. I always start
with the full Schur Complement factorization and accurate solves so
that it converges in 1 iterate. Then I
back off some of it.
2) That preconditioner is definitely not optimal
3) For questions about convergence, always send
-ksp_view -ksp_monitor_true_residual -ksp_converged_reason
and in this case turn on the monitors for the interior solves, e.g.
-fieldsplit_1_ksp_monitor_true_residual
Thanks,
Matt
> However I get really slow convergence of the outer GMRES method.
>
> Can someone give me some insight into this?
>
> Best regards
>
> Elias
>
>
> On 28.05.2015 18:02, Matthew Knepley wrote:
>
> On Thu, May 28, 2015 at 10:47 AM, Elias Karabelas <
> elias.karabelas at medunigraz.at> wrote:
>
>> Dear Members,
>>
>> I want to solve a Block System arising from the Discretization of a
>> stabilized Finite Element Formulation of the Stokes System.
>>
>> I have the following Block Structure
>>
>> A -B^T
>> B C
>>
>> The Preconditioner I intend to use is a block preconditioner of the Form
>>
>> A -B^T
>> S
>>
>> where S is an approximation of the Schur Complement. For applying the
>> inverse of the schur complement I want to use a Stabilized Least Squares
>> Commutator in the form
>>
>> S^-1 = (B diag(Q)^-1 B^T + C_1)^-1 (B diag(Q)^-1 A diag(Q)^-1 B^T + C_2)
>> (B diag(Q)^-1 B^T + C_1)^-1
>>
>> where Q is the mass matrix and C_1 and C_2 are some additional
>> stabilization matrices.
>>
>> I got from the Manual, that I can use the PCFieldSplit preconditioner for
>> generating the general Block preconditioner as indicated above. And I also
>> found that I can define some arbitrary PC with PCSHELL. My question is, if
>> it is possible to use PCSHELL to define the action of S^-1 as indicated
>> above.
>>
>
> 1) Use FieldSplit is the right PC to start with. Make sure you can do
> something simple like
>
> A -B^T
> C + B diag(A)^{-1} B^T
>
> with it before we do the more complicated thing.
>
> 2) You will want to implement a PC for the (1,1) block. You can use a
> PCSHELL, which is simpler to setup, but
> that means you will have to manually pull out the FieldSplit KSP and
> set it. If instead you define your own
> PC implementation, its more boilerplate code, but you could specify
> this PC from the command line without
> any FieldSplit specific code in your application.
>
> 3) Your PC will get two matrices, the MatSchurComplement, and the
> preconditioning matrix. If you set Q as the
> preconditioning matrix, or really if you set
>
> A 0
> 0 Q
>
> as the global preconditioning matrix, then the subsolve for (1,1) will get
> the Schur Complement and Q, and I think
> that is enough to build your Stabilized LSC PC.
>
> Let me know if this makes sense to you.
>
> Thanks,
>
> Matt
>
> Kind Regards
>> Elias Karabelas
>>
>> --
>> Elias Karabelas, Ph.D.
>>
>> Medical University of Graz
>> Institute of Biophysics
>> Harrachgasse 21/IV
>> 8010 Graz, Austria
>>
>> Phone: +43 316 380 7759
>> Email: elias.karabelas at medunigraz.at
>> Web : http://forschung.medunigraz.at/fodok/staff?name=EliasKarabelas
>>
>>
>
>
> --
> 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/20150529/4e8588b0/attachment.html>
More information about the petsc-users
mailing list