<html><head><meta http-equiv="Content-Type" content="text/html; charset=utf-8"></head><body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><br class=""><div><br class=""><blockquote type="cite" class=""><div class="">On Oct 16, 2020, at 11:33 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" class="">knepley@gmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><div dir="ltr" class=""><div dir="ltr" class="">On Fri, Oct 16, 2020 at 11:48 PM Alexey Kozlov <<a href="mailto:Alexey.V.Kozlov.2@nd.edu" class="">Alexey.V.Kozlov.2@nd.edu</a>> wrote:<br class=""></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 dir="ltr" class=""><p class="MsoNormal" style="margin:0in 0in 10pt;line-height:115%;font-size:11pt;font-family:Calibri,sans-serif">Thank you for your advice! My sparse matrix seems to be very
stiff so I have decided to concentrate on the direct solvers. I have very good
results with MUMPS. Due to a lack of time I haven’t got a good result with SuperLU_DIST
and haven’t compiled PETSc with Pastix yet but I have a feeling that MUMPS is
the best. I have run a sequential test case with built-in PETSc LU (-pc_type lu
-ksp_type preonly) and MUMPs (-pc_type lu -ksp_type preonly
-pc_factor_mat_solver_type mumps) with default settings and found that MUMPs was
about 50 times faster than the built-in LU and used about 3 times less RAM. Do
you have any idea why it could be?</p></div></blockquote><div class="">The numbers do not sound realistic, but of course we do not have your particular problem. In particular, the memory figure seems impossible. </div></div></div></div></blockquote><div><br class=""></div>   They are probably using a different ordering. Remember each external direct solver manages its own orderings and doesn't share even their names. (Not nice community behavior).</div><div><br class=""><blockquote type="cite" class=""><div class=""><div dir="ltr" class=""><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 dir="ltr" class=""><div style="margin: 0in 0in 10pt; line-height: 115%; font-size: 11pt; font-family: Calibri, sans-serif;" class=""> <br class="webkit-block-placeholder"></div><p class="MsoNormal" style="margin:0in 0in 10pt;line-height:115%;font-size:11pt;font-family:Calibri,sans-serif">My test case has about 100,000 complex equations with about 3,000,000
non-zeros. PETSc was compiled with the following options: ./configure
--with-blaslapack-dir=/opt/crc/i/intel/19.0/mkl --enable-g
--with-valgrind-dir=/opt/crc/v/valgrind/3.14/ompi --with-scalar-type=complex
--with-clanguage=c --with-openmp --with-debugging=0 COPTFLAGS='-mkl=parallel
-O2 -mavx -axCORE-AVX2 -no-prec-div -fp-model fast=2' FOPTFLAGS='-mkl=parallel
-O2 -mavx -axCORE-AVX2 -no-prec-div -fp-model fast=2' CXXOPTFLAGS='-mkl=parallel
-O2 -mavx -axCORE-AVX2 -no-prec-div -fp-model fast=2' --download-superlu_dist
--download-mumps --download-scalapack --download-metis --download-cmake
--download-parmetis --download-ptscotch. </p><p class="MsoNormal" style="margin:0in 0in 10pt;line-height:115%;font-size:11pt;font-family:Calibri,sans-serif">Running MUPMS in parallel using MPI also gave me a significant
gain in performance (about 10 times on a single cluster node).</p></div></blockquote><div class="">Again, this does not appear to make sense. The performance should be limited by memory bandwidth, and a single cluster node will not usually have</div><div class="">10x the bandwidth of a CPU, although it might be possible with a very old CPU.<br class=""></div><div class=""><br class=""></div><div class="">It would help to understand the performance if you would send the output of -log_view.</div><div class=""><br class=""></div><div class="">  Thanks,</div><div class=""><br class=""></div><div class="">    Matt</div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr" class=""><div style="margin: 0in 0in 10pt; line-height: 115%; font-size: 11pt; font-family: Calibri, sans-serif;" class=""> <br class="webkit-block-placeholder"></div>

<span style="font-size:11pt;line-height:115%;font-family:Calibri,sans-serif" class="">Could you, please, advise me whether I can adjust
some options for the direct solvers to improve performance? Should I try MUMPS
in OpenMP mode?</span>



