<div dir="ltr"><div dir="ltr">On Tue, May 2, 2023 at 9:07 AM Blauth, Sebastian <<a href="mailto:sebastian.blauth@itwm.fraunhofer.de">sebastian.blauth@itwm.fraunhofer.de</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div class="msg-3051768603283180690"><div lang="DE" style="overflow-wrap: break-word;"><div class="m_1708809103809953240WordSection1"><p class="MsoNormal">Hello,<u></u><u></u></p><p class="MsoNormal"><u></u> <u></u></p><p class="MsoNormal"><span lang="EN-US">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. <u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">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).<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US"><u></u> <u></u></span></p><p class="MsoNormal"><span lang="EN-US">To test the solver, I am using a simple 2D channel flow problem with do-nothing conditions at the outlet.<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US"><u></u> <u></u></span></p><p class="MsoNormal"><span lang="EN-US">I am facing the following difficulties at the moment:<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US"><u></u> <u></u></span></p><p class="MsoNormal"><span lang="EN-US">- 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)</span></p></div></div></div></blockquote><div><br></div><div>The diagonal of the momentum block is a bad preconditioner for the Schur complement, because S is spectrally equivalent to the mass matrix. You should build the mass matrix and use that as the preconditioning matrix for the Schur part. The FEniCS people can show you how to do that. This will provide mesh-independent convergence (you can see me doing this in SNES ex69).</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div class="msg-3051768603283180690"><div lang="DE" style="overflow-wrap: break-word;"><div class="m_1708809103809953240WordSection1"><p class="MsoNormal"><span lang="EN-US">    -ksp_type fgmres<u></u><u></u></span></p><p class="MsoNormal">    -ksp_rtol 1e-6<u></u><u></u></p><p class="MsoNormal">    <span lang="EN-US">-ksp_atol 1e-30<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -pc_type fieldsplit<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -pc_fieldsplit_type schur<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -pc_fieldsplit_schur_fact_type full<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -pc_fieldsplit_schur_precondition selfp<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -fieldsplit_0_ksp_type preonly<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -fieldsplit_0_pc_type lu<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -fieldsplit_1_ksp_type gmres<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -fieldsplit_1_ksp_pc_side right<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -fieldsplit_1_ksp_max_it 1000<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -fieldsplit_1_ksp_rtol 1e-1<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -fieldsplit_1_ksp_atol 1e-30<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -fieldsplit_1_pc_type lu<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -fieldsplit_1_ksp_converged_reason<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">    -ksp_converged_reason<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US"><u></u> <u></u></span></p><p class="MsoNormal"><span lang="EN-US">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. <u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">In particular, I get the following behavior:<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">For a 8x8 mesh, I need, on average, 25 iterations to solve fieldsplit_1<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">For a 16x16 mesh, I need 40 iterations<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">For a 32x32 mesh, I need 70 iterations<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">For a 64x64 mesh, I need 100 iterations<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US"><u></u> <u></u></span></p><p class="MsoNormal"><span lang="EN-US">However, the outer fgmres requires, as expected, always the same number of iterations to reach convergence (as expected).<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">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?<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US"><u></u> <u></u></span></p><p class="MsoNormal"><span lang="EN-US">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.<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US"><u></u> <u></u></span></p><p class="MsoNormal"><span lang="EN-US">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?<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US"><u></u> <u></u></span></p><p class="MsoNormal"><span lang="EN-US">Thanks a lot in advance and I am looking forward to your reply,<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US">Sebastian<u></u><u></u></span></p><p class="MsoNormal"><span lang="EN-US"><u></u> <u></u></span></p><p class="MsoNormal"><span lang="EN-US">--<u></u><u></u></span></p><p class="MsoNormal"><span>Dr. Sebastian Blauth<u></u><u></u></span></p><p class="MsoNormal"><span>Fraunhofer-Institut für<u></u><u></u></span></p><p class="MsoNormal"><span>Techno- und Wirtschaftsmathematik ITWM<u></u><u></u></span></p><p class="MsoNormal"><span>Abteilung Transportvorgänge<u></u><u></u></span></p><p class="MsoNormal"><span>Fraunhofer-Platz 1, 67663 Kaiserslautern<u></u><u></u></span></p><p class="MsoNormal"><span>Telefon: +49 631 31600-4968<u></u><u></u></span></p><p class="MsoNormal"><span><a href="mailto:sebastian.blauth@itwm.fraunhofer.de" target="_blank">sebastian.blauth@itwm.fraunhofer.de</a><u></u><u></u></span></p><p class="MsoNormal"><span><a href="https://www.itwm.fraunhofer.de" target="_blank">https://www.itwm.fraunhofer.de</a><u></u><u></u></span></p><p class="MsoNormal"><u></u> <u></u></p></div></div></div></blockquote></div><br clear="all"><div><br></div><span class="gmail_signature_prefix">-- </span><br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>