<div dir="ltr"><div>Hi, </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">We use a hash table to store the nonzeros on the fly, and then convert to packed storage on assembly.<br></blockquote><div><br></div><div>Maybe can you tell me which file implements this function? </div><div><br></div><div>Runfeng</div><div><br></div><div><br></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Runfeng Jin <<a href="mailto:jsfaraway@gmail.com">jsfaraway@gmail.com</a>> 于2023年7月3日周一 22:05写道:<br></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">Thank you for all your help!<div><br></div><div>Runfeng</div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> 于2023年7月3日周一 22:03写道:<br></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"><div dir="ltr"><div>On Mon, Jul 3, 2023 at 9:56 AM Runfeng Jin <<a href="mailto:jsfaraway@gmail.com" target="_blank">jsfaraway@gmail.com</a>> wrote:<br></div></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">Hi, impressive performance!<div>  I use the newest version of petsc(release branch), and it almost deletes all assembly and stash time in large processors (assembly time 64-4s/128-2s/256-0.2s, stash time all below 2s). For the zero programming cost, it really incredible. </div><div>  The order code has a regular arrangement of the number of nonzero-elements across rows, so I can have a good rough preallocation. And from the data, dedicatedly arrange data and roughly acquiring the max number of non-zero elements in rows can have a better performance than the new version without preallocation. However, in reality, I will use the newer version without preallocation for:1)less effort in programming and also nearly the same good performance 2) good memory usage(I see no unneeded memory after assembly) 3) dedicated preallocation is usually not very easy and cause extra time cost.</div><div>   Maybe it will be better that leave some space for the user to do a slight direction for the preallocation and thus acquire better performance. But have no idea how to direct it.</div><div>   And I am very curious about how petsc achieves this. How can it not know anything but achieve so good performance, and no wasted memory? May you have an explanation about this?</div></div></blockquote><div><br></div><div>We use a hash table to store the nonzeros on the fly, and then convert to packed storage on assembly.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> </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"><div>assemble time:</div><div>version\processors               4            8        16         32           64        128         256</div><div>     old                             14677s   4694s   1124s     572s        38s         8s          2s</div><div>     new                                50s      28s       15s        7.8s         4s          2s        0.4s</div><div>     older                              27s       24s        19s       12s         14s         -              -</div><div>stash time(max among all processors):</div><div><div>version\processors               4            8        16         32           64        128         256</div><div>     old                                 3145s   2554s   673s     329s       201s     142s     138s</div><div>     new                                2s         1s        ~0s        ~0s         ~0s          ~0s       ~0s</div><div>     older                              10s       73s        18s       5s            1s         -              -</div></div><div>old: my poor preallocation</div><div>new: newest version of petsc that do not preallocation</div><div>older: the best preallocation version of my code. </div><div><br></div><div><br></div><div>Runfeng</div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>> 于2023年7月3日周一 12:19写道:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div><br></div>   The main branch of PETSc now supports filling sparse matrices without providing any preallocation information.<div><br></div><div>   You can give it a try. Use your current fastest code but just remove ALL the preallocation calls. I would be interested in what kind of performance you get compared to your best current performance.</div><div><br></div><div>  Barry</div><div><br></div><div><br><div><br><blockquote type="cite"><div>On Jul 2, 2023, at 11:24 PM, Runfeng Jin <<a href="mailto:jsfaraway@gmail.com" target="_blank">jsfaraway@gmail.com</a>> wrote:</div><br><div><div dir="ltr"><div dir="ltr">Hi! Good advice!<div>    I set value with MatSetValues() API, which sets one part of a row at a time(I use a kind of tiling technology so I cannot get all values of a row at a time).<br>    I tested the number of malloc in these three cases.  The number of mallocs is decreasing with the increase of processors, and all these are very large(the matrix is 283234*283234, as can see in the following). This should be due to the unqualified preallocation. I use a rough preallocation, that every processor counts the number of nonzero elements for the first 10 rows, and uses the largest one to preallocate memory for all local rows. It seems that not work well. </div><blockquote style="margin:0px 0px 0px 40px;border:none;padding:0px">number_of_processors   number_of_max_mallocs_among_all_processors  <br>64                                     20000<br>128                                   17000<br>256                                   11000<br></blockquote>    I change my way to preallocate. I evenly take 100 rows in every local matrix and take the largest one to preallocate memory for all local rows. Now the assemble time is reduced to a very small time.<div><blockquote style="margin:0px 0px 0px 40px;border:none;padding:0px"><div>number_of_processors   number_of_max_mallocs_among_all_processors  <br>64                                     3000<br>128                                   700<br>256                                   500</div></blockquote><blockquote style="margin:0px 0px 0px 40px;border:none;padding:0px"><div><div dir="ltr"><div class="gmail_quote"><div>Event                Count          Time (sec)            Flop                                                              --- Global ---                  --- Stage ----              Total<br>                   Max Ratio        Max     Ratio       Max  Ratio  Mess            AvgLen  Reduct  %T %F %M %L %R        %T %F %M %L %R    Mflop/s<br></div><div>64                 1    1.0       3.8999e+01 1.0     0.00e+00 0.0 7.1e+03     2.9e+05 1.1e+01 15  0  1  8  3                     15  0  1  8  3                  0</div></div></div></div></blockquote><blockquote style="margin:0px 0px 0px 40px;border:none;padding:0px"><div>128               1    1.0       8.5714e+00 1.0     0.00e+00 0.0 2.6e+04     8.1e+04 1.1e+01  5  0  1  4  3                       5  0  1  4  3                   0<br></div><div>256               1    1.0        2.5512e+00 1.0    0.00e+00 0.0 1.0e+05     2.3e+04 1.1e+01  2  0  1  3  3                       2  0  1  3  3                   0<br></div></blockquote><div>So the reason "why assemble time is smaller with the increasing number of processors " may be because more processors divide the malloc job so that total time is reduced? </div> If so, I still have some questions:<div>    1. If preallocation is not accurate, will the performance of the assembly be affected? I mean, when processors receive the elements that should be stored in their local by MPI, then will the new mallocs  happen at this time point?</div><div>    2. I can not give an accurate preallocation for the large cost, so is there any better way to preallocate for my situation?<br><blockquote style="margin:0px 0px 0px 40px;border:none;padding:0px"><br></blockquote></div></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>> 于2023年7月2日周日 00:16写道:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div><br></div>   I see no reason not to trust the times below, they seem reasonable. You get more than 2 times speed from 64 to 128 and then about 1.38 from 128 to 256. <div><br></div><div>   The total amount of data moved (number of messages moved times average length) goes from 7.0e+03 * 2.8e+05  1.9600e+09 to 2.1060e+09 to 2.3000e+09. A pretty moderate amount of data increase, but note that each time you double the number of ranks, you also increase substantially the network's hardware to move data, so one would hope for a good speed up.</div><div><br></div><div>   Also, the load balance is very good, near 1. Often with assembly, we see very out-of-balance, and it is difficult to get good speedup when the balance is really off.</div><div><br></div><div>   It looks like over 90% of the entire run time is coming from setting and assembling the values? Also the setting values time dominates assembly time more with more ranks.  Are you setting a single value at a time or a collection of them? How big are the vectors?</div><div><br></div><div>   Run all three cases with -info :vec to see some information about how many mallocs where move to hold the stashed vector entries.</div><div><br></div><div><br></div><div><br><div><br><blockquote type="cite"><div>On Jun 30, 2023, at 10:25 PM, Runfeng Jin <<a href="mailto:jsfaraway@gmail.com" target="_blank">jsfaraway@gmail.com</a>> wrote:</div><br><div><div dir="ltr"><div class="gmail_quote"><div dir="ltr" class="gmail_attr"><br></div><br><div dir="ltr">Hi, <div>    Thanks for your reply. I try to use PetscLogEvent(), and the result shows same conclusion.</div><div>    What I have done is :</div><div>----------------</div><div>    PetscLogEvent Mat_assemble_event, Mat_setvalue_event, Mat_setAsse_event;<br>    PetscClassId classid;<br>    PetscLogDouble user_event_flops;<br>    PetscClassIdRegister("Test assemble and set value", &classid);<br>    PetscLogEventRegister("Test only assemble", classid, &Mat_assemble_event);<br>    PetscLogEventRegister("Test only set values", classid, &Mat_setvalue_event);<br>    PetscLogEventRegister("Test both assemble and set values", classid, &Mat_setAsse_event);<br>    PetscLogEventBegin(Mat_setAsse_event, 0, 0, 0, 0);<br>    PetscLogEventBegin(Mat_setvalue_event, 0, 0, 0, 0);<br></div><div>    ...compute elements and use MatSetValues. No call for assembly</div><div>    PetscLogEventEnd(Mat_setvalue_event, 0, 0, 0, 0);<br><br>    PetscLogEventBegin(Mat_assemble_event, 0, 0, 0, 0);<br>    MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);<br>    MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);<br>    PetscLogEventEnd(Mat_assemble_event, 0, 0, 0, 0);<br>    PetscLogEventEnd(Mat_setAsse_event, 0, 0, 0, 0);<br></div><div>----------------<br></div><div><br></div><div>    And the output as follows. By the way, dose petsc recorde all time between PetscLogEventBegin and PetscLogEventEnd? or just test the time of petsc API?</div></div></div></div></div></blockquote><div><br></div>   It is all of the time. </div><div><br><blockquote type="cite"><div><div dir="ltr"><div class="gmail_quote"><div dir="ltr"><div>----------------<br></div>Event                Count      Time (sec)     Flop                              --- Global ---  --- Stage ----  Total<br>                   Max Ratio  <b>Max</b>     Ratio   Max  Ratio  Mess   AvgLen  Reduct  %T %F %M %L %R  %T %F %M %L %R Mflop/s<br>64new               1 1.0 <b>2.3775e+02</b> 1.0 0.00e+00 0.0 6.2e+03 2.3e+04 9.0e+00 52  0  1  1  2  52  0  1  1  2     0<br>128new              1 1.0<b> 6.9945e+01</b> 1.0 0.00e+00 0.0 2.5e+04 1.1e+04 9.0e+00 30  0  1  1  2  30  0  1  1  2     0<br>256new              1 1.0 <b>1.7445e+01</b> 1.0 0.00e+00 0.0 9.9e+04 5.2e+03 9.0e+00 10  0  1  1  2  10  0  1  1  2     0<br><br>64:<br>only assemble       1 1.0 <b>2.6596e+02 </b>1.0 0.00e+00 0.0 7.0e+03 2.8e+05 1.1e+01 55  0  1  8  3  55  0  1  8  3     0<br>only setvalues      1 1.0 <b>1.9987e+02</b> 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 41  0  0  0  0  41  0  0  0  0     0<br>Test both           1 1.0 4.<b>6580e+02</b> 1.0 0.00e+00 0.0 7.0e+03 2.8e+05 1.5e+01 96  0  1  8  4  96  0  1  8  4     0<br><br>128:<br> only assemble      1 1.0 <b>6.9718e+01</b> 1.0 0.00e+00 0.0 2.6e+04 8.1e+04 1.1e+01 30  0  1  4  3  30  0  1  4  3     0<br>only setvalues      1 1.0 <b>1.4438e+02</b> 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 60  0  0  0  0  60  0  0  0  0     0<br>Test both           1 1.0 <b>2.1417e+02</b> 1.0 0.00e+00 0.0 2.6e+04 8.1e+04 1.5e+01 91  0  1  4  4  91  0  1  4  4     0<br><br>256:<br>only assemble       1 1.0 <b>1.7482e+01</b> 1.0 0.00e+00 0.0 1.0e+05 2.3e+04 1.1e+01 10  0  1  3  3  10  0  1  3  3     0<br>only setvalues      1 1.0 <b>1.3717e+02</b> 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 78  0  0  0  0  78  0  0  0  0     0<br><div>Test both           1 1.0 <b>1.5475e+02</b> 1.0 0.00e+00 0.0 1.0e+05 2.3e+04 1.5e+01 91  0  1  3  4  91  0  1  3  4     0 </div><div><br></div><div><br></div><div><br></div><div>Runfeng</div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>> 于2023年6月30日周五 23:35写道:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><br>
   You cannot look just at the VecAssemblyEnd() time, that will very likely give the wrong impression of the total time it takes to put the values in.<br>