</div></blockquote></div></div></div></blockquote><div><br class=""></div>    Look at the orderings and other options that MUMPs supports and try them out. Most can be accessed from the command line. You can run with -help to get a real brief summary of them but should read the MUMPs users manual for full details.</div><div><br class=""></div><div><br class=""><blockquote type="cite" class=""><div class=""><div dir="ltr" class=""><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"><br class=""><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Sat, Sep 19, 2020 at 7:40 AM Mark Adams <<a href="mailto:mfadams@lbl.gov" target="_blank" class="">mfadams@lbl.gov</a>> wrote:<br class=""></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr" class="">As Jed said high frequency is hard. AMG, as-is,  can be adapted (<a href="https://link.springer.com/article/10.1007/s00466-006-0047-8" target="_blank" class="">https://link.springer.com/article/10.1007/s00466-006-0047-8</a>) with parameters.<div class="">AMG for convection: use richardson/sor and not chebyshev smoothers and in smoothed aggregation (gamg) don't smooth (-pc_gamg_agg_nsmooths 0).</div><div class="">Mark</div></div><br class=""><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Sat, Sep 19, 2020 at 2:11 AM Alexey Kozlov <<a href="mailto:Alexey.V.Kozlov.2@nd.edu" target="_blank" class="">Alexey.V.Kozlov.2@nd.edu</a>> wrote:<br class=""></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr" class="">Thanks a lot! I'll check them out.<br class=""></div><br class=""><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Sat, Sep 19, 2020 at 1:41 AM Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank" class="">bsmith@petsc.dev</a>> wrote:<br class=""></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=""><div class=""><br class=""></div>  These are small enough that likely sparse direct solvers are the best use of your time and for general efficiency. <div class=""><br class=""></div><div class="">  PETSc supports 3 parallel direct solvers, SuperLU_DIST, MUMPs and Pastix. I recommend configuring PETSc for all three of them and then comparing them for problems of interest to you.</div><div class=""><br class=""></div><div class="">   --download-superlu_dist --download-mumps --download-pastix --download-scalapack (used by MUMPS) --download-metis --download-parmetis --download-ptscotch </div><div class=""><br class=""></div><div class="">  Barry</div><div class=""><br class=""><div class=""><br class=""><blockquote type="cite" class=""><div class="">On Sep 18, 2020, at 11:28 PM, Alexey Kozlov <<a href="mailto:Alexey.V.Kozlov.2@nd.edu" target="_blank" class="">Alexey.V.Kozlov.2@nd.edu</a>> wrote:</div><br class=""><div class=""><div dir="ltr" class="">Thanks for the tips! My matrix is complex and unsymmetric. My typical test case has of the order of one million equations. I use a 2nd-order finite-difference scheme with 19-point stencil, so my typical test case uses several GB of RAM.<br class=""></div><br class=""><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Fri, Sep 18, 2020 at 11:52 PM Jed Brown <<a href="mailto:jed@jedbrown.org" target="_blank" class="">jed@jedbrown.org</a>> wrote:<br class=""></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">Unfortunately, those are hard problems in which the "good" methods are technical and hard to make black-box.  There are "sweeping" methods that solve on 2D "slabs" with PML boundary conditions, H-matrix based methods, and fancy multigrid methods.  Attempting to solve with STRUMPACK is probably the easiest thing to try (--download-strumpack).<br class="">
<br class="">
<a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MATSOLVERSSTRUMPACK.html" rel="noreferrer" target="_blank" class="">https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MATSOLVERSSTRUMPACK.html</a><br class="">
<br class="">
Is the matrix complex symmetric?<br class="">
<br class="">
Note that you can use a direct solver (MUMPS, STRUMPACK, etc.) for a 3D problem like this if you have enough memory.  I'm assuming the memory or time is unacceptable and you want an iterative method with much lower setup costs.<br class="">
<br class="">
Alexey Kozlov <<a href="mailto:Alexey.V.Kozlov.2@nd.edu" target="_blank" class="">Alexey.V.Kozlov.2@nd.edu</a>> writes:<br class="">
<br class="">
> Dear all,<br class="">
><br class="">
> I am solving a convected wave equation in a frequency domain. This equation<br class="">
> is a 3D Helmholtz equation with added first-order derivatives and mixed<br class="">
> derivatives, and with complex coefficients. The discretized PDE results in<br class="">
> a sparse linear system (about 10^6 equations) which is solved in PETSc. I<br class="">
> am having difficulty with the code convergence at high frequency, skewed<br class="">
> grid, and high Mach number. I suspect it may be due to the preconditioner I<br class="">
> use. I am currently using the ILU preconditioner with the number of fill<br class="">
> levels 2 or 3, and BCGS or GMRES solvers. I suspect the state of the art<br class="">
> has evolved and there are better preconditioners for Helmholtz-like<br class="">
> problems. Could you, please, advise me on a better preconditioner?<br class="">
><br class="">
> Thanks,<br class="">
> Alexey<br class="">
><br class="">
> -- <br class="">
> Alexey V. Kozlov<br class="">
><br class="">
> Research Scientist<br class="">
> Department of Aerospace and Mechanical Engineering<br class="">
> University of Notre Dame<br class="">
><br class="">
> 117 Hessert Center<br class="">
> Notre Dame, IN 46556-5684<br class="">
> Phone: (574) 631-4335<br class="">
> Fax: (574) 631-8355<br class="">
> Email: <a href="mailto:akozlov@nd.edu" target="_blank" class="">akozlov@nd.edu</a><br class="">
</blockquote></div><br clear="all" class=""><br class="">-- <br class=""><div dir="ltr" class=""><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><span style="color:rgb(0,0,255)" class=""><font size="4" class=""><span style="font-family:"comic sans ms",sans-serif" class="">Alexey V. Kozlov</span></font><br class=""><br class="">Research Scientist<br class="">Department of Aerospace and Mechanical Engineering<br class="">University of Notre Dame<br class=""><br class="">117 Hessert Center<br class="">Notre Dame, IN 46556-5684<br class="">Phone: (574) 631-4335<br class="">Fax: (574) 631-8355<br class="">Email: <a href="mailto:akozlov@nd.edu" target="_blank" class="">akozlov@nd.edu</a></span><br class=""></div></div></div></div></div></div>
</div></blockquote></div><br class=""></div></div></blockquote></div><br clear="all" class=""><br class="">-- <br class=""><div dir="ltr" class=""><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><span style="color:rgb(0,0,255)" class=""><font size="4" class=""><span style="font-family:"comic sans ms",sans-serif" class="">Alexey V. Kozlov</span></font><br class=""><br class="">Research Scientist<br class="">Department of Aerospace and Mechanical Engineering<br class="">University of Notre Dame<br class=""><br class="">117 Hessert Center<br class="">Notre Dame, IN 46556-5684<br class="">Phone: (574) 631-4335<br class="">Fax: (574) 631-8355<br class="">Email: <a href="mailto:akozlov@nd.edu" target="_blank" class="">akozlov@nd.edu</a></span><br class=""></div></div></div></div></div></div>
</blockquote></div>
</blockquote></div><br clear="all" class=""><br class="">-- <br class=""><div dir="ltr" class=""><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><span style="color:rgb(0,0,255)" class=""><font size="4" class=""><span style="font-family:"comic sans ms",sans-serif" class="">Alexey V. Kozlov</span></font><br class=""><br class="">Research Scientist<br class="">Department of Aerospace and Mechanical Engineering<br class="">University of Notre Dame<br class=""><br class="">117 Hessert Center<br class="">Notre Dame, IN 46556-5684<br class="">Phone: (574) 631-4335<br class="">Fax: (574) 631-8355<br class="">Email: <a href="mailto:akozlov@nd.edu" target="_blank" class="">akozlov@nd.edu</a></span><br class=""></div></div></div></div></div></div>
</blockquote></div><br clear="all" class=""><div class=""><br class=""></div>-- <br class=""><div dir="ltr" class="gmail_signature"><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><div class="">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br class="">-- Norbert Wiener</div><div class=""><br class=""></div><div class=""><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank" class="">https://www.cse.buffalo.edu/~knepley/</a><br class=""></div></div></div></div></div></div></div></div>
</div></blockquote></div><br class=""></body></html>