[petsc-users] Scalable Solver for Incompressible Flow
    Blauth, Sebastian 
    sebastian.blauth at itwm.fraunhofer.de
       
    Tue May  2 03:52:16 CDT 2023
    
    
  
Hello,
 
I am having a problem using / configuring PETSc to obtain a scalable solver
for the incompressible Navier Stokes equations. I am discretizing the
equations using FEM (with the library fenics) and I am using the stable
P2-P1 Taylor-Hood elements. I have read and tried a lot regarding
preconditioners for incompressible Navier Stokes and I am aware that this is
very much an active research field, but maybe I can get some hints / tips. 
I am interested in solving large-scale 3D problems, but I cannot even set up
a scaleable 2D solver for the problems. All of my approaches at the moment
are trying to use a Schur Complement approach, but I cannot get a good
preconditioner for the Schur complement matrix. For the velocity block, I am
using the AMG provided by hypre (which seems to work fine and is most likely
not the problem).
 
To test the solver, I am using a simple 2D channel flow problem with
do-nothing conditions at the outlet.
 
I am facing the following difficulties at the moment:
 
- First, I am having trouble with using -pc_fieldsplit_schur_precondition
selfp. With this setup, the cost for solving the Schur complement part in
the fieldsplit preconditioner (approximately) increase when the mesh is
refined. I am using the following options for this setup (note that I am
using exact solves for the velocity part to debug, but using, e.g., gmres
with hypre boomeramg reaches a given tolerance with a number of iterations
that is independent of the mesh)
 
    -ksp_type fgmres
    -ksp_rtol 1e-6
    -ksp_atol 1e-30
    -pc_type fieldsplit
    -pc_fieldsplit_type schur
    -pc_fieldsplit_schur_fact_type full
    -pc_fieldsplit_schur_precondition selfp
    -fieldsplit_0_ksp_type preonly
    -fieldsplit_0_pc_type lu
    -fieldsplit_1_ksp_type gmres
    -fieldsplit_1_ksp_pc_side right
    -fieldsplit_1_ksp_max_it 1000
    -fieldsplit_1_ksp_rtol 1e-1
    -fieldsplit_1_ksp_atol 1e-30
    -fieldsplit_1_pc_type lu
    -fieldsplit_1_ksp_converged_reason
    -ksp_converged_reason
 
Note, that I use direct solvers for the subproblems to get an ideal
convergence. Even if I replace the direct solver with boomeramg, the
behavior is the same and the number of iterations does not change much. 
In particular, I get the following behavior:
For a 8x8 mesh, I need, on average, 25 iterations to solve fieldsplit_1
For a 16x16 mesh, I need 40 iterations
For a 32x32 mesh, I need 70 iterations
For a 64x64 mesh, I need 100 iterations
 
However, the outer fgmres requires, as expected, always the same number of
iterations to reach convergence (as expected).
I do understand that the selfp preconditioner for the Schur complement is
expected to deteriorate as the Reynolds number increases and the problem
becomes more convective in nature, but I had hoped that I can at least get a
scaleable preconditioner with respect to the mesh size out of it. Are there
any tips on how to achieve this?
 
My second problem is concerning the LSC preconditioner. When I am using
this, again both with exact solves of the linear problems or when using
boomeramg, I do not get a scalable solver with respect to the mesh size. On
the contrary, here the number of solves required for solving fieldsplit_1 to
a fixed relative tolerance seem to behave linearly w.r.t. the problem size.
For this problem, I suspect that the issue lies in the scaling of the LSC
preconditioner matrices (in the book of Elman, Sylvester and Wathen, the
matrices are scaled with the inverse of the diagonal velocity mass matrix).
Is it possible to achieve this with PETSc? I started experimenting with
supplying the velocity mass matrix as preconditioner matrix and using
use_amat, but I am not sure where / how to do it this way.
 
And finally, more of an observation and question: I noticed that the AMG
approximations for the velocity block became worse with increase of the
Reynolds number when using the default options. However, when using
-pc_hypre_boomeramg_relax_weight_all 0.0 I noticed that boomeramg performed
way more robustly w.r.t. the Reynolds number. Are there any other ways to
improve the AMG performance in this regard?
 
Thanks a lot in advance and I am looking forward to your reply,
Sebastian
 
--
Dr. Sebastian Blauth
Fraunhofer-Institut für
Techno- und Wirtschaftsmathematik ITWM
Abteilung Transportvorgänge
Fraunhofer-Platz 1, 67663 Kaiserslautern
Telefon: +49 631 31600-4968
sebastian.blauth at itwm.fraunhofer.de
https://www.itwm.fraunhofer.de
 
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20230502/79bbe7cf/attachment-0001.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: smime.p7s
Type: application/pkcs7-signature
Size: 7943 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20230502/79bbe7cf/attachment-0001.p7s>
    
    
More information about the petsc-users
mailing list