<br>
   You need to register a new Event and put a PetscLogEvent() just before you start generating the vector entries and calling VecSetValues() and put the PetscLogEventEnd() just after the VecAssemblyEnd() this is the only way to get an accurate accounting of the time.<br>
<br>
  Barry<br>
<br>
<br>
> On Jun 30, 2023, at 11:21 AM, Runfeng Jin <<a href="mailto:jsfaraway@gmail.com" target="_blank">jsfaraway@gmail.com</a>> wrote:<br>
> <br>
> Hello!<br>
> <br>
> When I use PETSc build a sbaij matrix, I find a strange thing. When I increase the number of processors, the assemble time become smaller. All these are totally same matrix. The assemble time mainly arouse from message passing, which because I use dynamic workload that it is random for which elements are computed by which processor.<br>
> But from instinct, if use more processors, then more possible that the processor computes elements storing in other processors. But from the output of log_view, It seems when use more processors, the processors compute more elements storing in its local(infer from that, with more processors, less total amount of passed messages).<br>
> <br>
> What could cause this happened? Thank you!<br>
> <br>
> <br>
>  Following is the output of log_view for 64\128\256 processors. Every row is time profiler of VecAssemblyEnd.<br>
> <br>
> ------------------------------------------------------------------------------------------------------------------------<br>
> processors                Count                      Time (sec)                                      Flop                                                               --- Global ---                               --- Stage ----                Total<br>
>                               Max    Ratio         Max                  Ratio                 Max  Ratio      Mess        AvgLen         Reduct               %T %F %M %L %R         %T %F %M %L %R       Mflop/s<br>
> 64                            1     1.0            2.3775e+02      1.0                   0.00e+00 0.0      6.2e+03    2.3e+04     9.0e+00                 52  0      1    1    2             52   0    1      1     2             0<br>
> 128                          1     1.0            6.9945e+01      1.0                   0.00e+00 0.0      2.5e+04    1.1e+04     9.0e+00                30   0      1     1  2              30   0    1       1    2             0<br>
> 256                          1     1.0           1.7445e+01        1.0                  0.00e+00 0.0      9.9e+04     5.2e+03    9.0e+00                10   0      1     1  2              10   0    1        1   2             0<br>
> <br>
> Runfeng Jin<br>
<br>
</blockquote></div>
</div></div>
</div></blockquote></div><br></div></div></blockquote></div></div>
</div></blockquote></div><br></div></div></blockquote></div>
</blockquote></div><br clear="all"><div><br></div><span class="gmail_signature_prefix">-- </span><br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>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><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>
</blockquote></div>
</blockquote></div></div>