<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Thu, Jun 25, 2015 at 5:51 AM, Lei Shi <span dir="ltr"><<a href="mailto:stoneszone@gmail.com" target="_blank">stoneszone@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div>Hello,</div></div></blockquote><div><br></div><div>1) In order to understand this, we have to disentagle the various effect. First, run the STREAMS benchmark</div><div><br></div><div>  make NPMAX=4 streams</div><div><br></div><div>This will tell you the maximum speedup you can expect on this machine.</div><div><br></div><div>2) For these test cases, also send the output of</div><div><br></div><div>  -ksp_view -ksp_converged_reason -ksp_monitor_true_residual</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div>I'm trying to improve the parallel efficiency of gmres solve in my. In my CFD solver, Petsc gmres is used to solve the linear system generated by the Newton's method. To test its efficiency, I started with a very simple inviscid subsonic 3D flow as the first testcase. The parallel efficiency of gmres solve with asm as the preconditioner is very bad. The results are from our latest cluster. Right now, I'm only looking at the wclock time of the ksp_solve.</div><div><ol><li>First I tested ASM with gmres and ilu 0 for the sub domain , the cpu time of 2 cores is almost the same as the serial run. Here is the options for this case<br></li></ol></div><blockquote style="margin:0px 0px 0px 40px;border:none;padding:0px"><div><div>-ksp_type gmres  -ksp_max_it 100 -ksp_rtol 1e-5 -ksp_atol 1e-50 </div></div><div><div>-ksp_gmres_restart 30 -ksp_pc_side right</div></div><div><div>-pc_type asm -sub_ksp_type gmres -sub_ksp_rtol 0.001 -sub_ksp_atol 1e-30</div></div><div><div>-sub_ksp_max_it 1000 -sub_pc_type ilu -sub_pc_factor_levels 0 </div></div><div><div>-sub_pc_factor_fill 1.9</div></div></blockquote><blockquote style="margin:0px 0px 0px 40px;border:none;padding:0px"><div>The iteration numbers increase a lot for parallel run.</div></blockquote><blockquote style="margin:0px 0px 0px 40px;border:none;padding:0px"><div><table cellspacing="0" cellpadding="0" dir="ltr" border="1" style="table-layout:fixed;font-size:13px;font-family:arial,sans,sans-serif;border-collapse:collapse;border:1px solid rgb(204,204,204)"><colgroup><col width="100"><col width="100"><col width="144"><col width="144"><col width="100"><col width="100"></colgroup><tbody><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom">cores</td><td style="padding:2px 3px;vertical-align:bottom">iterations</td><td style="padding:2px 3px;vertical-align:bottom">err</td><td style="padding:2px 3px;vertical-align:bottom">petsc solve wclock time</td><td style="padding:2px 3px;vertical-align:bottom">speedup</td><td style="padding:2px 3px;vertical-align:bottom">efficiency</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;text-align:right">1</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">2</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">1.15E-04</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">11.95</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">1</td><td></td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;text-align:right">2</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">5</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">2.05E-02</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">10.5</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">1.01</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">0.50</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;text-align:right">4</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">6</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">2.19E-02</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">7.64</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">1.39</td><td style="padding:2px 3px;vertical-align:bottom;text-align:right">0.34</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;text-align:right"><br></td><td style="padding:2px 3px;vertical-align:bottom;text-align:right"><br></td><td style="padding:2px 3px;vertical-align:bottom;text-align:right"><br></td><td style="padding:2px 3px;vertical-align:bottom;text-align:right"><br></td><td style="padding:2px 3px;vertical-align:bottom;text-align:right"><br></td><td style="padding:2px 3px;vertical-align:bottom;text-align:right"><br></td></tr></tbody></table></div></blockquote><div><div>     </div><div>      2.  Then I tested ASM with ilu 0 as the preconditoner only, the cpu time of 2 cores is better than the 1st test, but the speedup is still very bad. Here is the options i'm using<br></div></div><blockquote style="margin:0px 0px 0px 40px;border:none;padding:0px"><div><div>-ksp_type gmres  -ksp_max_it 100 -ksp_rtol 1e-5 -ksp_atol 1e-50 </div></div><div><div>-ksp_gmres_restart 30 -ksp_pc_side right</div></div><div><div>-pc_type asm -sub_pc_type ilu -sub_pc_factor_levels 0  -sub_pc_factor_fill 1.9</div></div></blockquote><blockquote style="margin:0px 0px 0px 40px;border:none;padding:0px"><table cellspacing="0" cellpadding="0" border="1" style="table-layout:fixed;font-size:13px;font-family:arial,sans,sans-serif;border-collapse:collapse;border:1px solid rgb(204,204,204)"><colgroup><col width="100"><col width="100"><col width="144"><col width="144"><col width="100"><col width="100"></colgroup><tbody><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom">cores</td><td style="padding:2px 3px;vertical-align:bottom">iterations</td><td style="padding:2px 3px;vertical-align:bottom">err</td><td style="padding:2px 3px;vertical-align:bottom">petsc solve cpu time</td><td style="padding:2px 3px;vertical-align:bottom">speedup</td><td style="padding:2px 3px;vertical-align:bottom">efficiency</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;text-align:left">1</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">10</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">4.54E-04</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">10.68</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">1</td><td></td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;text-align:left">2</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">11</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">9.55E-04</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">8.2</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">1.30</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">0.65</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;text-align:left">4</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">12</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">3.59E-04</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">5.26</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">2.03</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">0.50</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;text-align:left"><br></td><td style="padding:2px 3px;vertical-align:bottom;text-align:left"><br></td><td style="padding:2px 3px;vertical-align:bottom;text-align:left"><br></td><td style="padding:2px 3px;vertical-align:bottom;text-align:left"><br></td><td style="padding:2px 3px;vertical-align:bottom;text-align:left"><br></td><td style="padding:2px 3px;vertical-align:bottom;text-align:left"><br></td></tr></tbody></table></blockquote><div><br></div><div>   Those results are from a third order "DG" scheme with a very coarse 3D mesh (480 elements). I believe I should get some speedups for this test even on this coarse mesh. </div><div><br></div><div>  My question is why does the asm with a local solve take much longer time than the asm as a preconditioner only? Also the accuracy is very bad too I have tested changing the overlap of asm to 2, but make it even worse.</div><div><br></div><div>  If I used a larger mesh ~4000 elements, the 2nd case with asm as the preconditioner gives me a better speedup, but still not very good. </div><div><blockquote style="margin:0px 0px 0px 40px;border:none;padding:0px"><div><br></div><div><table cellspacing="0" cellpadding="0" border="1" style="table-layout:fixed;font-size:13px;font-family:arial,sans,sans-serif;border-collapse:collapse;border:1px solid rgb(204,204,204)"><colgroup><col width="100"><col width="100"><col width="100"><col width="100"><col width="100"><col width="100"></colgroup><tbody><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom">cores</td><td style="padding:2px 3px;vertical-align:bottom">iterations</td><td style="padding:2px 3px;vertical-align:bottom">err</td><td style="padding:2px 3px;vertical-align:bottom">petsc solve cpu time</td><td style="padding:2px 3px;vertical-align:bottom">speedup</td><td style="padding:2px 3px;vertical-align:bottom">efficiency</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;text-align:left">1</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">7</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">1.91E-02</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">97.32</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">1</td><td></td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;text-align:left">2</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">7</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">2.07E-02</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">64.94</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">1.5</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">0.74</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;text-align:left">4</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">7</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">2.61E-02</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">36.97</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">2.6</td><td style="padding:2px 3px;vertical-align:bottom;text-align:left">0.65</td></tr></tbody></table></div></blockquote></div><div><br></div><div><br></div><div>Attached are the log_summary dumped from petsc, any suggestions are welcome. I really appreciate it.</div><div><br></div><div><br clear="all"><div><div>Sincerely Yours,<br><br>Lei Shi <br>---------</div></div>
</div><img src="https://t.yesware.com/t/d1fcbaa1b12e0f6b1beef0b50d5ebbd873d1b8f9/e8040b69bc7ca00f2bce0335cfcc6ef8/spacer.gif" style="border:0px;width:0px;min-height:0px;overflow:hidden" width="0" height="0"><img src="http://t.yesware.com/t/d1fcbaa1b12e0f6b1beef0b50d5ebbd873d1b8f9/e8040b69bc7ca00f2bce0335cfcc6ef8/spacer.gif" style="border:0px;width:0px;min-height:0px;overflow:hidden" width="0" height="0"><font face="yw-d1fcbaa1b12e0f6b1beef0b50d5ebbd873d1b8f9-e8040b69bc7ca00f2bce0335cfcc6ef8--to"></font></div>
</blockquote></div><br><br clear="all"><div><br></div>-- <br><div class="gmail_signature">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></div>