<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Mon, Feb 8, 2016 at 5:37 AM, Dave May <span dir="ltr"><<a href="mailto:dave.mayhem23@gmail.com" target="_blank">dave.mayhem23@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr"><br><div class="gmail_extra"><br><div class="gmail_quote"><span class="">On 8 February 2016 at 12:31, Jacek Miloszewski <span dir="ltr"><<a href="mailto:jacek.miloszewski@gmail.com" target="_blank">jacek.miloszewski@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-style:solid;border-left-color:rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div><span style="font-size:12.8px">Dear PETSc users,</span></div><div><span style="font-size:12.8px"><br></span></div><div><span style="font-size:12.8px">I use PETSc to assembly a square matrix (in the attached example it is n = 4356) which has around 12% of non-zero entries. I timed my code using various number of process (data in table). Now I have 2 questions:</span></div><div><span style="font-size:12.8px"><br></span></div><div><span style="font-size:12.8px">1. Why with doubling number of processes the speed-up is 4x? I would expect 2x at the most.</span></div></div></blockquote><div><br></div></span><div>Your timing doesn't appear to include the time required to scatter off-processor values.<br></div><div>You should move the timer to be after the following calls<br><br></div><div>MatAssemblyBegin(a,MAT_FINAL_ASSEMBLY);<br></div><div><div>MatAssemblyEnd(a,MAT_FINAL_ASSEMBLY);</div></div></div></div></div></blockquote><div><br></div><div>Also, I would guess that you are not preallocating the matrix:</div><div><br></div><div> <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html#efficient-assembly">http://www.mcs.anl.gov/petsc/documentation/faq.html#efficient-assembly</a></div><div><br></div><div> Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr"><div class="gmail_extra"><div class="gmail_quote"><div> <br></div><div><div class="h5"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-style:solid;border-left-color:rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div><span style="font-size:12.8px"><br></span></div><div><span style="font-size:12.8px">2. Is there a way to speed-up matrix construction in general? I attach piece of my fortran code at the bottom. At compilation time I have the following knowledge about the matrix: the total number of non-zero matrix elements, all diagonal elements are non-zero.<br></span><br>Timing data:<br><br><table dir="ltr" style="table-layout:fixed;font-size:13px;font-family:arial,sans,sans-serif;border-collapse:collapse;border:1px solid rgb(204,204,204)" border="1" cellpadding="0" cellspacing="0"><colgroup><col width="100"><col width="100"><col width="100"></colgroup><tbody><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;border:1px solid rgb(0,0,0)">number of proc</td><td style="padding:2px 3px;vertical-align:bottom;border-top-width:1px;border-top-style:solid;border-top-color:rgb(0,0,0);border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0)">time [s]</td><td style="padding:2px 3px;vertical-align:bottom;border-top-width:1px;border-top-style:solid;border-top-color:rgb(0,0,0);border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0)">speedup ratio</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);border-left-width:1px;border-left-style:solid;border-left-color:rgb(0,0,0);text-align:right">1</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">2044.941</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0)"></td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);border-left-width:1px;border-left-style:solid;border-left-color:rgb(0,0,0);text-align:right">2</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">504.692</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">4.051859352</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);border-left-width:1px;border-left-style:solid;border-left-color:rgb(0,0,0);text-align:right">4</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">149.678</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">3.371851575</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);border-left-width:1px;border-left-style:solid;border-left-color:rgb(0,0,0);text-align:right">8</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">64.102</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">2.334997348</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);border-left-width:1px;border-left-style:solid;border-left-color:rgb(0,0,0);text-align:right">16</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">17.296</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">3.706174838</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);border-left-width:1px;border-left-style:solid;border-left-color:rgb(0,0,0);text-align:right">32</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">4.43</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">3.904288939</td></tr><tr style="height:21px"><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);border-left-width:1px;border-left-style:solid;border-left-color:rgb(0,0,0);text-align:right">64</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">1.096</td><td style="padding:2px 3px;vertical-align:bottom;border-right-width:1px;border-right-style:solid;border-right-color:rgb(0,0,0);border-bottom-width:1px;border-bottom-style:solid;border-bottom-color:rgb(0,0,0);text-align:right">4.041970803</td></tr></tbody></table></div><div><font color="#134f5c"><br></font></div><div><font color="#000000">Code:</font></div><div><font color="#134f5c"><br></font></div><div><font color="#000000">call MatCreate(PETSC_COMM_WORLD, a, ierr)</font></div><div><font color="#000000">call MatSetSizes(a, PETSC_DECIDE, PETSC_DECIDE, n, n, ierr)</font></div><div><font color="#000000">call MatSetFromOptions(a, ierr)</font></div><div><font color="#000000">call MatSetUp(a, ierr)</font></div><div><font color="#000000"><br></font></div><div><font color="#000000">call MatGetOwnershipRange(a, Istart, Iend, ierr)</font></div><div><font color="#000000"><br></font></div><div><font color="#000000">call system_clock(t1)</font></div><div><font color="#000000">t_dim = h_dim*e_dim</font></div><div><font color="#000000">do row = Istart, Iend - 1 ! row</font></div><div><font color="#000000"> do col = 0, t_dim - 1 ! col</font></div><div><font color="#000000"> call h_ij(row + 1, col + 1, n_h, n_e, b_h, b_e, h_dim, e_dim, e_sp, v, basis, ht, hs, info_1)</font></div><div><font color="#000000"> if (hs) then</font></div><div><font color="#000000"> hh = ht</font></div><div><font color="#000000"> call MatSetValues(a, 1, row, 1, col, hh, INSERT_VALUES, ierr)</font></div><div><font color="#000000"> end if</font></div><div><font color="#000000"> end do</font></div><div><font color="#000000">end do</font></div><div><font color="#000000">call system_clock(t2, ct)</font></div><div><font color="#000000">if (rank == 0) then</font></div><div><font color="#000000"> write(*, '(a, f0.3)') 'Matrix assembly time: ', real((t2 - t1), r8)/real(ct, r8)</font></div><div><font color="#000000">end if</font></div><span><font color="#888888"><div><br></div>-- <br><div><div dir="ltr"><div><div dir="ltr"><div>Best Wishes</div><div>Jacek Miloszewski</div></div></div></div></div>
</font></span></div>
</blockquote></div></div></div><br></div></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>