Hi Jie,

thanks a lot. The link to your repository was very helpful. I could significantly reduce the number of iterations using your implementation of computing the pressure laplacian through:

K_p = B diag(velocity_mass_matrix)^-1 B^T.

I have tested the dependence of the preconditioner on the time step and the mesh size and observed that it is likely to be mesh independent, i.e. if the mesh is refined the number of outer (GMRES) iterations does not change, which is great.

However there is a dependence on the size of the time step. Both variants with assembly and your implementation show an increase in the number of iterations.

time step
        iterations with B B^T
        iterations with assembly
1e-5
        6 to 7
        6 or 16 ringing
1e-4
        8 to 9
        18 to 19
1e-3
        10 to 11
        21 to 22
1e-2
        12 to 14
        24 to 27
1e-1
        23 to 24
        51 to 54


Did you make the same or similar observations? I have read some of the literature and as far as I understand it the preconditioner should be independent of the step size and viscosity.

Another point I would like to raise is related to  mesh independence, which is not quite true. In every step the solution of the pressure laplace matrix (K_p) requires by far the largest number of iterations. I am monitoring the iterations of all subsolvers:

     7 GMRES iterations for stokes system,  (A: 8, Kp: 682, Mp: 60)

For A (velocity block) one cycle of TrilinosWrappers::PreconditionAMG is used. Mp is preconditioned using PreconditionSSOR and for Kp I am using SparseILU. For the example quoted the pressure space has 4800 dofs. I am wondering why the SparseILU with roughly 100 iterations per precondition step does such a bad job in this case.

Of course, I don't expect you to answer all of this but it would great if you make some comments about your experience with this preconditioner. By the way I am also using an IMEX scheme in my solver, see https://github.com/sebglane/BoussinesqFluidSolver/tree/modular_version .

Best wishes,

Sebastian

On 01.11.18 17:30, Jie Cheng wrote:
Hi Sebastian,

I think what you described is correct. I could not see why it did not work out. But I recommend reading Timo's dissertation. Also, for the implementation, you can check out my code for the implicit scheme <https://github.com/OpenIFEM/OpenIFEM/blob/master/source/insim.cpp#L56> and explicit-implicit scheme <https://github.com/OpenIFEM/OpenIFEM/blob/master/source/insimex.cpp#L42>.

Best,
Jie

On Thu, Nov 1, 2018 at 5:58 PM SebG <seb_gl...@gmx.de <mailto:seb_gl...@gmx.de>> wrote:

    Dear Jie, dear deai.ii user,

    I am working on the Cahouet-Chabbard preconditioner in context
    buoyancy-driven flow problems. Somehow my preconditioner does not
    work quite well. The number of iterations depends on the time
    step, which should not be the case. With more than 50 iterations
    it is also quite large. I would like to ask if you or someone else
    could provide some details of your implementation or give tips.

    The velocity block of the system matrix is given by:
    alpha / timestep * M + gamma * c * K .
    M, K are the velocity mass and stiffness matrix. The scalars alpha
    and gamma are related to the time discretization and c is a
    dimensionless parameter. If I am not wrong, the Cahouet-Chabbard
    Schur complement approximation is given by
    S^-1 = alpha / timestep * K_p^-1 + gamma * c * M_p^-1 .

    I am assembling the pressure stiffness and pressure mass matrix
    explicitly. However my problem is a pure Dirichlet problem w.r.t.
    the velocity, which in contrast mean that preconditioner is using
    Neumann BCs. Therefore, I am constraining one pressure dof, which
    regularizes the pressure laplace matrix. This approach is
    discussed in another thread
    <https://groups.google.com/forum/#!msg/dealii/HBQD_WXuNOs/w0A56WAf2WcJ>.
    For this reason I have two ConstraintMatrix object one for the
    system matrix and one for the preconditioner.

    I also attached my code which is based step-32 but in serial.

    Best wishes,
    Sebastian

    Am Freitag, 20. Oktober 2017 09:36:07 UTC-7 schrieb Jie Cheng:

        Hi Martin and Wolfgang

          Thank you very much for the helpful comments and references.
        I'll start to read the works you mentioned.

            Jie -- I think this would be a very interesting program
            for others as well.
            Would you be interested in contributing it to the code
            gallery?


          I'd love to contribute after I work it out!

        Jie

-- The deal.II project is located at http://www.dealii.org/
    For mailing list/forum options, see
    https://groups.google.com/d/forum/dealii?hl=en
    ---
    You received this message because you are subscribed to a topic in
    the Google Groups "deal.II User Group" group.
    To unsubscribe from this topic, visit
    https://groups.google.com/d/topic/dealii/41VjIh5dzng/unsubscribe.
    To unsubscribe from this group and all its topics, send an email
    to dealii+unsubscr...@googlegroups.com
    <mailto:dealii+unsubscr...@googlegroups.com>.
    For more options, visit https://groups.google.com/d/optout.

--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en
---
You received this message because you are subscribed to a topic in the Google Groups "deal.II User Group" group. To unsubscribe from this topic, visit https://groups.google.com/d/topic/dealii/41VjIh5dzng/unsubscribe. To unsubscribe from this group and all its topics, send an email to dealii+unsubscr...@googlegroups.com <mailto:dealii+unsubscr...@googlegroups.com>.
For more options, visit https://groups.google.com/d/optout.

--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- You received this message because you are subscribed to the Google Groups "deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to dealii+unsubscr...@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.

Reply via email to