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, 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.
46 #include "gromacs/compat/make_unique.h"
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 chargeGroup, // Keep vsite in the same periodic image as the rest of it's charge group
427 none // No PBC treatment needed
430 /*! \brief Returns the PBC mode based on the system PBC and vsite properties
432 * \param[in] pbcPtr A pointer to a PBC struct or nullptr when no PBC treatment is required
433 * \param[in] vsite A pointer to the vsite struct, can be nullptr
435 static PbcMode getPbcMode(const t_pbc *pbcPtr,
436 const gmx_vsite_t *vsite)
438 if (pbcPtr == nullptr)
440 return PbcMode::none;
442 else if (vsite != nullptr && vsite->bHaveChargeGroups)
444 return PbcMode::chargeGroup;
452 static void construct_vsites_thread(const gmx_vsite_t *vsite,
455 const t_iparams ip[], const t_ilist ilist[],
456 const t_pbc *pbc_null)
468 const PbcMode pbcMode = getPbcMode(pbc_null, vsite);
469 /* We need another pbc pointer, as with charge groups we switch per vsite */
470 const t_pbc *pbc_null2 = pbc_null;
471 gmx::ArrayRef<const int> vsite_pbc;
473 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
475 if (ilist[ftype].nr == 0)
481 int nra = interaction_function[ftype].nratoms;
483 int nr = ilist[ftype].nr;
485 const t_iatom *ia = ilist[ftype].iatoms;
487 if (pbcMode == PbcMode::chargeGroup)
489 vsite_pbc = (*vsite->vsite_pbc_loc)[ftype - c_ftypeVsiteStart];
492 for (int i = 0; i < nr; )
495 /* The vsite and constructing atoms */
498 /* Constants for constructing vsites */
499 real a1 = ip[tp].vsite.a;
500 /* Check what kind of pbc we need to use */
503 if (pbcMode == PbcMode::all)
505 /* No charge groups, vsite follows its own pbc */
507 copy_rvec(x[avsite], xpbc);
509 else if (pbcMode == PbcMode::chargeGroup)
511 pbc_atom = vsite_pbc[i/(1 + nra)];
516 /* We need to copy the coordinates here,
517 * single for single atom cg's pbc_atom
518 * is the vsite itself.
520 copy_rvec(x[pbc_atom], xpbc);
522 pbc_null2 = pbc_null;
533 /* Copy the old position */
535 copy_rvec(x[avsite], xv);
537 /* Construct the vsite depending on type */
544 constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
550 constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
556 constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
562 constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
569 constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
577 constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
586 constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
590 inc = constr_vsiten(ia, ip, x, pbc_null2);
593 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
594 ftype, __FILE__, __LINE__);
599 /* Match the pbc of this vsite to the rest of its charge group */
601 int ishift = pbc_dx_aiuc(pbc_null, x[avsite], xpbc, dx);
602 if (ishift != CENTRAL)
604 rvec_add(xpbc, dx, x[avsite]);
609 /* Calculate velocity of vsite... */
611 rvec_sub(x[avsite], xv, vv);
612 svmul(inv_dt, vv, v[avsite]);
615 /* Increment loop variables */
623 void construct_vsites(const gmx_vsite_t *vsite,
626 const t_iparams ip[], const t_ilist ilist[],
627 int ePBC, gmx_bool bMolPBC,
631 const bool useDomdec = (vsite != nullptr && vsite->useDomdec);
632 GMX_ASSERT(!useDomdec || (cr != nullptr && DOMAINDECOMP(cr)), "When vsites are set up with domain decomposition, we need a valid commrec");
633 // TODO: Remove this assertion when we remove charge groups
634 GMX_ASSERT(vsite != nullptr || ePBC == epbcNONE, "Without a vsite struct we can not do PBC (in case we have charge groups)");
636 t_pbc pbc, *pbc_null;
638 /* We only need to do pbc when we have inter-cg vsites.
639 * Note that with domain decomposition we do not need to apply PBC here
640 * when we have at least 3 domains along each dimension. Currently we
641 * do not optimize this case.
643 if (ePBC != epbcNONE && (useDomdec || bMolPBC) &&
644 !(vsite != nullptr && vsite->n_intercg_vsite == 0))
646 /* This is wasting some CPU time as we now do this multiple times
650 clear_ivec(null_ivec);
651 pbc_null = set_pbc_dd(&pbc, ePBC,
652 useDomdec ? cr->dd->nc : null_ivec,
662 dd_move_x_vsites(cr->dd, box, x);
665 if (vsite == nullptr || vsite->nthreads == 1)
667 construct_vsites_thread(vsite,
674 #pragma omp parallel num_threads(vsite->nthreads)
678 const int th = gmx_omp_get_thread_num();
679 const VsiteThread &tData = *vsite->tData[th];
680 GMX_ASSERT(tData.rangeStart >= 0, "The thread data should be initialized before calling construct_vsites");
682 construct_vsites_thread(vsite,
686 if (tData.useInterdependentTask)
688 /* Here we don't need a barrier (unlike the spreading),
689 * since both tasks only construct vsites from particles,
690 * or local vsites, not from non-local vsites.
692 construct_vsites_thread(vsite,
694 ip, tData.idTask.ilist,
698 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
700 /* Now we can construct the vsites that might depend on other vsites */
701 construct_vsites_thread(vsite,
703 ip, vsite->tData[vsite->nthreads]->ilist,
708 static void spread_vsite2(const t_iatom ia[], real a,
710 rvec f[], rvec fshift[],
711 const t_pbc *pbc, const t_graph *g)
722 svmul(1 - a, f[av], fi);
723 svmul( a, f[av], fj);
732 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
734 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
739 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
740 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
748 if (fshift && (siv != CENTRAL || sij != CENTRAL))
750 rvec_inc(fshift[siv], f[av]);
751 rvec_dec(fshift[CENTRAL], fi);
752 rvec_dec(fshift[sij], fj);
755 /* TOTAL: 13 flops */
758 void constructVsitesGlobal(const gmx_mtop_t &mtop,
759 gmx::ArrayRef<gmx::RVec> x)
761 GMX_ASSERT(x.size() >= static_cast<gmx::index>(mtop.natoms), "x should contain the whole system");
762 GMX_ASSERT(!mtop.moleculeBlockIndices.empty(), "molblock indices are needed in constructVsitesGlobal");
764 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
766 const gmx_molblock_t &molb = mtop.molblock[mb];
767 const gmx_moltype_t &molt = mtop.moltype[molb.type];
768 if (vsiteIlistNrCount(molt.ilist) > 0)
770 int atomOffset = mtop.moleculeBlockIndices[mb].globalAtomStart;
771 for (int mol = 0; mol < molb.nmol; mol++)
773 t_ilist ilist[F_NRE];
774 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
776 ilist[ftype].nr = molt.ilist[ftype].size();
777 ilist[ftype].iatoms = const_cast<t_iatom *>(molt.ilist[ftype].iatoms.data());
780 construct_vsites(nullptr, as_rvec_array(x.data()) + atomOffset,
782 mtop.ffparams.iparams, ilist,
783 epbcNONE, TRUE, nullptr, nullptr);
784 atomOffset += molt.atoms.nr;
790 static void spread_vsite3(const t_iatom ia[], real a, real b,
792 rvec f[], rvec fshift[],
793 const t_pbc *pbc, const t_graph *g)
805 svmul(1 - a - b, f[av], fi);
806 svmul( a, f[av], fj);
807 svmul( b, f[av], fk);
817 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ia[1]), di);
819 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
821 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), di);
826 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
827 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
828 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
837 if (fshift && (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL))
839 rvec_inc(fshift[siv], f[av]);
840 rvec_dec(fshift[CENTRAL], fi);
841 rvec_dec(fshift[sij], fj);
842 rvec_dec(fshift[sik], fk);
845 /* TOTAL: 20 flops */
848 static void spread_vsite3FD(const t_iatom ia[], real a, real b,
850 rvec f[], rvec fshift[],
851 gmx_bool VirCorr, matrix dxdf,
852 const t_pbc *pbc, const t_graph *g)
854 real c, invl, fproj, a1;
855 rvec xvi, xij, xjk, xix, fv, temp;
856 t_iatom av, ai, aj, ak;
864 copy_rvec(f[av], fv);
866 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
867 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
870 /* xix goes from i to point x on the line jk */
871 xix[XX] = xij[XX]+a*xjk[XX];
872 xix[YY] = xij[YY]+a*xjk[YY];
873 xix[ZZ] = xij[ZZ]+a*xjk[ZZ];
876 invl = gmx::invsqrt(iprod(xix, xix));
880 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
882 temp[XX] = c*(fv[XX]-fproj*xix[XX]);
883 temp[YY] = c*(fv[YY]-fproj*xix[YY]);
884 temp[ZZ] = c*(fv[ZZ]-fproj*xix[ZZ]);
887 /* c is already calculated in constr_vsite3FD
888 storing c somewhere will save 26 flops! */
891 f[ai][XX] += fv[XX] - temp[XX];
892 f[ai][YY] += fv[YY] - temp[YY];
893 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
894 f[aj][XX] += a1*temp[XX];
895 f[aj][YY] += a1*temp[YY];
896 f[aj][ZZ] += a1*temp[ZZ];
897 f[ak][XX] += a*temp[XX];
898 f[ak][YY] += a*temp[YY];
899 f[ak][ZZ] += a*temp[ZZ];
904 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
906 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
908 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
913 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
920 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
922 rvec_dec(fshift[svi], fv);
923 fshift[CENTRAL][XX] += fv[XX] - (1 + a)*temp[XX];
924 fshift[CENTRAL][YY] += fv[YY] - (1 + a)*temp[YY];
925 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a)*temp[ZZ];
926 fshift[ sji][XX] += temp[XX];
927 fshift[ sji][YY] += temp[YY];
928 fshift[ sji][ZZ] += temp[ZZ];
929 fshift[ skj][XX] += a*temp[XX];
930 fshift[ skj][YY] += a*temp[YY];
931 fshift[ skj][ZZ] += a*temp[ZZ];
936 /* When VirCorr=TRUE, the virial for the current forces is not
937 * calculated from the redistributed forces. This means that
938 * the effect of non-linear virtual site constructions on the virial
939 * needs to be added separately. This contribution can be calculated
940 * in many ways, but the simplest and cheapest way is to use
941 * the first constructing atom ai as a reference position in space:
942 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
946 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
948 for (int i = 0; i < DIM; i++)
950 for (int j = 0; j < DIM; j++)
952 /* As xix is a linear combination of j and k, use that here */
953 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
958 /* TOTAL: 61 flops */
961 static void spread_vsite3FAD(const t_iatom ia[], real a, real b,
963 rvec f[], rvec fshift[],
964 gmx_bool VirCorr, matrix dxdf,
965 const t_pbc *pbc, const t_graph *g)
967 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
968 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
969 t_iatom av, ai, aj, ak;
970 int svi, sji, skj, d;
977 copy_rvec(f[ia[1]], fv);
979 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
980 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
983 invdij = gmx::invsqrt(iprod(xij, xij));
984 invdij2 = invdij * invdij;
985 c1 = iprod(xij, xjk) * invdij2;
986 xperp[XX] = xjk[XX] - c1*xij[XX];
987 xperp[YY] = xjk[YY] - c1*xij[YY];
988 xperp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
989 /* xperp in plane ijk, perp. to ij */
990 invdp = gmx::invsqrt(iprod(xperp, xperp));
995 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
996 storing them somewhere will save 45 flops! */
998 fproj = iprod(xij, fv)*invdij2;
999 svmul(fproj, xij, Fpij); /* proj. f on xij */
1000 svmul(iprod(xperp, fv)*invdp*invdp, xperp, Fppp); /* proj. f on xperp */
1001 svmul(b1*fproj, xperp, f3);
1004 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
1005 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
1006 for (d = 0; (d < DIM); d++)
1014 f[ai][XX] += fv[XX] - f1[XX] + c1*f2[XX] + f3[XX];
1015 f[ai][YY] += fv[YY] - f1[YY] + c1*f2[YY] + f3[YY];
1016 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1*f2[ZZ] + f3[ZZ];
1017 f[aj][XX] += f1[XX] - c2*f2[XX] - f3[XX];
1018 f[aj][YY] += f1[YY] - c2*f2[YY] - f3[YY];
1019 f[aj][ZZ] += f1[ZZ] - c2*f2[ZZ] - f3[ZZ];
1020 f[ak][XX] += f2[XX];
1021 f[ak][YY] += f2[YY];
1022 f[ak][ZZ] += f2[ZZ];
1027 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1029 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1031 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
1036 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1043 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
1045 rvec_dec(fshift[svi], fv);
1046 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1-c1)*f2[XX] + f3[XX];
1047 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1-c1)*f2[YY] + f3[YY];
1048 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1-c1)*f2[ZZ] + f3[ZZ];
1049 fshift[ sji][XX] += f1[XX] - c1 *f2[XX] - f3[XX];
1050 fshift[ sji][YY] += f1[YY] - c1 *f2[YY] - f3[YY];
1051 fshift[ sji][ZZ] += f1[ZZ] - c1 *f2[ZZ] - f3[ZZ];
1052 fshift[ skj][XX] += f2[XX];
1053 fshift[ skj][YY] += f2[YY];
1054 fshift[ skj][ZZ] += f2[ZZ];
1062 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1064 for (i = 0; i < DIM; i++)
1066 for (j = 0; j < DIM; j++)
1068 /* Note that xik=xij+xjk, so we have to add xij*f2 */
1071 + xij[i]*(f1[j] + (1 - c2)*f2[j] - f3[j])
1077 /* TOTAL: 113 flops */
1080 static void spread_vsite3OUT(const t_iatom ia[], real a, real b, real c,
1082 rvec f[], rvec fshift[],
1083 gmx_bool VirCorr, matrix dxdf,
1084 const t_pbc *pbc, const t_graph *g)
1086 rvec xvi, xij, xik, fv, fj, fk;
1097 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1098 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1101 copy_rvec(f[av], fv);
1108 fj[XX] = a*fv[XX] - xik[ZZ]*cfy + xik[YY]*cfz;
1109 fj[YY] = xik[ZZ]*cfx + a*fv[YY] - xik[XX]*cfz;
1110 fj[ZZ] = -xik[YY]*cfx + xik[XX]*cfy + a*fv[ZZ];
1112 fk[XX] = b*fv[XX] + xij[ZZ]*cfy - xij[YY]*cfz;
1113 fk[YY] = -xij[ZZ]*cfx + b*fv[YY] + xij[XX]*cfz;
1114 fk[ZZ] = xij[YY]*cfx - xij[XX]*cfy + b*fv[ZZ];
1117 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
1118 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
1119 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1120 rvec_inc(f[aj], fj);
1121 rvec_inc(f[ak], fk);
1126 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1128 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1130 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1135 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1142 if (fshift && (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL))
1144 rvec_dec(fshift[svi], fv);
1145 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
1146 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
1147 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1148 rvec_inc(fshift[sji], fj);
1149 rvec_inc(fshift[ski], fk);
1156 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1158 for (int i = 0; i < DIM; i++)
1160 for (int j = 0; j < DIM; j++)
1162 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j];
1167 /* TOTAL: 54 flops */
1170 static void spread_vsite4FD(const t_iatom ia[], real a, real b, real c,
1172 rvec f[], rvec fshift[],
1173 gmx_bool VirCorr, matrix dxdf,
1174 const t_pbc *pbc, const t_graph *g)
1176 real d, invl, fproj, a1;
1177 rvec xvi, xij, xjk, xjl, xix, fv, temp;
1178 int av, ai, aj, ak, al;
1180 int svi, sji, skj, slj, m;
1188 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1189 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1190 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1193 /* xix goes from i to point x on the plane jkl */
1194 for (m = 0; m < DIM; m++)
1196 xix[m] = xij[m] + a*xjk[m] + b*xjl[m];
1200 invl = gmx::invsqrt(iprod(xix, xix));
1202 /* 4 + ?10? flops */
1204 copy_rvec(f[av], fv);
1206 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
1208 for (m = 0; m < DIM; m++)
1210 temp[m] = d*(fv[m] - fproj*xix[m]);
1214 /* c is already calculated in constr_vsite3FD
1215 storing c somewhere will save 35 flops! */
1218 for (m = 0; m < DIM; m++)
1220 f[ai][m] += fv[m] - temp[m];
1221 f[aj][m] += a1*temp[m];
1222 f[ak][m] += a*temp[m];
1223 f[al][m] += b*temp[m];
1229 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1231 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1233 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
1235 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, aj), di);
1240 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1248 (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL))
1250 rvec_dec(fshift[svi], fv);
1251 for (m = 0; m < DIM; m++)
1253 fshift[CENTRAL][m] += fv[m] - (1 + a + b)*temp[m];
1254 fshift[ sji][m] += temp[m];
1255 fshift[ skj][m] += a*temp[m];
1256 fshift[ slj][m] += b*temp[m];
1265 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1267 for (i = 0; i < DIM; i++)
1269 for (j = 0; j < DIM; j++)
1271 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
1276 /* TOTAL: 77 flops */
1280 static void spread_vsite4FDN(const t_iatom ia[], real a, real b, real c,
1282 rvec f[], rvec fshift[],
1283 gmx_bool VirCorr, matrix dxdf,
1284 const t_pbc *pbc, const t_graph *g)
1286 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1287 rvec fv, fj, fk, fl;
1291 int av, ai, aj, ak, al;
1292 int svi, sij, sik, sil;
1294 /* DEBUG: check atom indices */
1301 copy_rvec(f[av], fv);
1303 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1304 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1305 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1318 rvec_sub(ra, xij, rja);
1319 rvec_sub(rb, xij, rjb);
1320 rvec_sub(rb, ra, rab);
1323 cprod(rja, rjb, rm);
1326 invrm = gmx::invsqrt(norm2(rm));
1327 denom = invrm*invrm;
1330 cfx = c*invrm*fv[XX];
1331 cfy = c*invrm*fv[YY];
1332 cfz = c*invrm*fv[ZZ];
1343 fj[XX] = ( -rm[XX]*rt[XX]) * cfx + ( rab[ZZ]-rm[YY]*rt[XX]) * cfy + (-rab[YY]-rm[ZZ]*rt[XX]) * cfz;
1344 fj[YY] = (-rab[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( rab[XX]-rm[ZZ]*rt[YY]) * cfz;
1345 fj[ZZ] = ( rab[YY]-rm[XX]*rt[ZZ]) * cfx + (-rab[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1356 fk[XX] = ( -rm[XX]*rt[XX]) * cfx + (-a*rjb[ZZ]-rm[YY]*rt[XX]) * cfy + ( a*rjb[YY]-rm[ZZ]*rt[XX]) * cfz;
1357 fk[YY] = ( a*rjb[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + (-a*rjb[XX]-rm[ZZ]*rt[YY]) * cfz;
1358 fk[ZZ] = (-a*rjb[YY]-rm[XX]*rt[ZZ]) * cfx + ( a*rjb[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1369 fl[XX] = ( -rm[XX]*rt[XX]) * cfx + ( b*rja[ZZ]-rm[YY]*rt[XX]) * cfy + (-b*rja[YY]-rm[ZZ]*rt[XX]) * cfz;
1370 fl[YY] = (-b*rja[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( b*rja[XX]-rm[ZZ]*rt[YY]) * cfz;
1371 fl[ZZ] = ( b*rja[YY]-rm[XX]*rt[ZZ]) * cfx + (-b*rja[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1374 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1375 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1376 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1377 rvec_inc(f[aj], fj);
1378 rvec_inc(f[ak], fk);
1379 rvec_inc(f[al], fl);
1384 ivec_sub(SHIFT_IVEC(g, av), SHIFT_IVEC(g, ai), di);
1386 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1388 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1390 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, ai), di);
1395 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1402 if (fshift && (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL))
1404 rvec_dec(fshift[svi], fv);
1405 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1406 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1407 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1408 rvec_inc(fshift[sij], fj);
1409 rvec_inc(fshift[sik], fk);
1410 rvec_inc(fshift[sil], fl);
1418 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1420 for (i = 0; i < DIM; i++)
1422 for (j = 0; j < DIM; j++)
1424 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j] + xil[i]*fl[j];
1429 /* Total: 207 flops (Yuck!) */
1433 static int spread_vsiten(const t_iatom ia[], const t_iparams ip[],
1435 rvec f[], rvec fshift[],
1436 const t_pbc *pbc, const t_graph *g)
1444 n3 = 3*ip[ia[0]].vsiten.n;
1446 copy_rvec(x[av], xv);
1448 for (i = 0; i < n3; i += 3)
1453 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
1458 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1464 a = ip[ia[i]].vsiten.a;
1465 svmul(a, f[av], fi);
1466 rvec_inc(f[ai], fi);
1467 if (fshift && siv != CENTRAL)
1469 rvec_inc(fshift[siv], fi);
1470 rvec_dec(fshift[CENTRAL], fi);
1479 static int vsite_count(const t_ilist *ilist, int ftype)
1481 if (ftype == F_VSITEN)
1483 return ilist[ftype].nr/3;
1487 return ilist[ftype].nr/(1 + interaction_function[ftype].nratoms);
1491 static void spread_vsite_f_thread(const gmx_vsite_t *vsite,
1493 rvec f[], rvec *fshift,
1494 gmx_bool VirCorr, matrix dxdf,
1495 t_iparams ip[], const t_ilist ilist[],
1496 const t_graph *g, const t_pbc *pbc_null)
1498 const PbcMode pbcMode = getPbcMode(pbc_null, vsite);
1499 /* We need another pbc pointer, as with charge groups we switch per vsite */
1500 const t_pbc *pbc_null2 = pbc_null;
1501 gmx::ArrayRef<const int> vsite_pbc;
1503 /* this loop goes backwards to be able to build *
1504 * higher type vsites from lower types */
1505 for (int ftype = c_ftypeVsiteEnd - 1; ftype >= c_ftypeVsiteStart; ftype--)
1507 if (ilist[ftype].nr == 0)
1513 int nra = interaction_function[ftype].nratoms;
1515 int nr = ilist[ftype].nr;
1517 const t_iatom *ia = ilist[ftype].iatoms;
1519 if (pbcMode == PbcMode::all)
1521 pbc_null2 = pbc_null;
1523 else if (pbcMode == PbcMode::chargeGroup)
1525 if (vsite->vsite_pbc_loc)
1527 vsite_pbc = (*vsite->vsite_pbc_loc)[ftype - c_ftypeVsiteStart];
1531 for (int i = 0; i < nr; )
1533 if (!vsite_pbc.empty())
1535 if (vsite_pbc[i/(1 + nra)] > -2)
1537 pbc_null2 = pbc_null;
1541 pbc_null2 = nullptr;
1547 /* Constants for constructing */
1549 a1 = ip[tp].vsite.a;
1550 /* Construct the vsite depending on type */
1554 spread_vsite2(ia, a1, x, f, fshift, pbc_null2, g);
1557 b1 = ip[tp].vsite.b;
1558 spread_vsite3(ia, a1, b1, x, f, fshift, pbc_null2, g);
1561 b1 = ip[tp].vsite.b;
1562 spread_vsite3FD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1565 b1 = ip[tp].vsite.b;
1566 spread_vsite3FAD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1569 b1 = ip[tp].vsite.b;
1570 c1 = ip[tp].vsite.c;
1571 spread_vsite3OUT(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1574 b1 = ip[tp].vsite.b;
1575 c1 = ip[tp].vsite.c;
1576 spread_vsite4FD(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1579 b1 = ip[tp].vsite.b;
1580 c1 = ip[tp].vsite.c;
1581 spread_vsite4FDN(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1584 inc = spread_vsiten(ia, ip, x, f, fshift, pbc_null2, g);
1587 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
1588 ftype, __FILE__, __LINE__);
1590 clear_rvec(f[ia[1]]);
1592 /* Increment loop variables */
1600 /*! \brief Clears the task force buffer elements that are written by task idTask */
1601 static void clearTaskForceBufferUsedElements(InterdependentTask *idTask)
1603 int ntask = idTask->spreadTask.size();
1604 for (int ti = 0; ti < ntask; ti++)
1606 const AtomIndex *atomList = &idTask->atomIndex[idTask->spreadTask[ti]];
1607 int natom = atomList->atom.size();
1608 RVec *force = idTask->force.data();
1609 for (int i = 0; i < natom; i++)
1611 clear_rvec(force[atomList->atom[i]]);
1616 void spread_vsite_f(const gmx_vsite_t *vsite,
1617 const rvec * gmx_restrict x,
1618 rvec * gmx_restrict f, rvec * gmx_restrict fshift,
1619 gmx_bool VirCorr, matrix vir,
1620 t_nrnb *nrnb, const t_idef *idef,
1621 int ePBC, gmx_bool bMolPBC, const t_graph *g, const matrix box,
1622 const t_commrec *cr, gmx_wallcycle *wcycle)
1624 wallcycle_start(wcycle, ewcVSITESPREAD);
1625 const bool useDomdec = vsite->useDomdec;
1626 GMX_ASSERT(!useDomdec || (cr != nullptr && DOMAINDECOMP(cr)), "When vsites are set up with domain decomposition, we need a valid commrec");
1628 t_pbc pbc, *pbc_null;
1630 /* We only need to do pbc when we have inter-cg vsites */
1631 if ((useDomdec || bMolPBC) && vsite->n_intercg_vsite)
1633 /* This is wasting some CPU time as we now do this multiple times
1636 pbc_null = set_pbc_dd(&pbc, ePBC, useDomdec ? cr->dd->nc : nullptr, FALSE, box);
1645 dd_clear_f_vsites(cr->dd, f);
1648 if (vsite->nthreads == 1)
1655 spread_vsite_f_thread(vsite,
1658 idef->iparams, idef->il,
1663 for (int i = 0; i < DIM; i++)
1665 for (int j = 0; j < DIM; j++)
1667 vir[i][j] += -0.5*dxdf[i][j];
1674 /* First spread the vsites that might depend on non-local vsites */
1677 clear_mat(vsite->tData[vsite->nthreads]->dxdf);
1679 spread_vsite_f_thread(vsite,
1681 VirCorr, vsite->tData[vsite->nthreads]->dxdf,
1683 vsite->tData[vsite->nthreads]->ilist,
1686 #pragma omp parallel num_threads(vsite->nthreads)
1690 int thread = gmx_omp_get_thread_num();
1691 VsiteThread &tData = *vsite->tData[thread];
1694 if (thread == 0 || fshift == nullptr)
1700 fshift_t = tData.fshift;
1702 for (int i = 0; i < SHIFTS; i++)
1704 clear_rvec(fshift_t[i]);
1709 clear_mat(tData.dxdf);
1712 if (tData.useInterdependentTask)
1714 /* Spread the vsites that spread outside our local range.
1715 * This is done using a thread-local force buffer force.
1716 * First we need to copy the input vsite forces to force.
1718 InterdependentTask *idTask = &tData.idTask;
1720 /* Clear the buffer elements set by our task during
1721 * the last call to spread_vsite_f.
1723 clearTaskForceBufferUsedElements(idTask);
1725 int nvsite = idTask->vsite.size();
1726 for (int i = 0; i < nvsite; i++)
1728 copy_rvec(f[idTask->vsite[i]],
1729 idTask->force[idTask->vsite[i]]);
1731 spread_vsite_f_thread(vsite,
1732 x, as_rvec_array(idTask->force.data()), fshift_t,
1733 VirCorr, tData.dxdf,
1738 /* We need a barrier before reducing forces below
1739 * that have been produced by a different thread above.
1743 /* Loop over all thread task and reduce forces they
1744 * produced on atoms that fall in our range.
1745 * Note that atomic reduction would be a simpler solution,
1746 * but that might not have good support on all platforms.
1748 int ntask = idTask->reduceTask.size();
1749 for (int ti = 0; ti < ntask; ti++)
1751 const InterdependentTask *idt_foreign = &vsite->tData[idTask->reduceTask[ti]]->idTask;
1752 const AtomIndex *atomList = &idt_foreign->atomIndex[thread];
1753 const RVec *f_foreign = idt_foreign->force.data();
1755 int natom = atomList->atom.size();
1756 for (int i = 0; i < natom; i++)
1758 int ind = atomList->atom[i];
1759 rvec_inc(f[ind], f_foreign[ind]);
1760 /* Clearing of f_foreign is done at the next step */
1763 /* Clear the vsite forces, both in f and force */
1764 for (int i = 0; i < nvsite; i++)
1766 int ind = tData.idTask.vsite[i];
1768 clear_rvec(tData.idTask.force[ind]);
1772 /* Spread the vsites that spread locally only */
1773 spread_vsite_f_thread(vsite,
1775 VirCorr, tData.dxdf,
1780 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1783 if (fshift != nullptr)
1785 for (int th = 1; th < vsite->nthreads; th++)
1787 for (int i = 0; i < SHIFTS; i++)
1789 rvec_inc(fshift[i], vsite->tData[th]->fshift[i]);
1796 for (int th = 0; th < vsite->nthreads + 1; th++)
1798 /* MSVC doesn't like matrix references, so we use a pointer */
1799 const matrix *dxdf = &vsite->tData[th]->dxdf;
1801 for (int i = 0; i < DIM; i++)
1803 for (int j = 0; j < DIM; j++)
1805 vir[i][j] += -0.5*(*dxdf)[i][j];
1814 dd_move_f_vsites(cr->dd, f, fshift);
1817 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(idef->il, F_VSITE2));
1818 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(idef->il, F_VSITE3));
1819 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(idef->il, F_VSITE3FD));
1820 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(idef->il, F_VSITE3FAD));
1821 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(idef->il, F_VSITE3OUT));
1822 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(idef->il, F_VSITE4FD));
1823 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(idef->il, F_VSITE4FDN));
1824 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(idef->il, F_VSITEN));
1826 wallcycle_stop(wcycle, ewcVSITESPREAD);
1829 /*! \brief Returns the an array with charge-group indices for each atom
1831 * \param[in] chargeGroups The charge group block struct
1833 static std::vector<int> atom2cg(const t_block &chargeGroups)
1835 std::vector<int> a2cg(chargeGroups.index[chargeGroups.nr], 0);
1837 for (int chargeGroup = 0; chargeGroup < chargeGroups.nr; chargeGroup++)
1839 std::fill(a2cg.begin() + chargeGroups.index[chargeGroup],
1840 a2cg.begin() + chargeGroups.index[chargeGroup + 1],
1847 int count_intercg_vsites(const gmx_mtop_t *mtop)
1849 int n_intercg_vsite = 0;
1850 for (const gmx_molblock_t &molb : mtop->molblock)
1852 const gmx_moltype_t &molt = mtop->moltype[molb.type];
1854 std::vector<int> a2cg = atom2cg(molt.cgs);
1855 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
1857 const int nral = NRAL(ftype);
1858 const InteractionList &il = molt.ilist[ftype];
1859 for (int i = 0; i < il.size(); i += 1 + nral)
1861 int cg = a2cg[il.iatoms[1 + i]];
1862 for (int a = 1; a < nral; a++)
1864 if (a2cg[il.iatoms[1 + a]] != cg)
1866 n_intercg_vsite += molb.nmol;
1874 return n_intercg_vsite;
1877 template <typename T>
1879 get_vsite_pbc(const t_iparams *iparams, const T *ilist,
1880 const t_atom *atom, const t_mdatoms *md,
1883 /* Make an atom to charge group index */
1884 std::vector<int> a2cg = atom2cg(cgs);
1886 /* Make an array that tells if the pbc of an atom is set */
1887 std::vector<bool> pbc_set(cgs.index[cgs.nr], false);
1888 /* PBC is set for all non vsites */
1889 for (int a = 0; a < cgs.index[cgs.nr]; a++)
1891 if ((atom && atom[a].ptype != eptVSite) ||
1892 (md && md->ptype[a] != eptVSite))
1900 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
1903 int nral = NRAL(ftype);
1904 const T &il = ilist[ftype];
1906 std::vector<int> &vsite_pbc_f = vsite_pbc[ftype - F_VSITE2];
1907 vsite_pbc_f.resize(il.size()/(1 + nral));
1910 while (i < il.size())
1912 int vsi = i/(1 + nral);
1913 t_iatom vsite = il.iatoms[i + 1];
1914 int cg_v = a2cg[vsite];
1915 /* A value of -2 signals that this vsite and its contructing
1916 * atoms are all within the same cg, so no pbc is required.
1918 vsite_pbc_f[vsi] = -2;
1919 /* Check if constructing atoms are outside the vsite's cg */
1921 if (ftype == F_VSITEN)
1923 nc3 = 3*iparams[il.iatoms[i]].vsiten.n;
1924 for (int j = 0; j < nc3; j += 3)
1926 if (a2cg[il.iatoms[i + j + 2]] != cg_v)
1928 vsite_pbc_f[vsi] = -1;
1934 for (int a = 1; a < nral; a++)
1936 if (a2cg[il.iatoms[i + 1 + a]] != cg_v)
1938 vsite_pbc_f[vsi] = -1;
1942 if (vsite_pbc_f[vsi] == -1)
1944 /* Check if this is the first processed atom of a vsite only cg */
1945 gmx_bool bViteOnlyCG_and_FirstAtom = TRUE;
1946 for (int a = cgs.index[cg_v]; a < cgs.index[cg_v + 1]; a++)
1948 /* Non-vsites already have pbc set, so simply check for pbc_set */
1951 bViteOnlyCG_and_FirstAtom = FALSE;
1955 if (bViteOnlyCG_and_FirstAtom)
1957 /* First processed atom of a vsite only charge group.
1958 * The pbc of the input coordinates to construct_vsites
1959 * should be preserved.
1961 vsite_pbc_f[vsi] = vsite;
1963 else if (cg_v != a2cg[il.iatoms[1 + i + 1]])
1965 /* This vsite has a different charge group index
1966 * than it's first constructing atom
1967 * and the charge group has more than one atom,
1968 * search for the first normal particle
1969 * or vsite that already had its pbc defined.
1970 * If nothing is found, use full pbc for this vsite.
1972 for (int a = cgs.index[cg_v]; a < cgs.index[cg_v + 1]; a++)
1974 if (a != vsite && pbc_set[a])
1976 vsite_pbc_f[vsi] = a;
1979 fprintf(debug, "vsite %d match pbc with atom %d\n",
1987 fprintf(debug, "vsite atom %d cg %d - %d pbc atom %d\n",
1988 vsite+1, cgs.index[cg_v] + 1, cgs.index[cg_v + 1],
1989 vsite_pbc_f[vsi] + 1);
1993 if (ftype == F_VSITEN)
1995 /* The other entries in vsite_pbc_f are not used for center vsites */
2003 /* This vsite now has its pbc defined */
2004 pbc_set[vsite] = true;
2013 std::unique_ptr<gmx_vsite_t>
2014 initVsite(const gmx_mtop_t &mtop,
2015 const t_commrec *cr)
2017 GMX_RELEASE_ASSERT(cr != nullptr, "We need a valid commrec");
2019 std::unique_ptr<gmx_vsite_t> vsite;
2021 /* check if there are vsites */
2023 for (int ftype = 0; ftype < F_NRE; ftype++)
2025 if (interaction_function[ftype].flags & IF_VSITE)
2027 GMX_ASSERT(ftype >= c_ftypeVsiteStart && ftype < c_ftypeVsiteEnd, "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
2029 nvsite += gmx_mtop_ftype_count(&mtop, ftype);
2033 GMX_ASSERT(ftype < c_ftypeVsiteStart || ftype >= c_ftypeVsiteEnd, "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
2042 vsite = gmx::compat::make_unique<gmx_vsite_t>();
2044 vsite->n_intercg_vsite = count_intercg_vsites(&mtop);
2046 vsite->bHaveChargeGroups = (ncg_mtop(&mtop) < mtop.natoms);
2048 vsite->useDomdec = (DOMAINDECOMP(cr) && cr->dd->nnodes > 1);
2050 /* If we don't have charge groups, the vsite follows its own pbc.
2052 * With charge groups, each vsite needs to follow the pbc of the charge
2053 * group. Thus we need to keep track of PBC. Here we assume that without
2054 * domain decomposition all molecules are whole (which will not be
2055 * the case with periodic molecules).
2057 if (vsite->bHaveChargeGroups &&
2058 vsite->n_intercg_vsite > 0 &&
2061 vsite->vsite_pbc_molt.resize(mtop.moltype.size());
2062 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
2064 const gmx_moltype_t &molt = mtop.moltype[mt];
2065 vsite->vsite_pbc_molt[mt] = get_vsite_pbc(mtop.ffparams.iparams,
2067 molt.atoms.atom, nullptr,
2071 vsite->vsite_pbc_loc = gmx::compat::make_unique<VsitePbc>();
2074 vsite->nthreads = gmx_omp_nthreads_get(emntVSITE);
2076 if (vsite->nthreads > 1)
2078 /* We need one extra thread data structure for the overlap vsites */
2079 vsite->tData.resize(vsite->nthreads + 1);
2080 #pragma omp parallel for num_threads(vsite->nthreads) schedule(static)
2081 for (int thread = 0; thread < vsite->nthreads; thread++)
2085 vsite->tData[thread] = gmx::compat::make_unique<VsiteThread>();
2087 InterdependentTask &idTask = vsite->tData[thread]->idTask;
2089 idTask.atomIndex.resize(vsite->nthreads);
2091 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2093 if (vsite->nthreads > 1)
2095 vsite->tData[vsite->nthreads] = gmx::compat::make_unique<VsiteThread>();
2102 gmx_vsite_t::gmx_vsite_t()
2106 gmx_vsite_t::~gmx_vsite_t()
2110 static inline void flagAtom(InterdependentTask *idTask, int atom,
2111 int thread, int nthread, int natperthread)
2113 if (!idTask->use[atom])
2115 idTask->use[atom] = true;
2116 thread = atom/natperthread;
2117 /* Assign all non-local atom force writes to thread 0 */
2118 if (thread >= nthread)
2122 idTask->atomIndex[thread].atom.push_back(atom);
2126 /*\brief Here we try to assign all vsites that are in our local range.
2128 * Our task local atom range is tData->rangeStart - tData->rangeEnd.
2129 * Vsites that depend only on local atoms, as indicated by taskIndex[]==thread,
2130 * are assigned to task tData->ilist. Vsites that depend on non-local atoms
2131 * but not on other vsites are assigned to task tData->id_task.ilist.
2132 * taskIndex[] is set for all vsites in our range, either to our local tasks
2133 * or to the single last task as taskIndex[]=2*nthreads.
2135 static void assignVsitesToThread(VsiteThread *tData,
2139 gmx::ArrayRef<int> taskIndex,
2140 const t_ilist *ilist,
2141 const t_iparams *ip,
2142 const unsigned short *ptype)
2144 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2146 tData->ilist[ftype].nr = 0;
2147 tData->idTask.ilist[ftype].nr = 0;
2149 int nral1 = 1 + NRAL(ftype);
2151 t_iatom *iat = ilist[ftype].iatoms;
2152 for (int i = 0; i < ilist[ftype].nr; )
2154 if (ftype == F_VSITEN)
2156 /* The 3 below is from 1+NRAL(ftype)=3 */
2157 inc = ip[iat[i]].vsiten.n*3;
2160 if (iat[1 + i] < tData->rangeStart ||
2161 iat[1 + i] >= tData->rangeEnd)
2163 /* This vsite belongs to a different thread */
2168 /* We would like to assign this vsite to task thread,
2169 * but it might depend on atoms outside the atom range of thread
2170 * or on another vsite not assigned to task thread.
2173 if (ftype != F_VSITEN)
2175 for (int j = i + 2; j < i + nral1; j++)
2177 /* Do a range check to avoid a harmless race on taskIndex */
2178 if (iat[j] < tData->rangeStart ||
2179 iat[j] >= tData->rangeEnd ||
2180 taskIndex[iat[j]] != thread)
2182 if (!tData->useInterdependentTask ||
2183 ptype[iat[j]] == eptVSite)
2185 /* At least one constructing atom is a vsite
2186 * that is not assigned to the same thread.
2187 * Put this vsite into a separate task.
2193 /* There are constructing atoms outside our range,
2194 * put this vsite into a second task to be executed
2195 * on the same thread. During construction no barrier
2196 * is needed between the two tasks on the same thread.
2197 * During spreading we need to run this task with
2198 * an additional thread-local intermediate force buffer
2199 * (or atomic reduction) and a barrier between the two
2202 task = nthread + thread;
2208 for (int j = i + 2; j < i + inc; j += 3)
2210 /* Do a range check to avoid a harmless race on taskIndex */
2211 if (iat[j] < tData->rangeStart ||
2212 iat[j] >= tData->rangeEnd ||
2213 taskIndex[iat[j]] != thread)
2215 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");
2217 task = nthread + thread;
2222 /* Update this vsite's thread index entry */
2223 taskIndex[iat[1+i]] = task;
2225 if (task == thread || task == nthread + thread)
2227 /* Copy this vsite to the thread data struct of thread */
2231 il_task = &tData->ilist[ftype];
2235 il_task = &tData->idTask.ilist[ftype];
2237 /* Ensure we have sufficient memory allocated */
2238 if (il_task->nr + inc > il_task->nalloc)
2240 il_task->nalloc = over_alloc_large(il_task->nr + inc);
2241 srenew(il_task->iatoms, il_task->nalloc);
2243 /* Copy the vsite data to the thread-task local array */
2244 for (int j = i; j < i + inc; j++)
2246 il_task->iatoms[il_task->nr++] = iat[j];
2248 if (task == nthread + thread)
2250 /* This vsite write outside our own task force block.
2251 * Put it into the interdependent task list and flag
2252 * the atoms involved for reduction.
2254 tData->idTask.vsite.push_back(iat[i + 1]);
2255 if (ftype != F_VSITEN)
2257 for (int j = i + 2; j < i + nral1; j++)
2259 flagAtom(&tData->idTask, iat[j],
2260 thread, nthread, natperthread);
2265 for (int j = i + 2; j < i + inc; j += 3)
2267 flagAtom(&tData->idTask, iat[j],
2268 thread, nthread, natperthread);
2279 /*! \brief Assign all vsites with taskIndex[]==task to task tData */
2280 static void assignVsitesToSingleTask(VsiteThread *tData,
2282 gmx::ArrayRef<const int> taskIndex,
2283 const t_ilist *ilist,
2284 const t_iparams *ip)
2286 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2288 tData->ilist[ftype].nr = 0;
2289 tData->idTask.ilist[ftype].nr = 0;
2291 int nral1 = 1 + NRAL(ftype);
2293 t_iatom *iat = ilist[ftype].iatoms;
2294 t_ilist *il_task = &tData->ilist[ftype];
2296 for (int i = 0; i < ilist[ftype].nr; )
2298 if (ftype == F_VSITEN)
2300 /* The 3 below is from 1+NRAL(ftype)=3 */
2301 inc = ip[iat[i]].vsiten.n*3;
2303 /* Check if the vsite is assigned to our task */
2304 if (taskIndex[iat[1 + i]] == task)
2306 /* Ensure we have sufficient memory allocated */
2307 if (il_task->nr + inc > il_task->nalloc)
2309 il_task->nalloc = over_alloc_large(il_task->nr + inc);
2310 srenew(il_task->iatoms, il_task->nalloc);
2312 /* Copy the vsite data to the thread-task local array */
2313 for (int j = i; j < i + inc; j++)
2315 il_task->iatoms[il_task->nr++] = iat[j];
2324 void split_vsites_over_threads(const t_ilist *ilist,
2325 const t_iparams *ip,
2326 const t_mdatoms *mdatoms,
2329 int vsite_atom_range, natperthread;
2331 if (vsite->nthreads == 1)
2337 /* The current way of distributing the vsites over threads in primitive.
2338 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
2339 * without taking into account how the vsites are distributed.
2340 * Without domain decomposition we at least tighten the upper bound
2341 * of the range (useful for common systems such as a vsite-protein
2343 * With domain decomposition, as long as the vsites are distributed
2344 * uniformly in each domain along the major dimension, usually x,
2345 * it will also perform well.
2347 if (!vsite->useDomdec)
2349 vsite_atom_range = -1;
2350 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2353 if (ftype != F_VSITEN)
2355 int nral1 = 1 + NRAL(ftype);
2356 const t_iatom *iat = ilist[ftype].iatoms;
2357 for (int i = 0; i < ilist[ftype].nr; i += nral1)
2359 for (int j = i + 1; j < i + nral1; j++)
2361 vsite_atom_range = std::max(vsite_atom_range, iat[j]);
2369 const t_iatom *iat = ilist[ftype].iatoms;
2372 while (i < ilist[ftype].nr)
2374 /* The 3 below is from 1+NRAL(ftype)=3 */
2375 vs_ind_end = i + ip[iat[i]].vsiten.n*3;
2377 vsite_atom_range = std::max(vsite_atom_range, iat[i+1]);
2378 while (i < vs_ind_end)
2380 vsite_atom_range = std::max(vsite_atom_range, iat[i+2]);
2388 natperthread = (vsite_atom_range + vsite->nthreads - 1)/vsite->nthreads;
2392 /* Any local or not local atom could be involved in virtual sites.
2393 * But since we usually have very few non-local virtual sites
2394 * (only non-local vsites that depend on local vsites),
2395 * we distribute the local atom range equally over the threads.
2396 * When assigning vsites to threads, we should take care that the last
2397 * threads also covers the non-local range.
2399 vsite_atom_range = mdatoms->nr;
2400 natperthread = (mdatoms->homenr + vsite->nthreads - 1)/vsite->nthreads;
2405 fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n", mdatoms->nr, vsite_atom_range, natperthread);
2408 /* To simplify the vsite assignment, we make an index which tells us
2409 * to which task particles, both non-vsites and vsites, are assigned.
2411 vsite->taskIndex.resize(mdatoms->nr);
2413 /* Initialize the task index array. Here we assign the non-vsite
2414 * particles to task=thread, so we easily figure out if vsites
2415 * depend on local and/or non-local particles in assignVsitesToThread.
2417 gmx::ArrayRef<int> taskIndex = vsite->taskIndex;
2420 for (int i = 0; i < mdatoms->nr; i++)
2422 if (mdatoms->ptype[i] == eptVSite)
2424 /* vsites are not assigned to a task yet */
2429 /* assign non-vsite particles to task thread */
2430 taskIndex[i] = thread;
2432 if (i == (thread + 1)*natperthread && thread < vsite->nthreads)
2439 #pragma omp parallel num_threads(vsite->nthreads)
2443 int thread = gmx_omp_get_thread_num();
2444 VsiteThread &tData = *vsite->tData[thread];
2446 /* Clear the buffer use flags that were set before */
2447 if (tData.useInterdependentTask)
2449 InterdependentTask &idTask = tData.idTask;
2451 /* To avoid an extra OpenMP barrier in spread_vsite_f,
2452 * we clear the force buffer at the next step,
2453 * so we need to do it here as well.
2455 clearTaskForceBufferUsedElements(&idTask);
2457 idTask.vsite.resize(0);
2458 for (int t = 0; t < vsite->nthreads; t++)
2460 AtomIndex &atomIndex = idTask.atomIndex[t];
2461 int natom = atomIndex.atom.size();
2462 for (int i = 0; i < natom; i++)
2464 idTask.use[atomIndex.atom[i]] = false;
2466 atomIndex.atom.resize(0);
2471 /* To avoid large f_buf allocations of #threads*vsite_atom_range
2472 * we don't use task2 with more than 200000 atoms. This doesn't
2473 * affect performance, since with such a large range relatively few
2474 * vsites will end up in the separate task.
2475 * Note that useTask2 should be the same for all threads.
2477 tData.useInterdependentTask = (vsite_atom_range <= 200000);
2478 if (tData.useInterdependentTask)
2480 size_t natoms_use_in_vsites = vsite_atom_range;
2481 InterdependentTask &idTask = tData.idTask;
2482 /* To avoid resizing and re-clearing every nstlist steps,
2483 * we never down size the force buffer.
2485 if (natoms_use_in_vsites > idTask.force.size() ||
2486 natoms_use_in_vsites > idTask.use.size())
2488 idTask.force.resize(natoms_use_in_vsites, { 0, 0, 0 });
2489 idTask.use.resize(natoms_use_in_vsites, false);
2493 /* Assign all vsites that can execute independently on threads */
2494 tData.rangeStart = thread *natperthread;
2495 if (thread < vsite->nthreads - 1)
2497 tData.rangeEnd = (thread + 1)*natperthread;
2501 /* The last thread should cover up to the end of the range */
2502 tData.rangeEnd = mdatoms->nr;
2504 assignVsitesToThread(&tData,
2505 thread, vsite->nthreads,
2508 ilist, ip, mdatoms->ptype);
2510 if (tData.useInterdependentTask)
2512 /* In the worst case, all tasks write to force ranges of
2513 * all other tasks, leading to #tasks^2 scaling (this is only
2514 * the overhead, the actual flops remain constant).
2515 * But in most cases there is far less coupling. To improve
2516 * scaling at high thread counts we therefore construct
2517 * an index to only loop over the actually affected tasks.
2519 InterdependentTask &idTask = tData.idTask;
2521 /* Ensure assignVsitesToThread finished on other threads */
2524 idTask.spreadTask.resize(0);
2525 idTask.reduceTask.resize(0);
2526 for (int t = 0; t < vsite->nthreads; t++)
2528 /* Do we write to the force buffer of task t? */
2529 if (!idTask.atomIndex[t].atom.empty())
2531 idTask.spreadTask.push_back(t);
2533 /* Does task t write to our force buffer? */
2534 if (!vsite->tData[t]->idTask.atomIndex[thread].atom.empty())
2536 idTask.reduceTask.push_back(t);
2541 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2543 /* Assign all remaining vsites, that will have taskIndex[]=2*vsite->nthreads,
2544 * to a single task that will not run in parallel with other tasks.
2546 assignVsitesToSingleTask(vsite->tData[vsite->nthreads].get(),
2551 if (debug && vsite->nthreads > 1)
2553 fprintf(debug, "virtual site useInterdependentTask %d, nuse:\n",
2554 static_cast<int>(vsite->tData[0]->useInterdependentTask));
2555 for (int th = 0; th < vsite->nthreads + 1; th++)
2557 fprintf(debug, " %4d", vsite->tData[th]->idTask.nuse);
2559 fprintf(debug, "\n");
2561 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2563 if (ilist[ftype].nr > 0)
2565 fprintf(debug, "%-20s thread dist:",
2566 interaction_function[ftype].longname);
2567 for (int th = 0; th < vsite->nthreads + 1; th++)
2569 fprintf(debug, " %4d %4d ",
2570 vsite->tData[th]->ilist[ftype].nr,
2571 vsite->tData[th]->idTask.ilist[ftype].nr);
2573 fprintf(debug, "\n");
2579 int nrOrig = vsiteIlistNrCount(ilist);
2581 for (int th = 0; th < vsite->nthreads + 1; th++)
2584 vsiteIlistNrCount(vsite->tData[th]->ilist) +
2585 vsiteIlistNrCount(vsite->tData[th]->idTask.ilist);
2587 GMX_ASSERT(nrThreaded == nrOrig, "The number of virtual sites assigned to all thread task has to match the total number of virtual sites");
2591 void set_vsite_top(gmx_vsite_t *vsite,
2592 const gmx_localtop_t *top,
2593 const t_mdatoms *md)
2595 if (vsite->n_intercg_vsite > 0 && vsite->bHaveChargeGroups)
2597 vsite->vsite_pbc_loc = gmx::compat::make_unique<VsitePbc>();
2598 *vsite->vsite_pbc_loc = get_vsite_pbc(top->idef.iparams,
2599 top->idef.il, nullptr, md,
2603 if (vsite->nthreads > 1)
2605 if (vsite->bHaveChargeGroups)
2607 gmx_fatal(FARGS, "The combination of threading, virtual sites and charge groups is not implemented");
2610 split_vsites_over_threads(top->idef.il, top->idef.iparams,