<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=""><div class=""><br class=""></div> You can use -mat_null_space_test to check it the null space you provide is within the null space of the operator. There is no practical way to test if the null space you provide is exactly the full null space of the operator but at least the check ensures that you are not providing something that is not in the null space.<div class=""><br class=""></div><div class=""> Barry</div><div class=""><br class=""></div><div class=""> Also if you run with GMRES and look at the norms of the residuals -ksp_monitor at each restart iteration (the default restart is 30), if they jump wildly at each restart this can indicate a problem with the nullspace. </div><div class=""><br class=""></div><div class=""><br class=""><div><br class=""><blockquote type="cite" class=""><div class="">On Sep 3, 2021, at 10:48 AM, Viktor Nazdrachev <<a href="mailto:numbersixvs@gmail.com" class="">numbersixvs@gmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><div dir="ltr" class=""><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">Hello Mark and Matthew!</span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">I attached log files for serial and parallel cases and corresponding information about GAMG preconditioner (using grep).
I have to notice, that assembling of global stiffness matrix in code was performed by MatSetValues subrotuine (not MatSetValuesBlocked) <br class="">
</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">!nnds – number of nodes</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">!dmn=3</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call MatCreate(Petsc_Comm_World,Mat_K,ierr)</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call MatSetFromOptions(Mat_K,ierr)</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call MatSetSizes(Mat_K,Petsc_Decide,Petsc_Decide,n,n,ierr_m)</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">…</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call MatMPIAIJSetPreallocation(Mat_K,0,dbw,0,obw,ierr)</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">…</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call MatSetOption(Mat_K,Mat_New_Nonzero_Allocation_Err,Petsc_False,ierr)</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">…<br class="">
do i=1,nels</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> call FormLocalK(i,k,indx,"Kp") ! find local stiffness matrix</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> indx=indxmap(indx,2) !find global indices for DOFs</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> call MatSetValues(Mat_K,ef_eldof,indx,ef_eldof,indx,k,Add_Values,ierr) </span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">end do</span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">But nullspace vector was created using VecSetBlockSize subroutine.</span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call VecCreate(Petsc_Comm_World,Vec_NullSpace,ierr)</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call VecSetBlockSize(Vec_NullSpace,dmn,ierr)</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call VecSetSizes(Vec_NullSpace,nnds*dmn,Petsc_Decide,ierr)</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call VecSetUp(Vec_NullSpace,ierr)</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call VecGetArrayF90(Vec_NullSpace,null_space,ierr)</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">…</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call VecRestoreArrayF90(Vec_NullSpace,null_space,ierr)</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call MatNullSpaceCreateRigidBody(Vec_NullSpace,matnull,ierr)</span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">call MatSetNearNullSpace(Mat_K,matnull,ierr)</span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">I suppose it can be one of the reasons of GAMG slow convergence.</span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">So I attached log files for parallel run with «pure» GAMG precondtioner.</span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">Kind regards,</span></div><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">Viktor Nazdrachev</span></div><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">R&D senior researcher</span></div><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">Geosteering Technologies LLC</span></div></div><br class=""><div class="gmail_quote"><div dir="ltr" class="gmail_attr">пт, 3 сент. 2021 г. в 15:11, Matthew Knepley <<a href="mailto:knepley@gmail.com" class="">knepley@gmail.com</a>>:<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=""><div dir="ltr" class="">On Fri, Sep 3, 2021 at 8:02 AM Mark Adams <<a href="mailto:mfadams@lbl.gov" target="_blank" class="">mfadams@lbl.gov</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=""><div dir="ltr" class=""><br class=""></div><br class=""><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Fri, Sep 3, 2021 at 1:57 AM Viktor Nazdrachev <<a href="mailto:numbersixvs@gmail.com" target="_blank" class="">numbersixvs@gmail.com</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=""><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">Hello, Lawrence!<br class="">
Thank you for your response!</span></div><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">I attached log files (txt files with convergence
behavior and RAM usage log in separate txt files) and resulting table with
convergence investigation data(xls). Data for main non-regular grid with 500K cells
and heterogeneous properties are in 500K folder, whereas data for simple
uniform 125K cells grid with constant properties are in 125K folder. </span></div><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""><br class=""></span></div><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>><i class=""> On 1 Sep 2021, at 09:42, </i></span><i class=""><span style="" class="">Наздрачёв</span></i><i class=""><span style="" class=""> </span><span style="" class="">Виктор</span></i><i class=""><span lang="EN-US" style="" class=""> <</span><span style="" class=""><a href="https://lists.mcs.anl.gov/mailman/listinfo/petsc-users" style="color:blue" target="_blank" class=""><span lang="EN-US" style="" class="">numbersixvs at gmail.com</span></a></span></i><i class=""><span lang="EN-US" style="" class="">> wrote:</span></i></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>><i class=""> </i></span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>><i class=""> I have a 3D elasticity problem with heterogeneous properties.</i></span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">> </span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>What does your coefficient variation look like? How large is the contrast?</span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">Young modulus varies from 1 to 10 GPa, Poisson ratio varies
from 0.3 to 0.44 and density – from 1700 to 2600 kg/m^3.</span></div></div></blockquote><div class=""><br class=""></div><div class="">That is not too bad. Poorly shaped elements are the next thing to worry about. Try to keep the aspect ratio below 10 if possible.</div><div 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=""><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>><i class=""> 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).</i></span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>><i class=""> </i></span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>><i class=""> The best performance and memory usage for single MPI process was obtained with HPDDM(BFBCG) solver 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.</i></span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">> </span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>How many iterations do you have in serial (and then in parallel)?</span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">Serial run is required 112 iterations to reach convergence (log_hpddm(bfbcg)_bjacobian_icc_1_mpi.txt), parallel run with 4 MPI – 680 iterations.</span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">I attached log files for all simulations (txt files
with convergence behavior and RAM usage log in separate txt files) and
resulting table with convergence/memory usage data(xls). Data for main
non-regular grid with 500K cells and heterogeneous properties are in 500K
folder, whereas data for simple uniform 125K cells grid with constant
properties are in 125K folder. </span></div><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>><i class=""> I`ve also tried PCGAMG (agg) preconditioner with IC</i></span><i class=""><span style="" class="">С</span></i><i class=""><span lang="EN-US" style="" 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></i></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">> </span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>Does the number of iterates increase in parallel? Again, how many iterations do you have?</span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="text-align:justify;margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">For case with 4 MPI processes and attached nullspace it is required 177 iterations to reach convergence (you may see detailed log in log_hpddm(bfbcg)_gamg_nearnullspace_4_mpi.txt). For comparison, 90 iterations are required for sequential run(log_hpddm(bfbcg)_gamg_nearnullspace_1_mpi.txt).</span></pre></div></blockquote><div class=""><br class=""></div><div class="">Again, do not use ICC. I am surprised to see such a large jump in iteration count, but get ICC off the table.</div><div class=""><br class=""></div><div class="">You will see variability in the iteration count with processor count with GAMG. As much as 10% +-. Maybe more (random) variability , but usually less.</div><div class=""><br class=""></div><div class="">You can decrease the memory a little, and the setup time a lot, by aggressively coarsening, at the expense of higher iteration counts. It's a balancing act.</div><div class=""><br class=""></div><div class="">You can run with the defaults, add '-info', grep on GAMG and send the ~30 lines of output if you want advice on parameters.</div></div></div></blockquote><div class=""><br class=""></div><div class="">Can you send the output of</div><div class=""><br class=""></div><div class=""> -ksp_view -ksp_monitor_true_residual -ksp_converged_reason</div><div class=""><br class=""></div><div class=""> Thanks,</div><div class=""><br class=""></div><div class=""> Matt</div><div 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=""><div class="gmail_quote"><div class="">Thanks,</div><div class="">Mark</div><div 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=""><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>><i 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)?</i></span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">> </span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>bjacobi is only a one-level method, so you would not expect process-independent convergence rate for this kind of problem. If the coefficient variation is not too extreme, then I would expect GAMG (or some other smoothed aggregation package, perhaps -pc_type ml (you need --download-ml)) would work well with some tuning.</span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">Thanks for idea, but, unfortunately, ML cannot be
compiled with 64bit integers (It is extremely necessary to perform computation
on mesh with more than 10M cells).</span></div><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>If you have extremely high contrast coefficients you might need something with stronger coarse grids. If you can assemble so-called Neumann matrices (</span><span style="" class=""><a href="https://petsc.org/release/docs/manualpages/Mat/MATIS.html#MATIS" style="color:blue" target="_blank" class=""><span lang="EN-US" style="" class="">https://petsc.org/release/docs/manualpages/Mat/MATIS.html#MATIS</span></a></span><span lang="EN-US" style="" class="">) then you could try the geneo scheme offered by PCHPDDM.</span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">I found strange convergence behavior for HPDDM
preconditioner. For 1 MPI process BFBCG solver did not converged
(log_hpddm(bfbcg)_pchpddm_1_mpi.txt), while for 4 MPI processes computation was
successful (1018 to reach convergence, log_hpddm(bfbcg)_pchpddm_4_mpi.txt).</span></div><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">But it should be mentioned that stiffness matrix was
created in AIJ format (our default matrix format in program).</span></div><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">Matrix conversion to MATIS format via MatConvert </span><span lang="EN-US" style="font-family: "Courier New";" class="">subroutine resulted in losing of convergence for both
serial and parallel run.</span><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""></span></div><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-family: "Courier New";" class=""><br class=""></span></div><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>><i class=""> Is this peak memory usage expected for gamg preconditioner? is there any way to reduce it?</i></span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">> </span></pre><pre style="margin:0cm 0cm 0.0001pt 35.4pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class="">>I think that peak memory usage comes from building the coarse grids. Can you run with `-info` and grep for GAMG, this will provide some output that more expert GAMG users can interpret.</span> </pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span></pre><pre style="margin:0cm 0cm 0.0001pt;font-size:10pt;font-family:"Courier New"" class=""><span lang="EN-US" style="" class=""> </span>Thanks, I`ll try to use a strong threshold only for coarse grids.</pre><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">Kind regards,</span></div><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">Viktor Nazdrachev</span></div><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">R&D senior researcher</span></div><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><div style="margin: 0cm; line-height: normal; font-size: 11pt; font-family: Calibri, sans-serif;" class=""><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">Geosteering Technologies LLC</span></div><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class="">
</span></p><p class="MsoNormal" style="margin:0cm;line-height:normal;font-size:11pt;font-family:Calibri,sans-serif"><span lang="EN-US" style="font-size: 10pt; font-family: "Courier New";" class=""> </span></p></div><br class=""><div class="gmail_quote"><div dir="ltr" class="gmail_attr">ср, 1 сент. 2021 г. в 12:02, Lawrence Mitchell <<a href="mailto:wence@gmx.li" target="_blank" class="">wence@gmx.li</a>>:<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"><br class="">
<br class="">
> On 1 Sep 2021, at 09:42, Наздрачёв Виктор <<a href="mailto:numbersixvs@gmail.com" target="_blank" class="">numbersixvs@gmail.com</a>> wrote:<br class="">
> <br class="">
> I have a 3D elasticity problem with heterogeneous properties.<br class="">
<br class="">
What does your coefficient variation look like? How large is the contrast?<br class="">
<br class="">
> 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).<br class="">
> <br class="">
> The best performance and memory usage for single MPI process was obtained with HPDDM(BFBCG) solver 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.<br class="">
<br class="">
How many iterations do you have in serial (and then in parallel)?<br class="">
<br class="">
> I`ve also tried PCGAMG (agg) preconditioner with ICС (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.<br class="">
<br class="">
Does the number of iterates increase in parallel? Again, how many iterations do you have?<br class="">
<br 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)?<br class="">
<br class="">
bjacobi is only a one-level method, so you would not expect process-independent convergence rate for this kind of problem. If the coefficient variation is not too extreme, then I would expect GAMG (or some other smoothed aggregation package, perhaps -pc_type ml (you need --download-ml)) would work well with some tuning.<br class="">
<br class="">
If you have extremely high contrast coefficients you might need something with stronger coarse grids. If you can assemble so-called Neumann matrices (<a href="https://petsc.org/release/docs/manualpages/Mat/MATIS.html#MATIS" rel="noreferrer" target="_blank" class="">https://petsc.org/release/docs/manualpages/Mat/MATIS.html#MATIS</a>) then you could try the geneo scheme offered by PCHPDDM.<br class="">
<br class="">
> Is this peak memory usage expected for gamg preconditioner? is there any way to reduce it?<br class="">
<br class="">
I think that peak memory usage comes from building the coarse grids. Can you run with `-info` and grep for GAMG, this will provide some output that more expert GAMG users can interpret.<br class="">
<br class="">
Lawrence<br class="">
<br class="">
</blockquote></div>
</blockquote></div></div>
</blockquote></div><br clear="all" class=""><div class=""><br class=""></div>-- <br class=""><div dir="ltr" class=""><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>
</blockquote></div>
<span id="cid:f_kt4h221g0"><true_residual_logs_and_greped.rar></span></div></blockquote></div><br class=""></div></body></html>