<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="">Dear Viktor,<br class=""><div><br class=""><blockquote type="cite" class=""><div class="">On 1 Sep 2021, at 10:42 AM, Наздрачёв Виктор <<a href="mailto:numbersixvs@gmail.com" class="">numbersixvs@gmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><div dir="ltr" class=""><p class="MsoNormal" style="margin:0cm 0cm 8pt;line-height:107%;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size:14pt;line-height:107%" class="">Dear all,</span></p><p class="MsoNormal" style="text-align:justify;margin:0cm 0cm 8pt;line-height:107%;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size:14pt;line-height:107%" class="">I have a 3D
elasticity problem with heterogeneous properties. There is unstructured grid
with aspect ratio varied from 4 to 25. Zero Dirichlet BCs  are imposed on bottom face of mesh. Also,
Neumann (traction) BCs are imposed on side faces. Gravity load is also
accounted for. The grid I use consists of 500k cells (which is approximately
1.6M of DOFs).</span></p><p class="MsoNormal" style="text-align:justify;margin:0cm 0cm 8pt;line-height:107%;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size:14pt;line-height:107%" class="">The best
performance and memory usage for single MPI process was obtained with
HPDDM(BFBCG) solver</span></p></div></div></blockquote><div>Block Krylov solvers are (most often) only useful if you have multiple right-hand sides, e.g., in the context of elasticity, multiple loadings.</div><div>Is that really the case? If not, you may as well stick to “standard” CG instead of the breakdown-free block (BFB) variant.</div><br class=""><blockquote type="cite" class=""><div class=""><div dir="ltr" class=""><p class="MsoNormal" style="text-align:justify;margin:0cm 0cm 8pt;line-height:107%;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size:14pt;line-height:107%" class=""> and bjacobian + ICC (1) in subdomains as preconditioner, it
took 1 m 45 s and RAM 5.0 GB. Parallel computation with 4 MPI processes took 2
m 46 s when using 5.6 GB of RAM. This because of number of iterations required
to achieve the same tolerance is significantly increased.</span></p><p class="MsoNormal" style="text-align:justify;margin:0cm 0cm 8pt;line-height:107%;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size:14pt;line-height:107%" class="">I`ve also
tried PCGAMG (agg) preconditioner with IC</span><span style="font-size:14pt;line-height:107%" class="">С</span><span lang="EN-US" style="font-size:14pt;line-height:107%" class=""> (1) sub-precondtioner. For single MPI process, the
calculation took 10 min and 3.4 GB of RAM. To improve the convergence rate, the
nullspace was attached using MatNullSpaceCreateRigidBody and MatSetNearNullSpace
subroutines.  This has reduced
calculation time to 3 m 58 s when using 4.3 GB of RAM. Also, there is peak
memory usage with 14.1 GB, which appears just before the start of the
iterations. Parallel computation with 4 MPI processes took 2 m 53 s when using
8.4 GB of RAM. In that case the peak memory usage is about 22 GB.</span></p></div></div></blockquote><div>I’m surprised that GAMG is converging so slowly. What do you mean by "ICC(1) sub-preconditioner"? Do you use that as a smoother or as a coarse level solver?</div><div>How many iterations are required to reach convergence?</div><div>Could you please maybe run the solver with -ksp_view -log_view and send us the output?</div><div>Most of the default parameters of GAMG should be good enough for 3D elasticity, provided that your MatNullSpace is correct.</div><div>One parameter that may need some adjustments though is the aggregation threshold -pc_gamg_threshold (you could try values in the [0.01; 0.1] range, that’s what I always use for elasticity problems).</div><div><br class=""></div><div>Thanks,</div><div>Pierre</div><div><br class=""></div><blockquote type="cite" class=""><div class=""><div dir="ltr" class=""><p class="MsoNormal" style="text-align:justify;margin:0cm 0cm 8pt;line-height:107%;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size:14pt;line-height:107%" class="">Are there
ways to avoid decreasing of the convergence rate for bjacobi precondtioner in
parallel mode? Does it make sense to use hierarchical or nested krylov methods
with a local gmres solver (sub_pc_type gmres) and some sub-precondtioner (for
example, sub_pc_type bjacobi)?</span></p><div style="text-align: justify; margin: 0cm 0cm 8pt; line-height: 107%; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size:14pt;line-height:107%" class=""> </span><br class="webkit-block-placeholder"></div><p class="MsoNormal" style="text-align:justify;margin:0cm 0cm 8pt;line-height:107%;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size:14pt;line-height:107%" class="">Is this peak
memory usage expected for gamg preconditioner? is there any way to reduce it?</span></p><div style="text-align: justify; margin: 0cm 0cm 8pt; line-height: 107%; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size:14pt;line-height:107%" class=""> </span><br class="webkit-block-placeholder"></div><p class="MsoNormal" style="text-align:justify;margin:0cm 0cm 8pt;line-height:107%;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size:14pt;line-height:107%" class="">What advice
would you give to improve the convergence rate with multiple MPI processes, but
keep memory consumption reasonable?</span></p><div style="margin: 0cm 0cm 8pt; line-height: 107%; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size:14pt;line-height:107%" class=""> </span><br class="webkit-block-placeholder"></div><p class="MsoNormal" style="margin:0cm 0cm 8pt;line-height:107%;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size:14pt;line-height:107%" class="">Kind regards,</span></p><p class="MsoNormal" style="margin:0cm 0cm 8pt;line-height:107%;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size:14pt;line-height:107%" class="">Viktor Nazdrachev</span></p><p class="MsoNormal" style="margin:0cm 0cm 8pt;line-height:107%;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size:14pt;line-height:107%" class="">R&D senior researcher</span></p><p class="MsoNormal" style="margin:0cm 0cm 8pt;line-height:107%;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size:14pt;line-height:107%" class="">Geosteering Technologies LLC</span></p></div>
</div></blockquote></div><br class=""></body></html>