[gmx-developers] std::vector<Rvec> slowness

Berk Hess hess at kth.se
Wed Mar 23 10:44:26 CET 2016


On 2016-03-23 10:42, Mark Abraham wrote:
> Hi,
>
> On Wed, Mar 23, 2016 at 9:44 AM Berk Hess <hess at kth.se 
> <mailto:hess at kth.se>> wrote:
>
>     Hi,
>
>     Luckily Szilard does thorough testing and noticed a performance
>     degradation in change set 25 of
>     https://gerrit.gromacs.org/#/c/5232/ The
>     only signifcant change with respect to previous sets is replacing C
>     pointers by std::vector. I traced the performance difference back to a
>     single loop, which must have become several factors slower to explain
>     the time difference. I get the performance back when replacing the
>     vector by a pointer extracted with .data(), see below. I looked at the
>     assembly code from gcc 5.3.1 and the vector case generated 200 extra
>     instructions, which makes it difficult to see what the actual
>     difference
>     is. The pointer case uses a lot of vmovss and vaddss, which the vector
>     one does much less, but this is only serial SIMD instruction. I
>     thought
>     that [] in vector might does bounds checking,
>
>
> Definitely not in release builds.
>
>     but it seems it does not.
>     Can anyone explain why the vector case can be so slow?
>
>     If this is a general issue (with RVec or more?), we need to always
>     extra
>     a pointer with .data() for use in all inner-loops. This is pretty
>     annoying and difficult to enforce.
>
>     Cheers,
>
>     Berk
>
>                             const std::vector<RVec> f_foreign =
>     idt_foreign->force
>
>
> This does a copy of the vector, and doesn't seem to be in any version 
> of this patch in gerrit. Is this what you meant to write?
I tried this. But my original "vectorized" patch set took a pointer from 
idt_foreign and did not copy the vector, that gives the same, slow, 
performance.

Berk
>
> Mark
>
>     or
>                             const RVec               *f_foreign  =
>     idt_foreign->force.data();
>
>                              int natom = atomList->atom.size();
>                              for (int i = 0; i < natom; i++)
>                              {
>                                  int ind = atomList->atom[i];
>                                  rvec_inc(f[ind], f_foreign[ind]);
>                              }
>     --
>     Gromacs Developers mailing list
>
>     * Please search the archive at
>     http://www.gromacs.org/Support/Mailing_Lists/GMX-developers_List
>     before posting!
>
>     * Can't post? Read http://www.gromacs.org/Support/Mailing_Lists
>
>     * For (un)subscribe requests visit
>     https://maillist.sys.kth.se/mailman/listinfo/gromacs.org_gmx-developers
>     or send a mail to gmx-developers-request at gromacs.org
>     <mailto:gmx-developers-request at gromacs.org>.
>
>
>

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://maillist.sys.kth.se/pipermail/gromacs.org_gmx-developers/attachments/20160323/fc2d9d5b/attachment.html>


More information about the gromacs.org_gmx-developers mailing list