2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
47 #include "gromacs/domdec/domdec.h"
48 #include "gromacs/domdec/domdec_struct.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdlib/gmx_omp_nthreads.h"
54 #include "gromacs/mdtypes/commrec.h"
55 #include "gromacs/mdtypes/mdatom.h"
56 #include "gromacs/pbcutil/ishift.h"
57 #include "gromacs/pbcutil/mshift.h"
58 #include "gromacs/pbcutil/pbc.h"
59 #include "gromacs/timing/wallcycle.h"
60 #include "gromacs/topology/ifunc.h"
61 #include "gromacs/topology/mtop_util.h"
62 #include "gromacs/topology/topology.h"
63 #include "gromacs/utility/exceptions.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/gmxassert.h"
66 #include "gromacs/utility/gmxomp.h"
68 /* The strategy used here for assigning virtual sites to (thread-)tasks
71 * We divide the atom range that vsites operate on (natoms_local with DD,
72 * 0 - last atom involved in vsites without DD) equally over all threads.
74 * Vsites in the local range constructed from atoms in the local range
75 * and/or other vsites that are fully local are assigned to a simple,
78 * Vsites that are not assigned after using the above criterion get assigned
79 * to a so called "interdependent" thread task when none of the constructing
80 * atoms is a vsite. These tasks are called interdependent, because one task
81 * accesses atoms assigned to a different task/thread.
82 * Note that this option is turned off with large (local) atom counts
83 * to avoid high memory usage.
85 * Any remaining vsites are assigned to a separate master thread task.
90 static void init_ilist(t_ilist *ilist)
92 for (int i = 0; i < F_NRE; i++)
96 ilist[i].iatoms = nullptr;
100 /*! \brief List of atom indices belonging to a task */
102 //! List of atom indices
103 std::vector<int> atom;
106 /*! \brief Data structure for thread tasks that use constructing atoms outside their own atom range */
107 struct InterdependentTask
109 //! The interaction lists, only vsite entries are used
110 t_ilist ilist[F_NRE];
111 //! Thread/task-local force buffer
112 std::vector<RVec> force;
113 //! The atom indices of the vsites of our task
114 std::vector<int> vsite;
115 //! Flags if elements in force are spread to or not
116 std::vector<bool> use;
117 //! The number of entries set to true in use
119 //! Array of atoms indices, size nthreads, covering all nuse set elements in use
120 std::vector<AtomIndex> atomIndex;
121 //! List of tasks (force blocks) this task spread forces to
122 std::vector<int> spreadTask;
123 //! List of tasks that write to this tasks force block range
124 std::vector<int> reduceTask;
133 /*! \brief Vsite thread task data structure */
136 //! Start of atom range of this task
138 //! End of atom range of this task
140 //! The interaction lists, only vsite entries are used
141 t_ilist ilist[F_NRE];
142 //! Local fshift accumulation buffer
144 //! Local virial dx*df accumulation buffer
146 //! Tells if interdependent task idTask should be used (in addition to the rest of this task), this bool has the same value on all threads
147 bool useInterdependentTask;
148 //! Data for vsites that involve constructing atoms in the atom range of other threads/tasks
149 InterdependentTask idTask;
151 /*! \brief Constructor */
157 clear_rvecs(SHIFTS, fshift);
159 useInterdependentTask = false;
163 /*! \brief Returns the sum of the vsite ilist sizes over all vsite types
165 * \param[in] ilist The interaction list
167 template <typename T>
168 static int vsiteIlistNrCount(const T *ilist)
171 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
173 nr += ilist[ftype].size();
179 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
183 return pbc_dx_aiuc(pbc, xi, xj, dx);
187 rvec_sub(xi, xj, dx);
192 /* Vsite construction routines */
194 static void constr_vsite2(const rvec xi, const rvec xj, rvec x, real a, const t_pbc *pbc)
202 pbc_dx_aiuc(pbc, xj, xi, dx);
203 x[XX] = xi[XX] + a*dx[XX];
204 x[YY] = xi[YY] + a*dx[YY];
205 x[ZZ] = xi[ZZ] + a*dx[ZZ];
209 x[XX] = b*xi[XX] + a*xj[XX];
210 x[YY] = b*xi[YY] + a*xj[YY];
211 x[ZZ] = b*xi[ZZ] + a*xj[ZZ];
215 /* TOTAL: 10 flops */
218 static void constr_vsite3(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b,
228 pbc_dx_aiuc(pbc, xj, xi, dxj);
229 pbc_dx_aiuc(pbc, xk, xi, dxk);
230 x[XX] = xi[XX] + a*dxj[XX] + b*dxk[XX];
231 x[YY] = xi[YY] + a*dxj[YY] + b*dxk[YY];
232 x[ZZ] = xi[ZZ] + a*dxj[ZZ] + b*dxk[ZZ];
236 x[XX] = c*xi[XX] + a*xj[XX] + b*xk[XX];
237 x[YY] = c*xi[YY] + a*xj[YY] + b*xk[YY];
238 x[ZZ] = c*xi[ZZ] + a*xj[ZZ] + b*xk[ZZ];
242 /* TOTAL: 17 flops */
245 static void constr_vsite3FD(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b,
251 pbc_rvec_sub(pbc, xj, xi, xij);
252 pbc_rvec_sub(pbc, xk, xj, xjk);
255 /* temp goes from i to a point on the line jk */
256 temp[XX] = xij[XX] + a*xjk[XX];
257 temp[YY] = xij[YY] + a*xjk[YY];
258 temp[ZZ] = xij[ZZ] + a*xjk[ZZ];
261 c = b*gmx::invsqrt(iprod(temp, temp));
264 x[XX] = xi[XX] + c*temp[XX];
265 x[YY] = xi[YY] + c*temp[YY];
266 x[ZZ] = xi[ZZ] + c*temp[ZZ];
269 /* TOTAL: 34 flops */
272 static void constr_vsite3FAD(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, const t_pbc *pbc)
275 real a1, b1, c1, invdij;
277 pbc_rvec_sub(pbc, xj, xi, xij);
278 pbc_rvec_sub(pbc, xk, xj, xjk);
281 invdij = gmx::invsqrt(iprod(xij, xij));
282 c1 = invdij * invdij * iprod(xij, xjk);
283 xp[XX] = xjk[XX] - c1*xij[XX];
284 xp[YY] = xjk[YY] - c1*xij[YY];
285 xp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
287 b1 = b*gmx::invsqrt(iprod(xp, xp));
290 x[XX] = xi[XX] + a1*xij[XX] + b1*xp[XX];
291 x[YY] = xi[YY] + a1*xij[YY] + b1*xp[YY];
292 x[ZZ] = xi[ZZ] + a1*xij[ZZ] + b1*xp[ZZ];
295 /* TOTAL: 63 flops */
298 static void constr_vsite3OUT(const rvec xi, const rvec xj, const rvec xk, rvec x,
299 real a, real b, real c, const t_pbc *pbc)
303 pbc_rvec_sub(pbc, xj, xi, xij);
304 pbc_rvec_sub(pbc, xk, xi, xik);
305 cprod(xij, xik, temp);
308 x[XX] = xi[XX] + a*xij[XX] + b*xik[XX] + c*temp[XX];
309 x[YY] = xi[YY] + a*xij[YY] + b*xik[YY] + c*temp[YY];
310 x[ZZ] = xi[ZZ] + a*xij[ZZ] + b*xik[ZZ] + c*temp[ZZ];
313 /* TOTAL: 33 flops */
316 static void constr_vsite4FD(const rvec xi, const rvec xj, const rvec xk, const rvec xl, rvec x,
317 real a, real b, real c, const t_pbc *pbc)
319 rvec xij, xjk, xjl, temp;
322 pbc_rvec_sub(pbc, xj, xi, xij);
323 pbc_rvec_sub(pbc, xk, xj, xjk);
324 pbc_rvec_sub(pbc, xl, xj, xjl);
327 /* temp goes from i to a point on the plane jkl */
328 temp[XX] = xij[XX] + a*xjk[XX] + b*xjl[XX];
329 temp[YY] = xij[YY] + a*xjk[YY] + b*xjl[YY];
330 temp[ZZ] = xij[ZZ] + a*xjk[ZZ] + b*xjl[ZZ];
333 d = c*gmx::invsqrt(iprod(temp, temp));
336 x[XX] = xi[XX] + d*temp[XX];
337 x[YY] = xi[YY] + d*temp[YY];
338 x[ZZ] = xi[ZZ] + d*temp[ZZ];
341 /* TOTAL: 43 flops */
344 static void constr_vsite4FDN(const rvec xi, const rvec xj, const rvec xk, const rvec xl, rvec x,
345 real a, real b, real c, const t_pbc *pbc)
347 rvec xij, xik, xil, ra, rb, rja, rjb, rm;
350 pbc_rvec_sub(pbc, xj, xi, xij);
351 pbc_rvec_sub(pbc, xk, xi, xik);
352 pbc_rvec_sub(pbc, xl, xi, xil);
365 rvec_sub(ra, xij, rja);
366 rvec_sub(rb, xij, rjb);
372 d = c*gmx::invsqrt(norm2(rm));
375 x[XX] = xi[XX] + d*rm[XX];
376 x[YY] = xi[YY] + d*rm[YY];
377 x[ZZ] = xi[ZZ] + d*rm[ZZ];
380 /* TOTAL: 47 flops */
384 static int constr_vsiten(const t_iatom *ia, const t_iparams ip[],
385 rvec *x, const t_pbc *pbc)
392 n3 = 3*ip[ia[0]].vsiten.n;
395 copy_rvec(x[ai], x1);
397 for (int i = 3; i < n3; i += 3)
400 a = ip[ia[i]].vsiten.a;
403 pbc_dx_aiuc(pbc, x[ai], x1, dx);
407 rvec_sub(x[ai], x1, dx);
409 dsum[XX] += a*dx[XX];
410 dsum[YY] += a*dx[YY];
411 dsum[ZZ] += a*dx[ZZ];
415 x[av][XX] = x1[XX] + dsum[XX];
416 x[av][YY] = x1[YY] + dsum[YY];
417 x[av][ZZ] = x1[ZZ] + dsum[ZZ];
422 /*! \brief PBC modes for vsite construction and spreading */
425 all, // Apply normal, simple PBC for all vsites
426 none // No PBC treatment needed
429 /*! \brief Returns the PBC mode based on the system PBC and vsite properties
431 * \param[in] pbcPtr A pointer to a PBC struct or nullptr when no PBC treatment is required
433 static PbcMode getPbcMode(const t_pbc *pbcPtr)
435 if (pbcPtr == nullptr)
437 return PbcMode::none;
445 static void construct_vsites_thread(rvec x[],
447 const t_iparams ip[], const t_ilist ilist[],
448 const t_pbc *pbc_null)
460 const PbcMode pbcMode = getPbcMode(pbc_null);
461 /* We need another pbc pointer, as with charge groups we switch per vsite */
462 const t_pbc *pbc_null2 = pbc_null;
464 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
466 if (ilist[ftype].nr == 0)
472 int nra = interaction_function[ftype].nratoms;
474 int nr = ilist[ftype].nr;
476 const t_iatom *ia = ilist[ftype].iatoms;
478 for (int i = 0; i < nr; )
481 /* The vsite and constructing atoms */
484 /* Constants for constructing vsites */
485 real a1 = ip[tp].vsite.a;
486 /* Copy the old position */
488 copy_rvec(x[avsite], xv);
490 /* Construct the vsite depending on type */
497 constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
503 constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
509 constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
515 constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
522 constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
530 constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
539 constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
543 inc = constr_vsiten(ia, ip, x, pbc_null2);
546 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
547 ftype, __FILE__, __LINE__);
550 if (pbcMode == PbcMode::all)
552 /* Keep the vsite in the same periodic image as before */
554 int ishift = pbc_dx_aiuc(pbc_null, x[avsite], xv, dx);
555 if (ishift != CENTRAL)
557 rvec_add(xv, dx, x[avsite]);
562 /* Calculate velocity of vsite... */
564 rvec_sub(x[avsite], xv, vv);
565 svmul(inv_dt, vv, v[avsite]);
568 /* Increment loop variables */
576 void construct_vsites(const gmx_vsite_t *vsite,
579 const t_iparams ip[], const t_ilist ilist[],
580 int ePBC, gmx_bool bMolPBC,
584 const bool useDomdec = (vsite != nullptr && vsite->useDomdec);
585 GMX_ASSERT(!useDomdec || (cr != nullptr && DOMAINDECOMP(cr)), "When vsites are set up with domain decomposition, we need a valid commrec");
586 // TODO: Remove this assertion when we remove charge groups
587 GMX_ASSERT(vsite != nullptr || ePBC == epbcNONE, "Without a vsite struct we can not do PBC (in case we have charge groups)");
589 t_pbc pbc, *pbc_null;
591 /* We only need to do pbc when we have inter-cg vsites.
592 * Note that with domain decomposition we do not need to apply PBC here
593 * when we have at least 3 domains along each dimension. Currently we
594 * do not optimize this case.
596 if (ePBC != epbcNONE && (useDomdec || bMolPBC) &&
597 !(vsite != nullptr && vsite->numInterUpdategroupVsites == 0))
599 /* This is wasting some CPU time as we now do this multiple times
603 clear_ivec(null_ivec);
604 pbc_null = set_pbc_dd(&pbc, ePBC,
605 useDomdec ? cr->dd->nc : null_ivec,
615 dd_move_x_vsites(cr->dd, box, x);
618 if (vsite == nullptr || vsite->nthreads == 1)
620 construct_vsites_thread(x, dt, v,
626 #pragma omp parallel num_threads(vsite->nthreads)
630 const int th = gmx_omp_get_thread_num();
631 const VsiteThread &tData = *vsite->tData[th];
632 GMX_ASSERT(tData.rangeStart >= 0, "The thread data should be initialized before calling construct_vsites");
634 construct_vsites_thread(x, dt, v,
637 if (tData.useInterdependentTask)
639 /* Here we don't need a barrier (unlike the spreading),
640 * since both tasks only construct vsites from particles,
641 * or local vsites, not from non-local vsites.
643 construct_vsites_thread(x, dt, v,
644 ip, tData.idTask.ilist,
648 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
650 /* Now we can construct the vsites that might depend on other vsites */
651 construct_vsites_thread(x, dt, v,
652 ip, vsite->tData[vsite->nthreads]->ilist,
657 static void spread_vsite2(const t_iatom ia[], real a,
659 rvec f[], rvec fshift[],
660 const t_pbc *pbc, const t_graph *g)
671 svmul(1 - a, f[av], fi);
672 svmul( a, f[av], fj);
681 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
683 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
688 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
689 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
697 if (fshift && (siv != CENTRAL || sij != CENTRAL))
699 rvec_inc(fshift[siv], f[av]);
700 rvec_dec(fshift[CENTRAL], fi);
701 rvec_dec(fshift[sij], fj);
704 /* TOTAL: 13 flops */
707 void constructVsitesGlobal(const gmx_mtop_t &mtop,
708 gmx::ArrayRef<gmx::RVec> x)
710 GMX_ASSERT(x.ssize() >= mtop.natoms, "x should contain the whole system");
711 GMX_ASSERT(!mtop.moleculeBlockIndices.empty(), "molblock indices are needed in constructVsitesGlobal");
713 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
715 const gmx_molblock_t &molb = mtop.molblock[mb];
716 const gmx_moltype_t &molt = mtop.moltype[molb.type];
717 if (vsiteIlistNrCount(molt.ilist.data()) > 0)
719 int atomOffset = mtop.moleculeBlockIndices[mb].globalAtomStart;
720 for (int mol = 0; mol < molb.nmol; mol++)
722 t_ilist ilist[F_NRE];
723 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
725 ilist[ftype].nr = molt.ilist[ftype].size();
726 ilist[ftype].iatoms = const_cast<t_iatom *>(molt.ilist[ftype].iatoms.data());
729 construct_vsites(nullptr, as_rvec_array(x.data()) + atomOffset,
731 mtop.ffparams.iparams.data(), ilist,
732 epbcNONE, TRUE, nullptr, nullptr);
733 atomOffset += molt.atoms.nr;
739 static void spread_vsite3(const t_iatom ia[], real a, real b,
741 rvec f[], rvec fshift[],
742 const t_pbc *pbc, const t_graph *g)
754 svmul(1 - a - b, f[av], fi);
755 svmul( a, f[av], fj);
756 svmul( b, f[av], fk);
766 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ia[1]), di);
768 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
770 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), di);
775 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
776 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
777 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
786 if (fshift && (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL))
788 rvec_inc(fshift[siv], f[av]);
789 rvec_dec(fshift[CENTRAL], fi);
790 rvec_dec(fshift[sij], fj);
791 rvec_dec(fshift[sik], fk);
794 /* TOTAL: 20 flops */
797 static void spread_vsite3FD(const t_iatom ia[], real a, real b,
799 rvec f[], rvec fshift[],
800 gmx_bool VirCorr, matrix dxdf,
801 const t_pbc *pbc, const t_graph *g)
803 real c, invl, fproj, a1;
804 rvec xvi, xij, xjk, xix, fv, temp;
805 t_iatom av, ai, aj, ak;
813 copy_rvec(f[av], fv);
815 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
816 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
819 /* xix goes from i to point x on the line jk */
820 xix[XX] = xij[XX]+a*xjk[XX];
821 xix[YY] = xij[YY]+a*xjk[YY];
822 xix[ZZ] = xij[ZZ]+a*xjk[ZZ];
825 invl = gmx::invsqrt(iprod(xix, xix));
829 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
831 temp[XX] = c*(fv[XX]-fproj*xix[XX]);
832 temp[YY] = c*(fv[YY]-fproj*xix[YY]);
833 temp[ZZ] = c*(fv[ZZ]-fproj*xix[ZZ]);
836 /* c is already calculated in constr_vsite3FD
837 storing c somewhere will save 26 flops! */
840 f[ai][XX] += fv[XX] - temp[XX];
841 f[ai][YY] += fv[YY] - temp[YY];
842 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
843 f[aj][XX] += a1*temp[XX];
844 f[aj][YY] += a1*temp[YY];
845 f[aj][ZZ] += a1*temp[ZZ];
846 f[ak][XX] += a*temp[XX];
847 f[ak][YY] += a*temp[YY];
848 f[ak][ZZ] += a*temp[ZZ];
853 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
855 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
857 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
862 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
869 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
871 rvec_dec(fshift[svi], fv);
872 fshift[CENTRAL][XX] += fv[XX] - (1 + a)*temp[XX];
873 fshift[CENTRAL][YY] += fv[YY] - (1 + a)*temp[YY];
874 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a)*temp[ZZ];
875 fshift[ sji][XX] += temp[XX];
876 fshift[ sji][YY] += temp[YY];
877 fshift[ sji][ZZ] += temp[ZZ];
878 fshift[ skj][XX] += a*temp[XX];
879 fshift[ skj][YY] += a*temp[YY];
880 fshift[ skj][ZZ] += a*temp[ZZ];
885 /* When VirCorr=TRUE, the virial for the current forces is not
886 * calculated from the redistributed forces. This means that
887 * the effect of non-linear virtual site constructions on the virial
888 * needs to be added separately. This contribution can be calculated
889 * in many ways, but the simplest and cheapest way is to use
890 * the first constructing atom ai as a reference position in space:
891 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
895 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
897 for (int i = 0; i < DIM; i++)
899 for (int j = 0; j < DIM; j++)
901 /* As xix is a linear combination of j and k, use that here */
902 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
907 /* TOTAL: 61 flops */
910 static void spread_vsite3FAD(const t_iatom ia[], real a, real b,
912 rvec f[], rvec fshift[],
913 gmx_bool VirCorr, matrix dxdf,
914 const t_pbc *pbc, const t_graph *g)
916 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
917 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
918 t_iatom av, ai, aj, ak;
919 int svi, sji, skj, d;
926 copy_rvec(f[ia[1]], fv);
928 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
929 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
932 invdij = gmx::invsqrt(iprod(xij, xij));
933 invdij2 = invdij * invdij;
934 c1 = iprod(xij, xjk) * invdij2;
935 xperp[XX] = xjk[XX] - c1*xij[XX];
936 xperp[YY] = xjk[YY] - c1*xij[YY];
937 xperp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
938 /* xperp in plane ijk, perp. to ij */
939 invdp = gmx::invsqrt(iprod(xperp, xperp));
944 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
945 storing them somewhere will save 45 flops! */
947 fproj = iprod(xij, fv)*invdij2;
948 svmul(fproj, xij, Fpij); /* proj. f on xij */
949 svmul(iprod(xperp, fv)*invdp*invdp, xperp, Fppp); /* proj. f on xperp */
950 svmul(b1*fproj, xperp, f3);
953 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
954 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
955 for (d = 0; (d < DIM); d++)
963 f[ai][XX] += fv[XX] - f1[XX] + c1*f2[XX] + f3[XX];
964 f[ai][YY] += fv[YY] - f1[YY] + c1*f2[YY] + f3[YY];
965 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1*f2[ZZ] + f3[ZZ];
966 f[aj][XX] += f1[XX] - c2*f2[XX] - f3[XX];
967 f[aj][YY] += f1[YY] - c2*f2[YY] - f3[YY];
968 f[aj][ZZ] += f1[ZZ] - c2*f2[ZZ] - f3[ZZ];
976 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
978 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
980 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
985 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
992 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
994 rvec_dec(fshift[svi], fv);
995 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1-c1)*f2[XX] + f3[XX];
996 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1-c1)*f2[YY] + f3[YY];
997 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1-c1)*f2[ZZ] + f3[ZZ];
998 fshift[ sji][XX] += f1[XX] - c1 *f2[XX] - f3[XX];
999 fshift[ sji][YY] += f1[YY] - c1 *f2[YY] - f3[YY];
1000 fshift[ sji][ZZ] += f1[ZZ] - c1 *f2[ZZ] - f3[ZZ];
1001 fshift[ skj][XX] += f2[XX];
1002 fshift[ skj][YY] += f2[YY];
1003 fshift[ skj][ZZ] += f2[ZZ];
1011 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1013 for (i = 0; i < DIM; i++)
1015 for (j = 0; j < DIM; j++)
1017 /* Note that xik=xij+xjk, so we have to add xij*f2 */
1020 + xij[i]*(f1[j] + (1 - c2)*f2[j] - f3[j])
1026 /* TOTAL: 113 flops */
1029 static void spread_vsite3OUT(const t_iatom ia[], real a, real b, real c,
1031 rvec f[], rvec fshift[],
1032 gmx_bool VirCorr, matrix dxdf,
1033 const t_pbc *pbc, const t_graph *g)
1035 rvec xvi, xij, xik, fv, fj, fk;
1046 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1047 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1050 copy_rvec(f[av], fv);
1057 fj[XX] = a*fv[XX] - xik[ZZ]*cfy + xik[YY]*cfz;
1058 fj[YY] = xik[ZZ]*cfx + a*fv[YY] - xik[XX]*cfz;
1059 fj[ZZ] = -xik[YY]*cfx + xik[XX]*cfy + a*fv[ZZ];
1061 fk[XX] = b*fv[XX] + xij[ZZ]*cfy - xij[YY]*cfz;
1062 fk[YY] = -xij[ZZ]*cfx + b*fv[YY] + xij[XX]*cfz;
1063 fk[ZZ] = xij[YY]*cfx - xij[XX]*cfy + b*fv[ZZ];
1066 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
1067 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
1068 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1069 rvec_inc(f[aj], fj);
1070 rvec_inc(f[ak], fk);
1075 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1077 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1079 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1084 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1091 if (fshift && (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL))
1093 rvec_dec(fshift[svi], fv);
1094 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
1095 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
1096 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1097 rvec_inc(fshift[sji], fj);
1098 rvec_inc(fshift[ski], fk);
1105 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1107 for (int i = 0; i < DIM; i++)
1109 for (int j = 0; j < DIM; j++)
1111 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j];
1116 /* TOTAL: 54 flops */
1119 static void spread_vsite4FD(const t_iatom ia[], real a, real b, real c,
1121 rvec f[], rvec fshift[],
1122 gmx_bool VirCorr, matrix dxdf,
1123 const t_pbc *pbc, const t_graph *g)
1125 real d, invl, fproj, a1;
1126 rvec xvi, xij, xjk, xjl, xix, fv, temp;
1127 int av, ai, aj, ak, al;
1129 int svi, sji, skj, slj, m;
1137 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1138 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1139 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1142 /* xix goes from i to point x on the plane jkl */
1143 for (m = 0; m < DIM; m++)
1145 xix[m] = xij[m] + a*xjk[m] + b*xjl[m];
1149 invl = gmx::invsqrt(iprod(xix, xix));
1151 /* 4 + ?10? flops */
1153 copy_rvec(f[av], fv);
1155 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
1157 for (m = 0; m < DIM; m++)
1159 temp[m] = d*(fv[m] - fproj*xix[m]);
1163 /* c is already calculated in constr_vsite3FD
1164 storing c somewhere will save 35 flops! */
1167 for (m = 0; m < DIM; m++)
1169 f[ai][m] += fv[m] - temp[m];
1170 f[aj][m] += a1*temp[m];
1171 f[ak][m] += a*temp[m];
1172 f[al][m] += b*temp[m];
1178 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1180 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1182 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
1184 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, aj), di);
1189 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1197 (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL))
1199 rvec_dec(fshift[svi], fv);
1200 for (m = 0; m < DIM; m++)
1202 fshift[CENTRAL][m] += fv[m] - (1 + a + b)*temp[m];
1203 fshift[ sji][m] += temp[m];
1204 fshift[ skj][m] += a*temp[m];
1205 fshift[ slj][m] += b*temp[m];
1214 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1216 for (i = 0; i < DIM; i++)
1218 for (j = 0; j < DIM; j++)
1220 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
1225 /* TOTAL: 77 flops */
1229 static void spread_vsite4FDN(const t_iatom ia[], real a, real b, real c,
1231 rvec f[], rvec fshift[],
1232 gmx_bool VirCorr, matrix dxdf,
1233 const t_pbc *pbc, const t_graph *g)
1235 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1236 rvec fv, fj, fk, fl;
1240 int av, ai, aj, ak, al;
1241 int svi, sij, sik, sil;
1243 /* DEBUG: check atom indices */
1250 copy_rvec(f[av], fv);
1252 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1253 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1254 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1267 rvec_sub(ra, xij, rja);
1268 rvec_sub(rb, xij, rjb);
1269 rvec_sub(rb, ra, rab);
1272 cprod(rja, rjb, rm);
1275 invrm = gmx::invsqrt(norm2(rm));
1276 denom = invrm*invrm;
1279 cfx = c*invrm*fv[XX];
1280 cfy = c*invrm*fv[YY];
1281 cfz = c*invrm*fv[ZZ];
1292 fj[XX] = ( -rm[XX]*rt[XX]) * cfx + ( rab[ZZ]-rm[YY]*rt[XX]) * cfy + (-rab[YY]-rm[ZZ]*rt[XX]) * cfz;
1293 fj[YY] = (-rab[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( rab[XX]-rm[ZZ]*rt[YY]) * cfz;
1294 fj[ZZ] = ( rab[YY]-rm[XX]*rt[ZZ]) * cfx + (-rab[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1305 fk[XX] = ( -rm[XX]*rt[XX]) * cfx + (-a*rjb[ZZ]-rm[YY]*rt[XX]) * cfy + ( a*rjb[YY]-rm[ZZ]*rt[XX]) * cfz;
1306 fk[YY] = ( a*rjb[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + (-a*rjb[XX]-rm[ZZ]*rt[YY]) * cfz;
1307 fk[ZZ] = (-a*rjb[YY]-rm[XX]*rt[ZZ]) * cfx + ( a*rjb[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1318 fl[XX] = ( -rm[XX]*rt[XX]) * cfx + ( b*rja[ZZ]-rm[YY]*rt[XX]) * cfy + (-b*rja[YY]-rm[ZZ]*rt[XX]) * cfz;
1319 fl[YY] = (-b*rja[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( b*rja[XX]-rm[ZZ]*rt[YY]) * cfz;
1320 fl[ZZ] = ( b*rja[YY]-rm[XX]*rt[ZZ]) * cfx + (-b*rja[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1323 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1324 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1325 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1326 rvec_inc(f[aj], fj);
1327 rvec_inc(f[ak], fk);
1328 rvec_inc(f[al], fl);
1333 ivec_sub(SHIFT_IVEC(g, av), SHIFT_IVEC(g, ai), di);
1335 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1337 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1339 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, ai), di);
1344 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1351 if (fshift && (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL))
1353 rvec_dec(fshift[svi], fv);
1354 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1355 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1356 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1357 rvec_inc(fshift[sij], fj);
1358 rvec_inc(fshift[sik], fk);
1359 rvec_inc(fshift[sil], fl);
1367 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1369 for (i = 0; i < DIM; i++)
1371 for (j = 0; j < DIM; j++)
1373 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j] + xil[i]*fl[j];
1378 /* Total: 207 flops (Yuck!) */
1382 static int spread_vsiten(const t_iatom ia[], const t_iparams ip[],
1384 rvec f[], rvec fshift[],
1385 const t_pbc *pbc, const t_graph *g)
1393 n3 = 3*ip[ia[0]].vsiten.n;
1395 copy_rvec(x[av], xv);
1397 for (i = 0; i < n3; i += 3)
1402 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
1407 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1413 a = ip[ia[i]].vsiten.a;
1414 svmul(a, f[av], fi);
1415 rvec_inc(f[ai], fi);
1416 if (fshift && siv != CENTRAL)
1418 rvec_inc(fshift[siv], fi);
1419 rvec_dec(fshift[CENTRAL], fi);
1428 static int vsite_count(const t_ilist *ilist, int ftype)
1430 if (ftype == F_VSITEN)
1432 return ilist[ftype].nr/3;
1436 return ilist[ftype].nr/(1 + interaction_function[ftype].nratoms);
1440 static void spread_vsite_f_thread(const rvec x[],
1441 rvec f[], rvec *fshift,
1442 gmx_bool VirCorr, matrix dxdf,
1443 t_iparams ip[], const t_ilist ilist[],
1444 const t_graph *g, const t_pbc *pbc_null)
1446 const PbcMode pbcMode = getPbcMode(pbc_null);
1447 /* We need another pbc pointer, as with charge groups we switch per vsite */
1448 const t_pbc *pbc_null2 = pbc_null;
1449 gmx::ArrayRef<const int> vsite_pbc;
1451 /* this loop goes backwards to be able to build *
1452 * higher type vsites from lower types */
1453 for (int ftype = c_ftypeVsiteEnd - 1; ftype >= c_ftypeVsiteStart; ftype--)
1455 if (ilist[ftype].nr == 0)
1461 int nra = interaction_function[ftype].nratoms;
1463 int nr = ilist[ftype].nr;
1465 const t_iatom *ia = ilist[ftype].iatoms;
1467 if (pbcMode == PbcMode::all)
1469 pbc_null2 = pbc_null;
1472 for (int i = 0; i < nr; )
1476 /* Constants for constructing */
1478 a1 = ip[tp].vsite.a;
1479 /* Construct the vsite depending on type */
1483 spread_vsite2(ia, a1, x, f, fshift, pbc_null2, g);
1486 b1 = ip[tp].vsite.b;
1487 spread_vsite3(ia, a1, b1, x, f, fshift, pbc_null2, g);
1490 b1 = ip[tp].vsite.b;
1491 spread_vsite3FD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1494 b1 = ip[tp].vsite.b;
1495 spread_vsite3FAD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1498 b1 = ip[tp].vsite.b;
1499 c1 = ip[tp].vsite.c;
1500 spread_vsite3OUT(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1503 b1 = ip[tp].vsite.b;
1504 c1 = ip[tp].vsite.c;
1505 spread_vsite4FD(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1508 b1 = ip[tp].vsite.b;
1509 c1 = ip[tp].vsite.c;
1510 spread_vsite4FDN(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1513 inc = spread_vsiten(ia, ip, x, f, fshift, pbc_null2, g);
1516 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
1517 ftype, __FILE__, __LINE__);
1519 clear_rvec(f[ia[1]]);
1521 /* Increment loop variables */
1529 /*! \brief Clears the task force buffer elements that are written by task idTask */
1530 static void clearTaskForceBufferUsedElements(InterdependentTask *idTask)
1532 int ntask = idTask->spreadTask.size();
1533 for (int ti = 0; ti < ntask; ti++)
1535 const AtomIndex *atomList = &idTask->atomIndex[idTask->spreadTask[ti]];
1536 int natom = atomList->atom.size();
1537 RVec *force = idTask->force.data();
1538 for (int i = 0; i < natom; i++)
1540 clear_rvec(force[atomList->atom[i]]);
1545 void spread_vsite_f(const gmx_vsite_t *vsite,
1546 const rvec * gmx_restrict x,
1547 rvec * gmx_restrict f, rvec * gmx_restrict fshift,
1548 gmx_bool VirCorr, matrix vir,
1549 t_nrnb *nrnb, const t_idef *idef,
1550 int ePBC, gmx_bool bMolPBC, const t_graph *g, const matrix box,
1551 const t_commrec *cr, gmx_wallcycle *wcycle)
1553 wallcycle_start(wcycle, ewcVSITESPREAD);
1554 const bool useDomdec = vsite->useDomdec;
1555 GMX_ASSERT(!useDomdec || (cr != nullptr && DOMAINDECOMP(cr)), "When vsites are set up with domain decomposition, we need a valid commrec");
1557 t_pbc pbc, *pbc_null;
1559 /* We only need to do pbc when we have inter-cg vsites */
1560 if ((useDomdec || bMolPBC) && vsite->numInterUpdategroupVsites)
1562 /* This is wasting some CPU time as we now do this multiple times
1565 pbc_null = set_pbc_dd(&pbc, ePBC, useDomdec ? cr->dd->nc : nullptr, FALSE, box);
1574 dd_clear_f_vsites(cr->dd, f);
1577 if (vsite->nthreads == 1)
1584 spread_vsite_f_thread(x, f, fshift,
1586 idef->iparams, idef->il,
1591 for (int i = 0; i < DIM; i++)
1593 for (int j = 0; j < DIM; j++)
1595 vir[i][j] += -0.5*dxdf[i][j];
1602 /* First spread the vsites that might depend on non-local vsites */
1605 clear_mat(vsite->tData[vsite->nthreads]->dxdf);
1607 spread_vsite_f_thread(x, f, fshift,
1608 VirCorr, vsite->tData[vsite->nthreads]->dxdf,
1610 vsite->tData[vsite->nthreads]->ilist,
1613 #pragma omp parallel num_threads(vsite->nthreads)
1617 int thread = gmx_omp_get_thread_num();
1618 VsiteThread &tData = *vsite->tData[thread];
1621 if (thread == 0 || fshift == nullptr)
1627 fshift_t = tData.fshift;
1629 for (int i = 0; i < SHIFTS; i++)
1631 clear_rvec(fshift_t[i]);
1636 clear_mat(tData.dxdf);
1639 if (tData.useInterdependentTask)
1641 /* Spread the vsites that spread outside our local range.
1642 * This is done using a thread-local force buffer force.
1643 * First we need to copy the input vsite forces to force.
1645 InterdependentTask *idTask = &tData.idTask;
1647 /* Clear the buffer elements set by our task during
1648 * the last call to spread_vsite_f.
1650 clearTaskForceBufferUsedElements(idTask);
1652 int nvsite = idTask->vsite.size();
1653 for (int i = 0; i < nvsite; i++)
1655 copy_rvec(f[idTask->vsite[i]],
1656 idTask->force[idTask->vsite[i]]);
1658 spread_vsite_f_thread(x, as_rvec_array(idTask->force.data()), fshift_t,
1659 VirCorr, tData.dxdf,
1664 /* We need a barrier before reducing forces below
1665 * that have been produced by a different thread above.
1669 /* Loop over all thread task and reduce forces they
1670 * produced on atoms that fall in our range.
1671 * Note that atomic reduction would be a simpler solution,
1672 * but that might not have good support on all platforms.
1674 int ntask = idTask->reduceTask.size();
1675 for (int ti = 0; ti < ntask; ti++)
1677 const InterdependentTask *idt_foreign = &vsite->tData[idTask->reduceTask[ti]]->idTask;
1678 const AtomIndex *atomList = &idt_foreign->atomIndex[thread];
1679 const RVec *f_foreign = idt_foreign->force.data();
1681 int natom = atomList->atom.size();
1682 for (int i = 0; i < natom; i++)
1684 int ind = atomList->atom[i];
1685 rvec_inc(f[ind], f_foreign[ind]);
1686 /* Clearing of f_foreign is done at the next step */
1689 /* Clear the vsite forces, both in f and force */
1690 for (int i = 0; i < nvsite; i++)
1692 int ind = tData.idTask.vsite[i];
1694 clear_rvec(tData.idTask.force[ind]);
1698 /* Spread the vsites that spread locally only */
1699 spread_vsite_f_thread(x, f, fshift_t,
1700 VirCorr, tData.dxdf,
1705 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1708 if (fshift != nullptr)
1710 for (int th = 1; th < vsite->nthreads; th++)
1712 for (int i = 0; i < SHIFTS; i++)
1714 rvec_inc(fshift[i], vsite->tData[th]->fshift[i]);
1721 for (int th = 0; th < vsite->nthreads + 1; th++)
1723 /* MSVC doesn't like matrix references, so we use a pointer */
1724 const matrix *dxdf = &vsite->tData[th]->dxdf;
1726 for (int i = 0; i < DIM; i++)
1728 for (int j = 0; j < DIM; j++)
1730 vir[i][j] += -0.5*(*dxdf)[i][j];
1739 dd_move_f_vsites(cr->dd, f, fshift);
1742 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(idef->il, F_VSITE2));
1743 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(idef->il, F_VSITE3));
1744 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(idef->il, F_VSITE3FD));
1745 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(idef->il, F_VSITE3FAD));
1746 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(idef->il, F_VSITE3OUT));
1747 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(idef->il, F_VSITE4FD));
1748 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(idef->il, F_VSITE4FDN));
1749 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(idef->il, F_VSITEN));
1751 wallcycle_stop(wcycle, ewcVSITESPREAD);
1754 /*! \brief Returns the an array with group indices for each atom
1756 * \param[in] grouping The paritioning of the atom range into atom groups
1758 static std::vector<int> makeAtomToGroupMapping(const gmx::RangePartitioning &grouping)
1760 std::vector<int> atomToGroup(grouping.fullRange().end(), 0);
1762 for (int group = 0; group < grouping.numBlocks(); group++)
1764 auto block = grouping.block(group);
1765 std::fill(atomToGroup.begin() + block.begin(),
1766 atomToGroup.begin() + block.end(),
1773 int countNonlinearVsites(const gmx_mtop_t &mtop)
1775 int numNonlinearVsites = 0;
1776 for (const gmx_molblock_t &molb : mtop.molblock)
1778 const gmx_moltype_t &molt = mtop.moltype[molb.type];
1780 for (const auto &ilist : extractILists(molt.ilist, IF_VSITE))
1782 if (ilist.functionType != F_VSITE2 &&
1783 ilist.functionType != F_VSITE3 &&
1784 ilist.functionType != F_VSITEN)
1786 numNonlinearVsites += molb.nmol*ilist.iatoms.size()/(1 + NRAL(ilist.functionType));
1791 return numNonlinearVsites;
1794 int countInterUpdategroupVsites(const gmx_mtop_t &mtop,
1795 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingPerMoleculetype)
1797 int n_intercg_vsite = 0;
1798 for (const gmx_molblock_t &molb : mtop.molblock)
1800 const gmx_moltype_t &molt = mtop.moltype[molb.type];
1802 std::vector<int> atomToGroup;
1803 if (!updateGroupingPerMoleculetype.empty())
1805 atomToGroup = makeAtomToGroupMapping(updateGroupingPerMoleculetype[molb.type]);
1807 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
1809 const int nral = NRAL(ftype);
1810 const InteractionList &il = molt.ilist[ftype];
1811 for (int i = 0; i < il.size(); i += 1 + nral)
1813 bool isInterGroup = atomToGroup.empty();
1816 const int group = atomToGroup[il.iatoms[1 + i]];
1817 for (int a = 1; a < nral; a++)
1819 if (atomToGroup[il.iatoms[1 + a]] != group)
1821 isInterGroup = true;
1828 n_intercg_vsite += molb.nmol;
1834 return n_intercg_vsite;
1837 std::unique_ptr<gmx_vsite_t>
1838 initVsite(const gmx_mtop_t &mtop,
1839 const t_commrec *cr)
1841 GMX_RELEASE_ASSERT(cr != nullptr, "We need a valid commrec");
1843 std::unique_ptr<gmx_vsite_t> vsite;
1845 /* check if there are vsites */
1847 for (int ftype = 0; ftype < F_NRE; ftype++)
1849 if (interaction_function[ftype].flags & IF_VSITE)
1851 GMX_ASSERT(ftype >= c_ftypeVsiteStart && ftype < c_ftypeVsiteEnd, "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
1853 nvsite += gmx_mtop_ftype_count(&mtop, ftype);
1857 GMX_ASSERT(ftype < c_ftypeVsiteStart || ftype >= c_ftypeVsiteEnd, "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
1866 vsite = std::make_unique<gmx_vsite_t>();
1868 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingPerMoleculetype;
1869 if (DOMAINDECOMP(cr))
1871 updateGroupingPerMoleculetype = getUpdateGroupingPerMoleculetype(*cr->dd);
1873 vsite->numInterUpdategroupVsites = countInterUpdategroupVsites(mtop, updateGroupingPerMoleculetype);
1875 vsite->useDomdec = (DOMAINDECOMP(cr) && cr->dd->nnodes > 1);
1877 vsite->nthreads = gmx_omp_nthreads_get(emntVSITE);
1879 if (vsite->nthreads > 1)
1881 /* We need one extra thread data structure for the overlap vsites */
1882 vsite->tData.resize(vsite->nthreads + 1);
1883 #pragma omp parallel for num_threads(vsite->nthreads) schedule(static)
1884 for (int thread = 0; thread < vsite->nthreads; thread++)
1888 vsite->tData[thread] = std::make_unique<VsiteThread>();
1890 InterdependentTask &idTask = vsite->tData[thread]->idTask;
1892 idTask.atomIndex.resize(vsite->nthreads);
1894 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1896 if (vsite->nthreads > 1)
1898 vsite->tData[vsite->nthreads] = std::make_unique<VsiteThread>();
1905 gmx_vsite_t::gmx_vsite_t()
1909 gmx_vsite_t::~gmx_vsite_t()
1913 static inline void flagAtom(InterdependentTask *idTask, int atom,
1914 int thread, int nthread, int natperthread)
1916 if (!idTask->use[atom])
1918 idTask->use[atom] = true;
1919 thread = atom/natperthread;
1920 /* Assign all non-local atom force writes to thread 0 */
1921 if (thread >= nthread)
1925 idTask->atomIndex[thread].atom.push_back(atom);
1929 /*\brief Here we try to assign all vsites that are in our local range.
1931 * Our task local atom range is tData->rangeStart - tData->rangeEnd.
1932 * Vsites that depend only on local atoms, as indicated by taskIndex[]==thread,
1933 * are assigned to task tData->ilist. Vsites that depend on non-local atoms
1934 * but not on other vsites are assigned to task tData->id_task.ilist.
1935 * taskIndex[] is set for all vsites in our range, either to our local tasks
1936 * or to the single last task as taskIndex[]=2*nthreads.
1938 static void assignVsitesToThread(VsiteThread *tData,
1942 gmx::ArrayRef<int> taskIndex,
1943 const t_ilist *ilist,
1944 const t_iparams *ip,
1945 const unsigned short *ptype)
1947 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
1949 tData->ilist[ftype].nr = 0;
1950 tData->idTask.ilist[ftype].nr = 0;
1952 int nral1 = 1 + NRAL(ftype);
1954 t_iatom *iat = ilist[ftype].iatoms;
1955 for (int i = 0; i < ilist[ftype].nr; )
1957 if (ftype == F_VSITEN)
1959 /* The 3 below is from 1+NRAL(ftype)=3 */
1960 inc = ip[iat[i]].vsiten.n*3;
1963 if (iat[1 + i] < tData->rangeStart ||
1964 iat[1 + i] >= tData->rangeEnd)
1966 /* This vsite belongs to a different thread */
1971 /* We would like to assign this vsite to task thread,
1972 * but it might depend on atoms outside the atom range of thread
1973 * or on another vsite not assigned to task thread.
1976 if (ftype != F_VSITEN)
1978 for (int j = i + 2; j < i + nral1; j++)
1980 /* Do a range check to avoid a harmless race on taskIndex */
1981 if (iat[j] < tData->rangeStart ||
1982 iat[j] >= tData->rangeEnd ||
1983 taskIndex[iat[j]] != thread)
1985 if (!tData->useInterdependentTask ||
1986 ptype[iat[j]] == eptVSite)
1988 /* At least one constructing atom is a vsite
1989 * that is not assigned to the same thread.
1990 * Put this vsite into a separate task.
1996 /* There are constructing atoms outside our range,
1997 * put this vsite into a second task to be executed
1998 * on the same thread. During construction no barrier
1999 * is needed between the two tasks on the same thread.
2000 * During spreading we need to run this task with
2001 * an additional thread-local intermediate force buffer
2002 * (or atomic reduction) and a barrier between the two
2005 task = nthread + thread;
2011 for (int j = i + 2; j < i + inc; j += 3)
2013 /* Do a range check to avoid a harmless race on taskIndex */
2014 if (iat[j] < tData->rangeStart ||
2015 iat[j] >= tData->rangeEnd ||
2016 taskIndex[iat[j]] != thread)
2018 GMX_ASSERT(ptype[iat[j]] != eptVSite, "A vsite to be assigned in assignVsitesToThread has a vsite as a constructing atom that does not belong to our task, such vsites should be assigned to the single 'master' task");
2020 task = nthread + thread;
2025 /* Update this vsite's thread index entry */
2026 taskIndex[iat[1+i]] = task;
2028 if (task == thread || task == nthread + thread)
2030 /* Copy this vsite to the thread data struct of thread */
2034 il_task = &tData->ilist[ftype];
2038 il_task = &tData->idTask.ilist[ftype];
2040 /* Ensure we have sufficient memory allocated */
2041 if (il_task->nr + inc > il_task->nalloc)
2043 il_task->nalloc = over_alloc_large(il_task->nr + inc);
2044 srenew(il_task->iatoms, il_task->nalloc);
2046 /* Copy the vsite data to the thread-task local array */
2047 for (int j = i; j < i + inc; j++)
2049 il_task->iatoms[il_task->nr++] = iat[j];
2051 if (task == nthread + thread)
2053 /* This vsite write outside our own task force block.
2054 * Put it into the interdependent task list and flag
2055 * the atoms involved for reduction.
2057 tData->idTask.vsite.push_back(iat[i + 1]);
2058 if (ftype != F_VSITEN)
2060 for (int j = i + 2; j < i + nral1; j++)
2062 flagAtom(&tData->idTask, iat[j],
2063 thread, nthread, natperthread);
2068 for (int j = i + 2; j < i + inc; j += 3)
2070 flagAtom(&tData->idTask, iat[j],
2071 thread, nthread, natperthread);
2082 /*! \brief Assign all vsites with taskIndex[]==task to task tData */
2083 static void assignVsitesToSingleTask(VsiteThread *tData,
2085 gmx::ArrayRef<const int> taskIndex,
2086 const t_ilist *ilist,
2087 const t_iparams *ip)
2089 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2091 tData->ilist[ftype].nr = 0;
2092 tData->idTask.ilist[ftype].nr = 0;
2094 int nral1 = 1 + NRAL(ftype);
2096 t_iatom *iat = ilist[ftype].iatoms;
2097 t_ilist *il_task = &tData->ilist[ftype];
2099 for (int i = 0; i < ilist[ftype].nr; )
2101 if (ftype == F_VSITEN)
2103 /* The 3 below is from 1+NRAL(ftype)=3 */
2104 inc = ip[iat[i]].vsiten.n*3;
2106 /* Check if the vsite is assigned to our task */
2107 if (taskIndex[iat[1 + i]] == task)
2109 /* Ensure we have sufficient memory allocated */
2110 if (il_task->nr + inc > il_task->nalloc)
2112 il_task->nalloc = over_alloc_large(il_task->nr + inc);
2113 srenew(il_task->iatoms, il_task->nalloc);
2115 /* Copy the vsite data to the thread-task local array */
2116 for (int j = i; j < i + inc; j++)
2118 il_task->iatoms[il_task->nr++] = iat[j];
2127 void split_vsites_over_threads(const t_ilist *ilist,
2128 const t_iparams *ip,
2129 const t_mdatoms *mdatoms,
2132 int vsite_atom_range, natperthread;
2134 if (vsite->nthreads == 1)
2140 /* The current way of distributing the vsites over threads in primitive.
2141 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
2142 * without taking into account how the vsites are distributed.
2143 * Without domain decomposition we at least tighten the upper bound
2144 * of the range (useful for common systems such as a vsite-protein
2146 * With domain decomposition, as long as the vsites are distributed
2147 * uniformly in each domain along the major dimension, usually x,
2148 * it will also perform well.
2150 if (!vsite->useDomdec)
2152 vsite_atom_range = -1;
2153 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2156 if (ftype != F_VSITEN)
2158 int nral1 = 1 + NRAL(ftype);
2159 const t_iatom *iat = ilist[ftype].iatoms;
2160 for (int i = 0; i < ilist[ftype].nr; i += nral1)
2162 for (int j = i + 1; j < i + nral1; j++)
2164 vsite_atom_range = std::max(vsite_atom_range, iat[j]);
2172 const t_iatom *iat = ilist[ftype].iatoms;
2175 while (i < ilist[ftype].nr)
2177 /* The 3 below is from 1+NRAL(ftype)=3 */
2178 vs_ind_end = i + ip[iat[i]].vsiten.n*3;
2180 vsite_atom_range = std::max(vsite_atom_range, iat[i+1]);
2181 while (i < vs_ind_end)
2183 vsite_atom_range = std::max(vsite_atom_range, iat[i+2]);
2191 natperthread = (vsite_atom_range + vsite->nthreads - 1)/vsite->nthreads;
2195 /* Any local or not local atom could be involved in virtual sites.
2196 * But since we usually have very few non-local virtual sites
2197 * (only non-local vsites that depend on local vsites),
2198 * we distribute the local atom range equally over the threads.
2199 * When assigning vsites to threads, we should take care that the last
2200 * threads also covers the non-local range.
2202 vsite_atom_range = mdatoms->nr;
2203 natperthread = (mdatoms->homenr + vsite->nthreads - 1)/vsite->nthreads;
2208 fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n", mdatoms->nr, vsite_atom_range, natperthread);
2211 /* To simplify the vsite assignment, we make an index which tells us
2212 * to which task particles, both non-vsites and vsites, are assigned.
2214 vsite->taskIndex.resize(mdatoms->nr);
2216 /* Initialize the task index array. Here we assign the non-vsite
2217 * particles to task=thread, so we easily figure out if vsites
2218 * depend on local and/or non-local particles in assignVsitesToThread.
2220 gmx::ArrayRef<int> taskIndex = vsite->taskIndex;
2223 for (int i = 0; i < mdatoms->nr; i++)
2225 if (mdatoms->ptype[i] == eptVSite)
2227 /* vsites are not assigned to a task yet */
2232 /* assign non-vsite particles to task thread */
2233 taskIndex[i] = thread;
2235 if (i == (thread + 1)*natperthread && thread < vsite->nthreads)
2242 #pragma omp parallel num_threads(vsite->nthreads)
2246 int thread = gmx_omp_get_thread_num();
2247 VsiteThread &tData = *vsite->tData[thread];
2249 /* Clear the buffer use flags that were set before */
2250 if (tData.useInterdependentTask)
2252 InterdependentTask &idTask = tData.idTask;
2254 /* To avoid an extra OpenMP barrier in spread_vsite_f,
2255 * we clear the force buffer at the next step,
2256 * so we need to do it here as well.
2258 clearTaskForceBufferUsedElements(&idTask);
2260 idTask.vsite.resize(0);
2261 for (int t = 0; t < vsite->nthreads; t++)
2263 AtomIndex &atomIndex = idTask.atomIndex[t];
2264 int natom = atomIndex.atom.size();
2265 for (int i = 0; i < natom; i++)
2267 idTask.use[atomIndex.atom[i]] = false;
2269 atomIndex.atom.resize(0);
2274 /* To avoid large f_buf allocations of #threads*vsite_atom_range
2275 * we don't use task2 with more than 200000 atoms. This doesn't
2276 * affect performance, since with such a large range relatively few
2277 * vsites will end up in the separate task.
2278 * Note that useTask2 should be the same for all threads.
2280 tData.useInterdependentTask = (vsite_atom_range <= 200000);
2281 if (tData.useInterdependentTask)
2283 size_t natoms_use_in_vsites = vsite_atom_range;
2284 InterdependentTask &idTask = tData.idTask;
2285 /* To avoid resizing and re-clearing every nstlist steps,
2286 * we never down size the force buffer.
2288 if (natoms_use_in_vsites > idTask.force.size() ||
2289 natoms_use_in_vsites > idTask.use.size())
2291 idTask.force.resize(natoms_use_in_vsites, { 0, 0, 0 });
2292 idTask.use.resize(natoms_use_in_vsites, false);
2296 /* Assign all vsites that can execute independently on threads */
2297 tData.rangeStart = thread *natperthread;
2298 if (thread < vsite->nthreads - 1)
2300 tData.rangeEnd = (thread + 1)*natperthread;
2304 /* The last thread should cover up to the end of the range */
2305 tData.rangeEnd = mdatoms->nr;
2307 assignVsitesToThread(&tData,
2308 thread, vsite->nthreads,
2311 ilist, ip, mdatoms->ptype);
2313 if (tData.useInterdependentTask)
2315 /* In the worst case, all tasks write to force ranges of
2316 * all other tasks, leading to #tasks^2 scaling (this is only
2317 * the overhead, the actual flops remain constant).
2318 * But in most cases there is far less coupling. To improve
2319 * scaling at high thread counts we therefore construct
2320 * an index to only loop over the actually affected tasks.
2322 InterdependentTask &idTask = tData.idTask;
2324 /* Ensure assignVsitesToThread finished on other threads */
2327 idTask.spreadTask.resize(0);
2328 idTask.reduceTask.resize(0);
2329 for (int t = 0; t < vsite->nthreads; t++)
2331 /* Do we write to the force buffer of task t? */
2332 if (!idTask.atomIndex[t].atom.empty())
2334 idTask.spreadTask.push_back(t);
2336 /* Does task t write to our force buffer? */
2337 if (!vsite->tData[t]->idTask.atomIndex[thread].atom.empty())
2339 idTask.reduceTask.push_back(t);
2344 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2346 /* Assign all remaining vsites, that will have taskIndex[]=2*vsite->nthreads,
2347 * to a single task that will not run in parallel with other tasks.
2349 assignVsitesToSingleTask(vsite->tData[vsite->nthreads].get(),
2354 if (debug && vsite->nthreads > 1)
2356 fprintf(debug, "virtual site useInterdependentTask %d, nuse:\n",
2357 static_cast<int>(vsite->tData[0]->useInterdependentTask));
2358 for (int th = 0; th < vsite->nthreads + 1; th++)
2360 fprintf(debug, " %4d", vsite->tData[th]->idTask.nuse);
2362 fprintf(debug, "\n");
2364 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2366 if (ilist[ftype].nr > 0)
2368 fprintf(debug, "%-20s thread dist:",
2369 interaction_function[ftype].longname);
2370 for (int th = 0; th < vsite->nthreads + 1; th++)
2372 fprintf(debug, " %4d %4d ",
2373 vsite->tData[th]->ilist[ftype].nr,
2374 vsite->tData[th]->idTask.ilist[ftype].nr);
2376 fprintf(debug, "\n");
2382 int nrOrig = vsiteIlistNrCount(ilist);
2384 for (int th = 0; th < vsite->nthreads + 1; th++)
2387 vsiteIlistNrCount(vsite->tData[th]->ilist) +
2388 vsiteIlistNrCount(vsite->tData[th]->idTask.ilist);
2390 GMX_ASSERT(nrThreaded == nrOrig, "The number of virtual sites assigned to all thread task has to match the total number of virtual sites");
2394 void set_vsite_top(gmx_vsite_t *vsite,
2395 const gmx_localtop_t *top,
2396 const t_mdatoms *md)
2398 if (vsite->nthreads > 1)
2400 split_vsites_over_threads(top->idef.il, top->idef.iparams,