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 The GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Implements the VirtualSitesHandler class and vsite standalone functions
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_mdlib
55 #include "gromacs/domdec/domdec.h"
56 #include "gromacs/domdec/domdec_struct.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/gmxlib/nrnb.h"
59 #include "gromacs/math/functions.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/mdlib/gmx_omp_nthreads.h"
62 #include "gromacs/mdtypes/commrec.h"
63 #include "gromacs/mdtypes/mdatom.h"
64 #include "gromacs/pbcutil/ishift.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/timing/wallcycle.h"
67 #include "gromacs/topology/ifunc.h"
68 #include "gromacs/topology/mtop_util.h"
69 #include "gromacs/topology/topology.h"
70 #include "gromacs/utility/exceptions.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/gmxassert.h"
73 #include "gromacs/utility/gmxomp.h"
75 /* The strategy used here for assigning virtual sites to (thread-)tasks
78 * We divide the atom range that vsites operate on (natoms_local with DD,
79 * 0 - last atom involved in vsites without DD) equally over all threads.
81 * Vsites in the local range constructed from atoms in the local range
82 * and/or other vsites that are fully local are assigned to a simple,
85 * Vsites that are not assigned after using the above criterion get assigned
86 * to a so called "interdependent" thread task when none of the constructing
87 * atoms is a vsite. These tasks are called interdependent, because one task
88 * accesses atoms assigned to a different task/thread.
89 * Note that this option is turned off with large (local) atom counts
90 * to avoid high memory usage.
92 * Any remaining vsites are assigned to a separate master thread task.
98 //! VirialHandling is often used outside VirtualSitesHandler class members
99 using VirialHandling = VirtualSitesHandler::VirialHandling;
102 * \brief Information on PBC and domain decomposition for virtual sites
107 //! Constructs without PBC and DD
108 DomainInfo() = default;
110 //! Constructs with PBC and DD, if !=nullptr
111 DomainInfo(PbcType pbcType, bool haveInterUpdateGroupVirtualSites, gmx_domdec_t* domdec) :
113 useMolPbc_(pbcType != PbcType::No && haveInterUpdateGroupVirtualSites),
118 //! Returns whether we are using domain decomposition with more than 1 DD rank
119 bool useDomdec() const { return (domdec_ != nullptr); }
122 const PbcType pbcType_ = PbcType::No;
123 //! Whether molecules are broken over PBC
124 const bool useMolPbc_ = false;
125 //! Pointer to the domain decomposition struct, nullptr without PP DD
126 const gmx_domdec_t* domdec_ = nullptr;
130 * \brief List of atom indices belonging to a task
134 //! List of atom indices
135 std::vector<int> atom;
139 * \brief Data structure for thread tasks that use constructing atoms outside their own atom range
141 struct InterdependentTask
143 //! The interaction lists, only vsite entries are used
144 InteractionLists ilist;
145 //! Thread/task-local force buffer
146 std::vector<RVec> force;
147 //! The atom indices of the vsites of our task
148 std::vector<int> vsite;
149 //! Flags if elements in force are spread to or not
150 std::vector<bool> use;
151 //! The number of entries set to true in use
153 //! Array of atoms indices, size nthreads, covering all nuse set elements in use
154 std::vector<AtomIndex> atomIndex;
155 //! List of tasks (force blocks) this task spread forces to
156 std::vector<int> spreadTask;
157 //! List of tasks that write to this tasks force block range
158 std::vector<int> reduceTask;
162 * \brief Vsite thread task data structure
166 //! Start of atom range of this task
168 //! End of atom range of this task
170 //! The interaction lists, only vsite entries are used
171 std::array<InteractionList, F_NRE> ilist;
172 //! Local fshift accumulation buffer
173 std::array<RVec, SHIFTS> fshift;
174 //! Local virial dx*df accumulation buffer
176 //! 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
177 bool useInterdependentTask;
178 //! Data for vsites that involve constructing atoms in the atom range of other threads/tasks
179 InterdependentTask idTask;
181 /*! \brief Constructor */
186 for (auto& elem : fshift)
188 elem = { 0.0_real, 0.0_real, 0.0_real };
191 useInterdependentTask = false;
197 * \brief Information on how the virtual site work is divided over thread tasks
202 //! Constructor, retrieves the number of threads to use from gmx_omp_nthreads.h
205 //! Returns the number of threads to use for vsite operations
206 int numThreads() const { return numThreads_; }
208 //! Returns the thread data for the given thread
209 const VsiteThread& threadData(int threadIndex) const { return *tData_[threadIndex]; }
211 //! Returns the thread data for the given thread
212 VsiteThread& threadData(int threadIndex) { return *tData_[threadIndex]; }
214 //! Returns the thread data for vsites that depend on non-local vsites
215 const VsiteThread& threadDataNonLocalDependent() const { return *tData_[numThreads_]; }
217 //! Returns the thread data for vsites that depend on non-local vsites
218 VsiteThread& threadDataNonLocalDependent() { return *tData_[numThreads_]; }
220 //! Set VSites and distribute VSite work over threads, should be called after DD partitioning
221 void setVirtualSites(ArrayRef<const InteractionList> ilist,
222 ArrayRef<const t_iparams> iparams,
223 const t_mdatoms& mdatoms,
227 //! Number of threads used for vsite operations
228 const int numThreads_;
229 //! Thread local vsites and work structs
230 std::vector<std::unique_ptr<VsiteThread>> tData_;
231 //! Work array for dividing vsites over threads
232 std::vector<int> taskIndex_;
236 * \brief Impl class for VirtualSitesHandler
238 class VirtualSitesHandler::Impl
241 //! Constructor, domdec should be nullptr without DD
242 Impl(const gmx_mtop_t& mtop, gmx_domdec_t* domdec, PbcType pbcType);
244 //! Returns the number of virtual sites acting over multiple update groups
245 int numInterUpdategroupVirtualSites() const { return numInterUpdategroupVirtualSites_; }
247 //! Set VSites and distribute VSite work over threads, should be called after DD partitioning
248 void setVirtualSites(ArrayRef<const InteractionList> ilist, const t_mdatoms& mdatoms);
250 /*! \brief Create positions of vsite atoms based for the local system
252 * \param[in,out] x The coordinates
253 * \param[in] dt The time step
254 * \param[in,out] v When != nullptr, velocities for vsites are set as displacement/dt
255 * \param[in] box The box
257 void construct(ArrayRef<RVec> x, real dt, ArrayRef<RVec> v, const matrix box) const;
259 /*! \brief Spread the force operating on the vsite atoms on the surrounding atoms.
261 * vsite should point to a valid object.
262 * The virialHandling parameter determines how virial contributions are handled.
263 * If this is set to Linear, shift forces are accumulated into fshift.
264 * If this is set to NonLinear, non-linear contributions are added to virial.
265 * This non-linear correction is required when the virial is not calculated
266 * afterwards from the particle position and forces, but in a different way,
267 * as for instance for the PME mesh contribution.
269 void spreadForces(ArrayRef<const RVec> x,
271 VirialHandling virialHandling,
272 ArrayRef<RVec> fshift,
276 gmx_wallcycle* wcycle);
279 //! The number of vsites that cross update groups, when =0 no PBC treatment is needed
280 const int numInterUpdategroupVirtualSites_;
281 //! PBC and DD information
282 const DomainInfo domainInfo_;
283 //! The interaction parameters
284 const ArrayRef<const t_iparams> iparams_;
285 //! The interaction lists
286 ArrayRef<const InteractionList> ilists_;
287 //! Information for handling vsite threading
288 ThreadingInfo threadingInfo_;
291 VirtualSitesHandler::~VirtualSitesHandler() = default;
293 int VirtualSitesHandler::numInterUpdategroupVirtualSites() const
295 return impl_->numInterUpdategroupVirtualSites();
298 /*! \brief Returns the sum of the vsite ilist sizes over all vsite types
300 * \param[in] ilist The interaction list
302 static int vsiteIlistNrCount(ArrayRef<const InteractionList> ilist)
305 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
307 nr += ilist[ftype].size();
313 //! Computes the distance between xi and xj, pbc is used when pbc!=nullptr
314 static int pbc_rvec_sub(const t_pbc* pbc, const rvec xi, const rvec xj, rvec dx)
318 return pbc_dx_aiuc(pbc, xi, xj, dx);
322 rvec_sub(xi, xj, dx);
327 //! Returns the 1/norm(x)
328 static inline real inverseNorm(const rvec x)
330 return gmx::invsqrt(iprod(x, x));
334 /* Vsite construction routines */
336 static void constr_vsite1(const rvec xi, rvec x)
343 static void constr_vsite2(const rvec xi, const rvec xj, rvec x, real a, const t_pbc* pbc)
351 pbc_dx_aiuc(pbc, xj, xi, dx);
352 x[XX] = xi[XX] + a * dx[XX];
353 x[YY] = xi[YY] + a * dx[YY];
354 x[ZZ] = xi[ZZ] + a * dx[ZZ];
358 x[XX] = b * xi[XX] + a * xj[XX];
359 x[YY] = b * xi[YY] + a * xj[YY];
360 x[ZZ] = b * xi[ZZ] + a * xj[ZZ];
364 /* TOTAL: 10 flops */
367 static void constr_vsite2FD(const rvec xi, const rvec xj, rvec x, real a, const t_pbc* pbc)
370 pbc_rvec_sub(pbc, xj, xi, xij);
373 const real b = a * inverseNorm(xij);
376 x[XX] = xi[XX] + b * xij[XX];
377 x[YY] = xi[YY] + b * xij[YY];
378 x[ZZ] = xi[ZZ] + b * xij[ZZ];
381 /* TOTAL: 25 flops */
384 static void constr_vsite3(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, const t_pbc* pbc)
393 pbc_dx_aiuc(pbc, xj, xi, dxj);
394 pbc_dx_aiuc(pbc, xk, xi, dxk);
395 x[XX] = xi[XX] + a * dxj[XX] + b * dxk[XX];
396 x[YY] = xi[YY] + a * dxj[YY] + b * dxk[YY];
397 x[ZZ] = xi[ZZ] + a * dxj[ZZ] + b * dxk[ZZ];
401 x[XX] = c * xi[XX] + a * xj[XX] + b * xk[XX];
402 x[YY] = c * xi[YY] + a * xj[YY] + b * xk[YY];
403 x[ZZ] = c * xi[ZZ] + a * xj[ZZ] + b * xk[ZZ];
407 /* TOTAL: 17 flops */
410 static void constr_vsite3FD(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, const t_pbc* pbc)
415 pbc_rvec_sub(pbc, xj, xi, xij);
416 pbc_rvec_sub(pbc, xk, xj, xjk);
419 /* temp goes from i to a point on the line jk */
420 temp[XX] = xij[XX] + a * xjk[XX];
421 temp[YY] = xij[YY] + a * xjk[YY];
422 temp[ZZ] = xij[ZZ] + a * xjk[ZZ];
425 c = b * inverseNorm(temp);
428 x[XX] = xi[XX] + c * temp[XX];
429 x[YY] = xi[YY] + c * temp[YY];
430 x[ZZ] = xi[ZZ] + c * temp[ZZ];
433 /* TOTAL: 34 flops */
436 static void constr_vsite3FAD(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, const t_pbc* pbc)
439 real a1, b1, c1, invdij;
441 pbc_rvec_sub(pbc, xj, xi, xij);
442 pbc_rvec_sub(pbc, xk, xj, xjk);
445 invdij = inverseNorm(xij);
446 c1 = invdij * invdij * iprod(xij, xjk);
447 xp[XX] = xjk[XX] - c1 * xij[XX];
448 xp[YY] = xjk[YY] - c1 * xij[YY];
449 xp[ZZ] = xjk[ZZ] - c1 * xij[ZZ];
451 b1 = b * inverseNorm(xp);
454 x[XX] = xi[XX] + a1 * xij[XX] + b1 * xp[XX];
455 x[YY] = xi[YY] + a1 * xij[YY] + b1 * xp[YY];
456 x[ZZ] = xi[ZZ] + a1 * xij[ZZ] + b1 * xp[ZZ];
459 /* TOTAL: 63 flops */
463 constr_vsite3OUT(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, real c, const t_pbc* pbc)
467 pbc_rvec_sub(pbc, xj, xi, xij);
468 pbc_rvec_sub(pbc, xk, xi, xik);
469 cprod(xij, xik, temp);
472 x[XX] = xi[XX] + a * xij[XX] + b * xik[XX] + c * temp[XX];
473 x[YY] = xi[YY] + a * xij[YY] + b * xik[YY] + c * temp[YY];
474 x[ZZ] = xi[ZZ] + a * xij[ZZ] + b * xik[ZZ] + c * temp[ZZ];
477 /* TOTAL: 33 flops */
480 static void constr_vsite4FD(const rvec xi,
490 rvec xij, xjk, xjl, temp;
493 pbc_rvec_sub(pbc, xj, xi, xij);
494 pbc_rvec_sub(pbc, xk, xj, xjk);
495 pbc_rvec_sub(pbc, xl, xj, xjl);
498 /* temp goes from i to a point on the plane jkl */
499 temp[XX] = xij[XX] + a * xjk[XX] + b * xjl[XX];
500 temp[YY] = xij[YY] + a * xjk[YY] + b * xjl[YY];
501 temp[ZZ] = xij[ZZ] + a * xjk[ZZ] + b * xjl[ZZ];
504 d = c * inverseNorm(temp);
507 x[XX] = xi[XX] + d * temp[XX];
508 x[YY] = xi[YY] + d * temp[YY];
509 x[ZZ] = xi[ZZ] + d * temp[ZZ];
512 /* TOTAL: 43 flops */
515 static void constr_vsite4FDN(const rvec xi,
525 rvec xij, xik, xil, ra, rb, rja, rjb, rm;
528 pbc_rvec_sub(pbc, xj, xi, xij);
529 pbc_rvec_sub(pbc, xk, xi, xik);
530 pbc_rvec_sub(pbc, xl, xi, xil);
533 ra[XX] = a * xik[XX];
534 ra[YY] = a * xik[YY];
535 ra[ZZ] = a * xik[ZZ];
537 rb[XX] = b * xil[XX];
538 rb[YY] = b * xil[YY];
539 rb[ZZ] = b * xil[ZZ];
543 rvec_sub(ra, xij, rja);
544 rvec_sub(rb, xij, rjb);
550 d = c * inverseNorm(rm);
553 x[XX] = xi[XX] + d * rm[XX];
554 x[YY] = xi[YY] + d * rm[YY];
555 x[ZZ] = xi[ZZ] + d * rm[ZZ];
558 /* TOTAL: 47 flops */
561 static int constr_vsiten(const t_iatom* ia, ArrayRef<const t_iparams> ip, ArrayRef<RVec> x, const t_pbc* pbc)
568 n3 = 3 * ip[ia[0]].vsiten.n;
571 copy_rvec(x[ai], x1);
573 for (int i = 3; i < n3; i += 3)
576 a = ip[ia[i]].vsiten.a;
579 pbc_dx_aiuc(pbc, x[ai], x1, dx);
583 rvec_sub(x[ai], x1, dx);
585 dsum[XX] += a * dx[XX];
586 dsum[YY] += a * dx[YY];
587 dsum[ZZ] += a * dx[ZZ];
591 x[av][XX] = x1[XX] + dsum[XX];
592 x[av][YY] = x1[YY] + dsum[YY];
593 x[av][ZZ] = x1[ZZ] + dsum[ZZ];
600 //! PBC modes for vsite construction and spreading
603 all, //!< Apply normal, simple PBC for all vsites
604 none //!< No PBC treatment needed
607 /*! \brief Returns the PBC mode based on the system PBC and vsite properties
609 * \param[in] pbcPtr A pointer to a PBC struct or nullptr when no PBC treatment is required
611 static PbcMode getPbcMode(const t_pbc* pbcPtr)
613 if (pbcPtr == nullptr)
615 return PbcMode::none;
623 /*! \brief Executes the vsite construction task for a single thread
625 * \param[in,out] x Coordinates to construct vsites for
626 * \param[in] dt Time step, needed when v is not empty
627 * \param[in,out] v When not empty, velocities are generated for virtual sites
628 * \param[in] ip Interaction parameters for all interaction, only vsite parameters are used
629 * \param[in] ilist The interaction lists, only vsites are usesd
630 * \param[in] pbc_null PBC struct, used for PBC distance calculations when !=nullptr
632 static void construct_vsites_thread(ArrayRef<RVec> x,
635 ArrayRef<const t_iparams> ip,
636 ArrayRef<const InteractionList> ilist,
637 const t_pbc* pbc_null)
649 const PbcMode pbcMode = getPbcMode(pbc_null);
650 /* We need another pbc pointer, as with charge groups we switch per vsite */
651 const t_pbc* pbc_null2 = pbc_null;
653 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
655 if (ilist[ftype].empty())
661 int nra = interaction_function[ftype].nratoms;
663 int nr = ilist[ftype].size();
665 const t_iatom* ia = ilist[ftype].iatoms.data();
667 for (int i = 0; i < nr;)
670 /* The vsite and constructing atoms */
673 /* Constants for constructing vsites */
674 real a1 = ip[tp].vsite.a;
675 /* Copy the old position */
677 copy_rvec(x[avsite], xv);
679 /* Construct the vsite depending on type */
684 case F_VSITE1: constr_vsite1(x[ai], x[avsite]); break;
687 constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
691 constr_vsite2FD(x[ai], x[aj], x[avsite], a1, pbc_null2);
697 constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
703 constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
709 constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
716 constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
724 constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1, pbc_null2);
732 constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1, pbc_null2);
734 case F_VSITEN: inc = constr_vsiten(ia, ip, x, pbc_null2); break;
736 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d", ftype, __FILE__, __LINE__);
739 if (pbcMode == PbcMode::all)
741 /* Keep the vsite in the same periodic image as before */
743 int ishift = pbc_dx_aiuc(pbc_null, x[avsite], xv, dx);
744 if (ishift != CENTRAL)
746 rvec_add(xv, dx, x[avsite]);
751 /* Calculate velocity of vsite... */
753 rvec_sub(x[avsite], xv, vv);
754 svmul(inv_dt, vv, v[avsite]);
757 /* Increment loop variables */
765 /*! \brief Dispatch the vsite construction tasks for all threads
767 * \param[in] threadingInfo Used to divide work over threads when != nullptr
768 * \param[in,out] x Coordinates to construct vsites for
769 * \param[in] dt Time step, needed when v is not empty
770 * \param[in,out] v When not empty, velocities are generated for virtual sites
771 * \param[in] ip Interaction parameters for all interaction, only vsite parameters are used
772 * \param[in] ilist The interaction lists, only vsites are usesd
773 * \param[in] domainInfo Information about PBC and DD
774 * \param[in] box Used for PBC when PBC is set in domainInfo
776 static void construct_vsites(const ThreadingInfo* threadingInfo,
780 ArrayRef<const t_iparams> ip,
781 ArrayRef<const InteractionList> ilist,
782 const DomainInfo& domainInfo,
785 const bool useDomdec = domainInfo.useDomdec();
787 t_pbc pbc, *pbc_null;
789 /* We only need to do pbc when we have inter update-group vsites.
790 * Note that with domain decomposition we do not need to apply PBC here
791 * when we have at least 3 domains along each dimension. Currently we
792 * do not optimize this case.
794 if (domainInfo.pbcType_ != PbcType::No && domainInfo.useMolPbc_)
796 /* This is wasting some CPU time as we now do this multiple times
800 clear_ivec(null_ivec);
801 pbc_null = set_pbc_dd(&pbc, domainInfo.pbcType_,
802 useDomdec ? domainInfo.domdec_->numCells : null_ivec, FALSE, box);
811 dd_move_x_vsites(*domainInfo.domdec_, box, as_rvec_array(x.data()));
814 if (threadingInfo == nullptr || threadingInfo->numThreads() == 1)
816 construct_vsites_thread(x, dt, v, ip, ilist, pbc_null);
820 #pragma omp parallel num_threads(threadingInfo->numThreads())
824 const int th = gmx_omp_get_thread_num();
825 const VsiteThread& tData = threadingInfo->threadData(th);
826 GMX_ASSERT(tData.rangeStart >= 0,
827 "The thread data should be initialized before calling construct_vsites");
829 construct_vsites_thread(x, dt, v, ip, tData.ilist, pbc_null);
830 if (tData.useInterdependentTask)
832 /* Here we don't need a barrier (unlike the spreading),
833 * since both tasks only construct vsites from particles,
834 * or local vsites, not from non-local vsites.
836 construct_vsites_thread(x, dt, v, ip, tData.idTask.ilist, pbc_null);
839 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
841 /* Now we can construct the vsites that might depend on other vsites */
842 construct_vsites_thread(x, dt, v, ip, threadingInfo->threadDataNonLocalDependent().ilist, pbc_null);
846 void VirtualSitesHandler::Impl::construct(ArrayRef<RVec> x, real dt, ArrayRef<RVec> v, const matrix box) const
848 construct_vsites(&threadingInfo_, x, dt, v, iparams_, ilists_, domainInfo_, box);
851 void VirtualSitesHandler::construct(ArrayRef<RVec> x, real dt, ArrayRef<RVec> v, const matrix box) const
853 impl_->construct(x, dt, v, box);
856 void constructVirtualSites(ArrayRef<RVec> x, ArrayRef<const t_iparams> ip, ArrayRef<const InteractionList> ilist)
860 const DomainInfo domainInfo;
861 construct_vsites(nullptr, x, 0, {}, ip, ilist, domainInfo, nullptr);
865 /* Force spreading routines */
867 static void spread_vsite1(const t_iatom ia[], ArrayRef<RVec> f)
869 const int av = ia[1];
870 const int ai = ia[2];
875 template<VirialHandling virialHandling>
876 static void spread_vsite2(const t_iatom ia[],
878 ArrayRef<const RVec> x,
880 ArrayRef<RVec> fshift,
890 svmul(1 - a, f[av], fi);
898 if (virialHandling == VirialHandling::Pbc)
904 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
905 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
913 if (siv != CENTRAL || sij != CENTRAL)
915 rvec_inc(fshift[siv], f[av]);
916 rvec_dec(fshift[CENTRAL], fi);
917 rvec_dec(fshift[sij], fj);
921 /* TOTAL: 13 flops */
924 void constructVirtualSitesGlobal(const gmx_mtop_t& mtop, gmx::ArrayRef<gmx::RVec> x)
926 GMX_ASSERT(x.ssize() >= mtop.natoms, "x should contain the whole system");
927 GMX_ASSERT(!mtop.moleculeBlockIndices.empty(),
928 "molblock indices are needed in constructVsitesGlobal");
930 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
932 const gmx_molblock_t& molb = mtop.molblock[mb];
933 const gmx_moltype_t& molt = mtop.moltype[molb.type];
934 if (vsiteIlistNrCount(molt.ilist) > 0)
936 int atomOffset = mtop.moleculeBlockIndices[mb].globalAtomStart;
937 for (int mol = 0; mol < molb.nmol; mol++)
939 constructVirtualSites(x.subArray(atomOffset, molt.atoms.nr), mtop.ffparams.iparams,
941 atomOffset += molt.atoms.nr;
947 template<VirialHandling virialHandling>
948 static void spread_vsite2FD(const t_iatom ia[],
950 ArrayRef<const RVec> x,
952 ArrayRef<RVec> fshift,
956 const int av = ia[1];
957 const int ai = ia[2];
958 const int aj = ia[3];
960 copy_rvec(f[av], fv);
963 int sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
966 const real invDistance = inverseNorm(xij);
967 const real b = a * invDistance;
970 const real fproj = iprod(xij, fv) * invDistance * invDistance;
973 fj[XX] = b * (fv[XX] - fproj * xij[XX]);
974 fj[YY] = b * (fv[YY] - fproj * xij[YY]);
975 fj[ZZ] = b * (fv[ZZ] - fproj * xij[ZZ]);
978 /* b is already calculated in constr_vsite2FD
979 storing b somewhere will save flops. */
981 f[ai][XX] += fv[XX] - fj[XX];
982 f[ai][YY] += fv[YY] - fj[YY];
983 f[ai][ZZ] += fv[ZZ] - fj[ZZ];
989 if (virialHandling == VirialHandling::Pbc)
995 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1002 if (svi != CENTRAL || sji != CENTRAL)
1004 rvec_dec(fshift[svi], fv);
1005 fshift[CENTRAL][XX] += fv[XX] - fj[XX];
1006 fshift[CENTRAL][YY] += fv[YY] - fj[YY];
1007 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ];
1008 fshift[sji][XX] += fj[XX];
1009 fshift[sji][YY] += fj[YY];
1010 fshift[sji][ZZ] += fj[ZZ];
1014 if (virialHandling == VirialHandling::NonLinear)
1016 /* Under this condition, the virial for the current forces is not
1017 * calculated from the redistributed forces. This means that
1018 * the effect of non-linear virtual site constructions on the virial
1019 * needs to be added separately. This contribution can be calculated
1020 * in many ways, but the simplest and cheapest way is to use
1021 * the first constructing atom ai as a reference position in space:
1022 * subtract (xv-xi)*fv and add (xj-xi)*fj.
1026 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1028 for (int i = 0; i < DIM; i++)
1030 for (int j = 0; j < DIM; j++)
1032 /* As xix is a linear combination of j and k, use that here */
1033 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * fj[j];
1038 /* TOTAL: 38 flops */
1041 template<VirialHandling virialHandling>
1042 static void spread_vsite3(const t_iatom ia[],
1045 ArrayRef<const RVec> x,
1047 ArrayRef<RVec> fshift,
1050 rvec fi, fj, fk, dx;
1058 svmul(1 - a - b, f[av], fi);
1059 svmul(a, f[av], fj);
1060 svmul(b, f[av], fk);
1063 rvec_inc(f[ai], fi);
1064 rvec_inc(f[aj], fj);
1065 rvec_inc(f[ak], fk);
1068 if (virialHandling == VirialHandling::Pbc)
1075 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
1076 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
1077 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
1086 if (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL)
1088 rvec_inc(fshift[siv], f[av]);
1089 rvec_dec(fshift[CENTRAL], fi);
1090 rvec_dec(fshift[sij], fj);
1091 rvec_dec(fshift[sik], fk);
1095 /* TOTAL: 20 flops */
1098 template<VirialHandling virialHandling>
1099 static void spread_vsite3FD(const t_iatom ia[],
1102 ArrayRef<const RVec> x,
1104 ArrayRef<RVec> fshift,
1109 rvec xvi, xij, xjk, xix, fv, temp;
1110 t_iatom av, ai, aj, ak;
1117 copy_rvec(f[av], fv);
1119 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1120 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1123 /* xix goes from i to point x on the line jk */
1124 xix[XX] = xij[XX] + a * xjk[XX];
1125 xix[YY] = xij[YY] + a * xjk[YY];
1126 xix[ZZ] = xij[ZZ] + a * xjk[ZZ];
1129 const real invDistance = inverseNorm(xix);
1130 const real c = b * invDistance;
1131 /* 4 + ?10? flops */
1133 fproj = iprod(xix, fv) * invDistance * invDistance; /* = (xix . f)/(xix . xix) */
1135 temp[XX] = c * (fv[XX] - fproj * xix[XX]);
1136 temp[YY] = c * (fv[YY] - fproj * xix[YY]);
1137 temp[ZZ] = c * (fv[ZZ] - fproj * xix[ZZ]);
1140 /* c is already calculated in constr_vsite3FD
1141 storing c somewhere will save 26 flops! */
1144 f[ai][XX] += fv[XX] - temp[XX];
1145 f[ai][YY] += fv[YY] - temp[YY];
1146 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
1147 f[aj][XX] += a1 * temp[XX];
1148 f[aj][YY] += a1 * temp[YY];
1149 f[aj][ZZ] += a1 * temp[ZZ];
1150 f[ak][XX] += a * temp[XX];
1151 f[ak][YY] += a * temp[YY];
1152 f[ak][ZZ] += a * temp[ZZ];
1155 if (virialHandling == VirialHandling::Pbc)
1160 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1167 if (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL)
1169 rvec_dec(fshift[svi], fv);
1170 fshift[CENTRAL][XX] += fv[XX] - (1 + a) * temp[XX];
1171 fshift[CENTRAL][YY] += fv[YY] - (1 + a) * temp[YY];
1172 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a) * temp[ZZ];
1173 fshift[sji][XX] += temp[XX];
1174 fshift[sji][YY] += temp[YY];
1175 fshift[sji][ZZ] += temp[ZZ];
1176 fshift[skj][XX] += a * temp[XX];
1177 fshift[skj][YY] += a * temp[YY];
1178 fshift[skj][ZZ] += a * temp[ZZ];
1182 if (virialHandling == VirialHandling::NonLinear)
1184 /* Under this condition, the virial for the current forces is not
1185 * calculated from the redistributed forces. This means that
1186 * the effect of non-linear virtual site constructions on the virial
1187 * needs to be added separately. This contribution can be calculated
1188 * in many ways, but the simplest and cheapest way is to use
1189 * the first constructing atom ai as a reference position in space:
1190 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
1194 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1196 for (int i = 0; i < DIM; i++)
1198 for (int j = 0; j < DIM; j++)
1200 /* As xix is a linear combination of j and k, use that here */
1201 dxdf[i][j] += -xiv[i] * fv[j] + xix[i] * temp[j];
1206 /* TOTAL: 61 flops */
1209 template<VirialHandling virialHandling>
1210 static void spread_vsite3FAD(const t_iatom ia[],
1213 ArrayRef<const RVec> x,
1215 ArrayRef<RVec> fshift,
1219 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
1220 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
1221 t_iatom av, ai, aj, ak;
1228 copy_rvec(f[ia[1]], fv);
1230 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1231 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1234 invdij = inverseNorm(xij);
1235 invdij2 = invdij * invdij;
1236 c1 = iprod(xij, xjk) * invdij2;
1237 xperp[XX] = xjk[XX] - c1 * xij[XX];
1238 xperp[YY] = xjk[YY] - c1 * xij[YY];
1239 xperp[ZZ] = xjk[ZZ] - c1 * xij[ZZ];
1240 /* xperp in plane ijk, perp. to ij */
1241 invdp = inverseNorm(xperp);
1246 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
1247 storing them somewhere will save 45 flops! */
1249 fproj = iprod(xij, fv) * invdij2;
1250 svmul(fproj, xij, Fpij); /* proj. f on xij */
1251 svmul(iprod(xperp, fv) * invdp * invdp, xperp, Fppp); /* proj. f on xperp */
1252 svmul(b1 * fproj, xperp, f3);
1255 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
1256 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
1257 for (int d = 0; d < DIM; d++)
1265 f[ai][XX] += fv[XX] - f1[XX] + c1 * f2[XX] + f3[XX];
1266 f[ai][YY] += fv[YY] - f1[YY] + c1 * f2[YY] + f3[YY];
1267 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1 * f2[ZZ] + f3[ZZ];
1268 f[aj][XX] += f1[XX] - c2 * f2[XX] - f3[XX];
1269 f[aj][YY] += f1[YY] - c2 * f2[YY] - f3[YY];
1270 f[aj][ZZ] += f1[ZZ] - c2 * f2[ZZ] - f3[ZZ];
1271 f[ak][XX] += f2[XX];
1272 f[ak][YY] += f2[YY];
1273 f[ak][ZZ] += f2[ZZ];
1276 if (virialHandling == VirialHandling::Pbc)
1282 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1289 if (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL)
1291 rvec_dec(fshift[svi], fv);
1292 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1 - c1) * f2[XX] + f3[XX];
1293 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1 - c1) * f2[YY] + f3[YY];
1294 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1 - c1) * f2[ZZ] + f3[ZZ];
1295 fshift[sji][XX] += f1[XX] - c1 * f2[XX] - f3[XX];
1296 fshift[sji][YY] += f1[YY] - c1 * f2[YY] - f3[YY];
1297 fshift[sji][ZZ] += f1[ZZ] - c1 * f2[ZZ] - f3[ZZ];
1298 fshift[skj][XX] += f2[XX];
1299 fshift[skj][YY] += f2[YY];
1300 fshift[skj][ZZ] += f2[ZZ];
1304 if (virialHandling == VirialHandling::NonLinear)
1307 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1309 for (int i = 0; i < DIM; i++)
1311 for (int j = 0; j < DIM; j++)
1313 /* Note that xik=xij+xjk, so we have to add xij*f2 */
1314 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * (f1[j] + (1 - c2) * f2[j] - f3[j])
1320 /* TOTAL: 113 flops */
1323 template<VirialHandling virialHandling>
1324 static void spread_vsite3OUT(const t_iatom ia[],
1328 ArrayRef<const RVec> x,
1330 ArrayRef<RVec> fshift,
1334 rvec xvi, xij, xik, fv, fj, fk;
1344 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1345 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1348 copy_rvec(f[av], fv);
1355 fj[XX] = a * fv[XX] - xik[ZZ] * cfy + xik[YY] * cfz;
1356 fj[YY] = xik[ZZ] * cfx + a * fv[YY] - xik[XX] * cfz;
1357 fj[ZZ] = -xik[YY] * cfx + xik[XX] * cfy + a * fv[ZZ];
1359 fk[XX] = b * fv[XX] + xij[ZZ] * cfy - xij[YY] * cfz;
1360 fk[YY] = -xij[ZZ] * cfx + b * fv[YY] + xij[XX] * cfz;
1361 fk[ZZ] = xij[YY] * cfx - xij[XX] * cfy + b * fv[ZZ];
1364 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
1365 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
1366 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1367 rvec_inc(f[aj], fj);
1368 rvec_inc(f[ak], fk);
1371 if (virialHandling == VirialHandling::Pbc)
1376 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1383 if (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL)
1385 rvec_dec(fshift[svi], fv);
1386 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
1387 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
1388 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1389 rvec_inc(fshift[sji], fj);
1390 rvec_inc(fshift[ski], fk);
1394 if (virialHandling == VirialHandling::NonLinear)
1398 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1400 for (int i = 0; i < DIM; i++)
1402 for (int j = 0; j < DIM; j++)
1404 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * fj[j] + xik[i] * fk[j];
1409 /* TOTAL: 54 flops */
1412 template<VirialHandling virialHandling>
1413 static void spread_vsite4FD(const t_iatom ia[],
1417 ArrayRef<const RVec> x,
1419 ArrayRef<RVec> fshift,
1424 rvec xvi, xij, xjk, xjl, xix, fv, temp;
1425 int av, ai, aj, ak, al;
1426 int sji, skj, slj, m;
1434 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1435 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1436 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1439 /* xix goes from i to point x on the plane jkl */
1440 for (m = 0; m < DIM; m++)
1442 xix[m] = xij[m] + a * xjk[m] + b * xjl[m];
1446 const real invDistance = inverseNorm(xix);
1447 const real d = c * invDistance;
1448 /* 4 + ?10? flops */
1450 copy_rvec(f[av], fv);
1452 fproj = iprod(xix, fv) * invDistance * invDistance; /* = (xix . f)/(xix . xix) */
1454 for (m = 0; m < DIM; m++)
1456 temp[m] = d * (fv[m] - fproj * xix[m]);
1460 /* c is already calculated in constr_vsite3FD
1461 storing c somewhere will save 35 flops! */
1464 for (m = 0; m < DIM; m++)
1466 f[ai][m] += fv[m] - temp[m];
1467 f[aj][m] += a1 * temp[m];
1468 f[ak][m] += a * temp[m];
1469 f[al][m] += b * temp[m];
1473 if (virialHandling == VirialHandling::Pbc)
1478 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1485 if (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL)
1487 rvec_dec(fshift[svi], fv);
1488 for (m = 0; m < DIM; m++)
1490 fshift[CENTRAL][m] += fv[m] - (1 + a + b) * temp[m];
1491 fshift[sji][m] += temp[m];
1492 fshift[skj][m] += a * temp[m];
1493 fshift[slj][m] += b * temp[m];
1498 if (virialHandling == VirialHandling::NonLinear)
1503 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1505 for (i = 0; i < DIM; i++)
1507 for (j = 0; j < DIM; j++)
1509 dxdf[i][j] += -xiv[i] * fv[j] + xix[i] * temp[j];
1514 /* TOTAL: 77 flops */
1517 template<VirialHandling virialHandling>
1518 static void spread_vsite4FDN(const t_iatom ia[],
1522 ArrayRef<const RVec> x,
1524 ArrayRef<RVec> fshift,
1528 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1529 rvec fv, fj, fk, fl;
1532 int av, ai, aj, ak, al;
1535 /* DEBUG: check atom indices */
1542 copy_rvec(f[av], fv);
1544 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1545 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1546 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1549 ra[XX] = a * xik[XX];
1550 ra[YY] = a * xik[YY];
1551 ra[ZZ] = a * xik[ZZ];
1553 rb[XX] = b * xil[XX];
1554 rb[YY] = b * xil[YY];
1555 rb[ZZ] = b * xil[ZZ];
1559 rvec_sub(ra, xij, rja);
1560 rvec_sub(rb, xij, rjb);
1561 rvec_sub(rb, ra, rab);
1564 cprod(rja, rjb, rm);
1567 invrm = inverseNorm(rm);
1568 denom = invrm * invrm;
1571 cfx = c * invrm * fv[XX];
1572 cfy = c * invrm * fv[YY];
1573 cfz = c * invrm * fv[ZZ];
1584 fj[XX] = (-rm[XX] * rt[XX]) * cfx + (rab[ZZ] - rm[YY] * rt[XX]) * cfy
1585 + (-rab[YY] - rm[ZZ] * rt[XX]) * cfz;
1586 fj[YY] = (-rab[ZZ] - rm[XX] * rt[YY]) * cfx + (-rm[YY] * rt[YY]) * cfy
1587 + (rab[XX] - rm[ZZ] * rt[YY]) * cfz;
1588 fj[ZZ] = (rab[YY] - rm[XX] * rt[ZZ]) * cfx + (-rab[XX] - rm[YY] * rt[ZZ]) * cfy
1589 + (-rm[ZZ] * rt[ZZ]) * cfz;
1595 rt[XX] *= denom * a;
1596 rt[YY] *= denom * a;
1597 rt[ZZ] *= denom * a;
1600 fk[XX] = (-rm[XX] * rt[XX]) * cfx + (-a * rjb[ZZ] - rm[YY] * rt[XX]) * cfy
1601 + (a * rjb[YY] - rm[ZZ] * rt[XX]) * cfz;
1602 fk[YY] = (a * rjb[ZZ] - rm[XX] * rt[YY]) * cfx + (-rm[YY] * rt[YY]) * cfy
1603 + (-a * rjb[XX] - rm[ZZ] * rt[YY]) * cfz;
1604 fk[ZZ] = (-a * rjb[YY] - rm[XX] * rt[ZZ]) * cfx + (a * rjb[XX] - rm[YY] * rt[ZZ]) * cfy
1605 + (-rm[ZZ] * rt[ZZ]) * cfz;
1611 rt[XX] *= denom * b;
1612 rt[YY] *= denom * b;
1613 rt[ZZ] *= denom * b;
1616 fl[XX] = (-rm[XX] * rt[XX]) * cfx + (b * rja[ZZ] - rm[YY] * rt[XX]) * cfy
1617 + (-b * rja[YY] - rm[ZZ] * rt[XX]) * cfz;
1618 fl[YY] = (-b * rja[ZZ] - rm[XX] * rt[YY]) * cfx + (-rm[YY] * rt[YY]) * cfy
1619 + (b * rja[XX] - rm[ZZ] * rt[YY]) * cfz;
1620 fl[ZZ] = (b * rja[YY] - rm[XX] * rt[ZZ]) * cfx + (-b * rja[XX] - rm[YY] * rt[ZZ]) * cfy
1621 + (-rm[ZZ] * rt[ZZ]) * cfz;
1624 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1625 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1626 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1627 rvec_inc(f[aj], fj);
1628 rvec_inc(f[ak], fk);
1629 rvec_inc(f[al], fl);
1632 if (virialHandling == VirialHandling::Pbc)
1637 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1644 if (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL)
1646 rvec_dec(fshift[svi], fv);
1647 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1648 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1649 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1650 rvec_inc(fshift[sij], fj);
1651 rvec_inc(fshift[sik], fk);
1652 rvec_inc(fshift[sil], fl);
1656 if (virialHandling == VirialHandling::NonLinear)
1661 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1663 for (i = 0; i < DIM; i++)
1665 for (j = 0; j < DIM; j++)
1667 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * fj[j] + xik[i] * fk[j] + xil[i] * fl[j];
1672 /* Total: 207 flops (Yuck!) */
1675 template<VirialHandling virialHandling>
1676 static int spread_vsiten(const t_iatom ia[],
1677 ArrayRef<const t_iparams> ip,
1678 ArrayRef<const RVec> x,
1680 ArrayRef<RVec> fshift,
1688 n3 = 3 * ip[ia[0]].vsiten.n;
1690 copy_rvec(x[av], xv);
1692 for (i = 0; i < n3; i += 3)
1697 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1703 a = ip[ia[i]].vsiten.a;
1704 svmul(a, f[av], fi);
1705 rvec_inc(f[ai], fi);
1707 if (virialHandling == VirialHandling::Pbc && siv != CENTRAL)
1709 rvec_inc(fshift[siv], fi);
1710 rvec_dec(fshift[CENTRAL], fi);
1720 //! Returns the number of virtual sites in the interaction list, for VSITEN the number of atoms
1721 static int vsite_count(ArrayRef<const InteractionList> ilist, int ftype)
1723 if (ftype == F_VSITEN)
1725 return ilist[ftype].size() / 3;
1729 return ilist[ftype].size() / (1 + interaction_function[ftype].nratoms);
1733 //! Executes the force spreading task for a single thread
1734 template<VirialHandling virialHandling>
1735 static void spreadForceForThread(ArrayRef<const RVec> x,
1737 ArrayRef<RVec> fshift,
1739 ArrayRef<const t_iparams> ip,
1740 ArrayRef<const InteractionList> ilist,
1741 const t_pbc* pbc_null)
1743 const PbcMode pbcMode = getPbcMode(pbc_null);
1744 /* We need another pbc pointer, as with charge groups we switch per vsite */
1745 const t_pbc* pbc_null2 = pbc_null;
1746 gmx::ArrayRef<const int> vsite_pbc;
1748 /* this loop goes backwards to be able to build *
1749 * higher type vsites from lower types */
1750 for (int ftype = c_ftypeVsiteEnd - 1; ftype >= c_ftypeVsiteStart; ftype--)
1752 if (ilist[ftype].empty())
1758 int nra = interaction_function[ftype].nratoms;
1760 int nr = ilist[ftype].size();
1762 const t_iatom* ia = ilist[ftype].iatoms.data();
1764 if (pbcMode == PbcMode::all)
1766 pbc_null2 = pbc_null;
1769 for (int i = 0; i < nr;)
1773 /* Constants for constructing */
1775 a1 = ip[tp].vsite.a;
1776 /* Construct the vsite depending on type */
1779 case F_VSITE1: spread_vsite1(ia, f); break;
1781 spread_vsite2<virialHandling>(ia, a1, x, f, fshift, pbc_null2);
1784 spread_vsite2FD<virialHandling>(ia, a1, x, f, fshift, dxdf, pbc_null2);
1787 b1 = ip[tp].vsite.b;
1788 spread_vsite3<virialHandling>(ia, a1, b1, x, f, fshift, pbc_null2);
1791 b1 = ip[tp].vsite.b;
1792 spread_vsite3FD<virialHandling>(ia, a1, b1, x, f, fshift, dxdf, pbc_null2);
1795 b1 = ip[tp].vsite.b;
1796 spread_vsite3FAD<virialHandling>(ia, a1, b1, x, f, fshift, dxdf, pbc_null2);
1799 b1 = ip[tp].vsite.b;
1800 c1 = ip[tp].vsite.c;
1801 spread_vsite3OUT<virialHandling>(ia, a1, b1, c1, x, f, fshift, dxdf, pbc_null2);
1804 b1 = ip[tp].vsite.b;
1805 c1 = ip[tp].vsite.c;
1806 spread_vsite4FD<virialHandling>(ia, a1, b1, c1, x, f, fshift, dxdf, pbc_null2);
1809 b1 = ip[tp].vsite.b;
1810 c1 = ip[tp].vsite.c;
1811 spread_vsite4FDN<virialHandling>(ia, a1, b1, c1, x, f, fshift, dxdf, pbc_null2);
1814 inc = spread_vsiten<virialHandling>(ia, ip, x, f, fshift, pbc_null2);
1817 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d", ftype, __FILE__, __LINE__);
1819 clear_rvec(f[ia[1]]);
1821 /* Increment loop variables */
1829 //! Wrapper function for calling the templated thread-local spread function
1830 static void spreadForceWrapper(ArrayRef<const RVec> x,
1832 const VirialHandling virialHandling,
1833 ArrayRef<RVec> fshift,
1835 const bool clearDxdf,
1836 ArrayRef<const t_iparams> ip,
1837 ArrayRef<const InteractionList> ilist,
1838 const t_pbc* pbc_null)
1840 if (virialHandling == VirialHandling::NonLinear && clearDxdf)
1845 switch (virialHandling)
1847 case VirialHandling::None:
1848 spreadForceForThread<VirialHandling::None>(x, f, fshift, dxdf, ip, ilist, pbc_null);
1850 case VirialHandling::Pbc:
1851 spreadForceForThread<VirialHandling::Pbc>(x, f, fshift, dxdf, ip, ilist, pbc_null);
1853 case VirialHandling::NonLinear:
1854 spreadForceForThread<VirialHandling::NonLinear>(x, f, fshift, dxdf, ip, ilist, pbc_null);
1859 //! Clears the task force buffer elements that are written by task idTask
1860 static void clearTaskForceBufferUsedElements(InterdependentTask* idTask)
1862 int ntask = idTask->spreadTask.size();
1863 for (int ti = 0; ti < ntask; ti++)
1865 const AtomIndex* atomList = &idTask->atomIndex[idTask->spreadTask[ti]];
1866 int natom = atomList->atom.size();
1867 RVec* force = idTask->force.data();
1868 for (int i = 0; i < natom; i++)
1870 clear_rvec(force[atomList->atom[i]]);
1875 void VirtualSitesHandler::Impl::spreadForces(ArrayRef<const RVec> x,
1877 const VirialHandling virialHandling,
1878 ArrayRef<RVec> fshift,
1882 gmx_wallcycle* wcycle)
1884 wallcycle_start(wcycle, ewcVSITESPREAD);
1886 const bool useDomdec = domainInfo_.useDomdec();
1888 t_pbc pbc, *pbc_null;
1890 if (domainInfo_.useMolPbc_)
1892 /* This is wasting some CPU time as we now do this multiple times
1895 pbc_null = set_pbc_dd(&pbc, domainInfo_.pbcType_,
1896 useDomdec ? domainInfo_.domdec_->numCells : nullptr, FALSE, box);
1905 dd_clear_f_vsites(*domainInfo_.domdec_, f);
1908 const int numThreads = threadingInfo_.numThreads();
1910 if (numThreads == 1)
1913 spreadForceWrapper(x, f, virialHandling, fshift, dxdf, true, iparams_, ilists_, pbc_null);
1915 if (virialHandling == VirialHandling::NonLinear)
1917 for (int i = 0; i < DIM; i++)
1919 for (int j = 0; j < DIM; j++)
1921 virial[i][j] += -0.5 * dxdf[i][j];
1928 /* First spread the vsites that might depend on non-local vsites */
1929 auto& nlDependentVSites = threadingInfo_.threadDataNonLocalDependent();
1930 spreadForceWrapper(x, f, virialHandling, fshift, nlDependentVSites.dxdf, true, iparams_,
1931 nlDependentVSites.ilist, pbc_null);
1933 #pragma omp parallel num_threads(numThreads)
1937 int thread = gmx_omp_get_thread_num();
1938 VsiteThread& tData = threadingInfo_.threadData(thread);
1940 ArrayRef<RVec> fshift_t;
1941 if (virialHandling == VirialHandling::Pbc)
1949 fshift_t = tData.fshift;
1951 for (int i = 0; i < SHIFTS; i++)
1953 clear_rvec(fshift_t[i]);
1958 if (tData.useInterdependentTask)
1960 /* Spread the vsites that spread outside our local range.
1961 * This is done using a thread-local force buffer force.
1962 * First we need to copy the input vsite forces to force.
1964 InterdependentTask* idTask = &tData.idTask;
1966 /* Clear the buffer elements set by our task during
1967 * the last call to spread_vsite_f.
1969 clearTaskForceBufferUsedElements(idTask);
1971 int nvsite = idTask->vsite.size();
1972 for (int i = 0; i < nvsite; i++)
1974 copy_rvec(f[idTask->vsite[i]], idTask->force[idTask->vsite[i]]);
1976 spreadForceWrapper(x, idTask->force, virialHandling, fshift_t, tData.dxdf, true,
1977 iparams_, tData.idTask.ilist, pbc_null);
1979 /* We need a barrier before reducing forces below
1980 * that have been produced by a different thread above.
1984 /* Loop over all thread task and reduce forces they
1985 * produced on atoms that fall in our range.
1986 * Note that atomic reduction would be a simpler solution,
1987 * but that might not have good support on all platforms.
1989 int ntask = idTask->reduceTask.size();
1990 for (int ti = 0; ti < ntask; ti++)
1992 const InterdependentTask& idt_foreign =
1993 threadingInfo_.threadData(idTask->reduceTask[ti]).idTask;
1994 const AtomIndex& atomList = idt_foreign.atomIndex[thread];
1995 const RVec* f_foreign = idt_foreign.force.data();
1997 for (int ind : atomList.atom)
1999 rvec_inc(f[ind], f_foreign[ind]);
2000 /* Clearing of f_foreign is done at the next step */
2003 /* Clear the vsite forces, both in f and force */
2004 for (int i = 0; i < nvsite; i++)
2006 int ind = tData.idTask.vsite[i];
2008 clear_rvec(tData.idTask.force[ind]);
2012 /* Spread the vsites that spread locally only */
2013 spreadForceWrapper(x, f, virialHandling, fshift_t, tData.dxdf, false, iparams_,
2014 tData.ilist, pbc_null);
2016 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2019 if (virialHandling == VirialHandling::Pbc)
2021 for (int th = 1; th < numThreads; th++)
2023 for (int i = 0; i < SHIFTS; i++)
2025 rvec_inc(fshift[i], threadingInfo_.threadData(th).fshift[i]);
2030 if (virialHandling == VirialHandling::NonLinear)
2032 for (int th = 0; th < numThreads + 1; th++)
2034 /* MSVC doesn't like matrix references, so we use a pointer */
2035 const matrix& dxdf = threadingInfo_.threadData(th).dxdf;
2037 for (int i = 0; i < DIM; i++)
2039 for (int j = 0; j < DIM; j++)
2041 virial[i][j] += -0.5 * dxdf[i][j];
2050 dd_move_f_vsites(*domainInfo_.domdec_, f, fshift);
2053 inc_nrnb(nrnb, eNR_VSITE1, vsite_count(ilists_, F_VSITE1));
2054 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(ilists_, F_VSITE2));
2055 inc_nrnb(nrnb, eNR_VSITE2FD, vsite_count(ilists_, F_VSITE2FD));
2056 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(ilists_, F_VSITE3));
2057 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(ilists_, F_VSITE3FD));
2058 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(ilists_, F_VSITE3FAD));
2059 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(ilists_, F_VSITE3OUT));
2060 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(ilists_, F_VSITE4FD));
2061 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(ilists_, F_VSITE4FDN));
2062 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(ilists_, F_VSITEN));
2064 wallcycle_stop(wcycle, ewcVSITESPREAD);
2067 /*! \brief Returns the an array with group indices for each atom
2069 * \param[in] grouping The paritioning of the atom range into atom groups
2071 static std::vector<int> makeAtomToGroupMapping(const gmx::RangePartitioning& grouping)
2073 std::vector<int> atomToGroup(grouping.fullRange().end(), 0);
2075 for (int group = 0; group < grouping.numBlocks(); group++)
2077 auto block = grouping.block(group);
2078 std::fill(atomToGroup.begin() + block.begin(), atomToGroup.begin() + block.end(), group);
2084 int countNonlinearVsites(const gmx_mtop_t& mtop)
2086 int numNonlinearVsites = 0;
2087 for (const gmx_molblock_t& molb : mtop.molblock)
2089 const gmx_moltype_t& molt = mtop.moltype[molb.type];
2091 for (const auto& ilist : extractILists(molt.ilist, IF_VSITE))
2093 if (ilist.functionType != F_VSITE2 && ilist.functionType != F_VSITE3
2094 && ilist.functionType != F_VSITEN)
2096 numNonlinearVsites += molb.nmol * ilist.iatoms.size() / (1 + NRAL(ilist.functionType));
2101 return numNonlinearVsites;
2104 void VirtualSitesHandler::spreadForces(ArrayRef<const RVec> x,
2106 const VirialHandling virialHandling,
2107 ArrayRef<RVec> fshift,
2111 gmx_wallcycle* wcycle)
2113 impl_->spreadForces(x, f, virialHandling, fshift, virial, nrnb, box, wcycle);
2116 int countInterUpdategroupVsites(const gmx_mtop_t& mtop,
2117 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingPerMoleculetype)
2119 int n_intercg_vsite = 0;
2120 for (const gmx_molblock_t& molb : mtop.molblock)
2122 const gmx_moltype_t& molt = mtop.moltype[molb.type];
2124 std::vector<int> atomToGroup;
2125 if (!updateGroupingPerMoleculetype.empty())
2127 atomToGroup = makeAtomToGroupMapping(updateGroupingPerMoleculetype[molb.type]);
2129 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2131 const int nral = NRAL(ftype);
2132 const InteractionList& il = molt.ilist[ftype];
2133 for (int i = 0; i < il.size(); i += 1 + nral)
2135 bool isInterGroup = atomToGroup.empty();
2138 const int group = atomToGroup[il.iatoms[1 + i]];
2139 for (int a = 1; a < nral; a++)
2141 if (atomToGroup[il.iatoms[1 + a]] != group)
2143 isInterGroup = true;
2150 n_intercg_vsite += molb.nmol;
2156 return n_intercg_vsite;
2159 std::unique_ptr<VirtualSitesHandler> makeVirtualSitesHandler(const gmx_mtop_t& mtop,
2160 const t_commrec* cr,
2163 GMX_RELEASE_ASSERT(cr != nullptr, "We need a valid commrec");
2165 std::unique_ptr<VirtualSitesHandler> vsite;
2167 /* check if there are vsites */
2169 for (int ftype = 0; ftype < F_NRE; ftype++)
2171 if (interaction_function[ftype].flags & IF_VSITE)
2173 GMX_ASSERT(ftype >= c_ftypeVsiteStart && ftype < c_ftypeVsiteEnd,
2174 "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
2176 nvsite += gmx_mtop_ftype_count(&mtop, ftype);
2180 GMX_ASSERT(ftype < c_ftypeVsiteStart || ftype >= c_ftypeVsiteEnd,
2181 "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
2190 return std::make_unique<VirtualSitesHandler>(mtop, cr->dd, pbcType);
2193 ThreadingInfo::ThreadingInfo() : numThreads_(gmx_omp_nthreads_get(emntVSITE))
2195 if (numThreads_ > 1)
2197 /* We need one extra thread data structure for the overlap vsites */
2198 tData_.resize(numThreads_ + 1);
2199 #pragma omp parallel for num_threads(numThreads_) schedule(static)
2200 for (int thread = 0; thread < numThreads_; thread++)
2204 tData_[thread] = std::make_unique<VsiteThread>();
2206 InterdependentTask& idTask = tData_[thread]->idTask;
2208 idTask.atomIndex.resize(numThreads_);
2210 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2212 if (numThreads_ > 1)
2214 tData_[numThreads_] = std::make_unique<VsiteThread>();
2219 //! Returns the number of inter update-group vsites
2220 static int getNumInterUpdategroupVsites(const gmx_mtop_t& mtop, const gmx_domdec_t* domdec)
2222 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingPerMoleculetype;
2225 updateGroupingPerMoleculetype = getUpdateGroupingPerMoleculetype(*domdec);
2228 return countInterUpdategroupVsites(mtop, updateGroupingPerMoleculetype);
2231 VirtualSitesHandler::Impl::Impl(const gmx_mtop_t& mtop, gmx_domdec_t* domdec, const PbcType pbcType) :
2232 numInterUpdategroupVirtualSites_(getNumInterUpdategroupVsites(mtop, domdec)),
2233 domainInfo_({ pbcType, pbcType != PbcType::No && numInterUpdategroupVirtualSites_ > 0, domdec }),
2234 iparams_(mtop.ffparams.iparams)
2238 VirtualSitesHandler::VirtualSitesHandler(const gmx_mtop_t& mtop, gmx_domdec_t* domdec, const PbcType pbcType) :
2239 impl_(new Impl(mtop, domdec, pbcType))
2243 //! Flag that atom \p atom which is home in another task, if it has not already been added before
2244 static inline void flagAtom(InterdependentTask* idTask, const int atom, const int numThreads, const int numAtomsPerThread)
2246 if (!idTask->use[atom])
2248 idTask->use[atom] = true;
2249 int thread = atom / numAtomsPerThread;
2250 /* Assign all non-local atom force writes to thread 0 */
2251 if (thread >= numThreads)
2255 idTask->atomIndex[thread].atom.push_back(atom);
2259 /*! \brief Here we try to assign all vsites that are in our local range.
2261 * Our task local atom range is tData->rangeStart - tData->rangeEnd.
2262 * Vsites that depend only on local atoms, as indicated by taskIndex[]==thread,
2263 * are assigned to task tData->ilist. Vsites that depend on non-local atoms
2264 * but not on other vsites are assigned to task tData->id_task.ilist.
2265 * taskIndex[] is set for all vsites in our range, either to our local tasks
2266 * or to the single last task as taskIndex[]=2*nthreads.
2268 static void assignVsitesToThread(VsiteThread* tData,
2272 gmx::ArrayRef<int> taskIndex,
2273 ArrayRef<const InteractionList> ilist,
2274 ArrayRef<const t_iparams> ip,
2275 const unsigned short* ptype)
2277 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2279 tData->ilist[ftype].clear();
2280 tData->idTask.ilist[ftype].clear();
2282 int nral1 = 1 + NRAL(ftype);
2284 const int* iat = ilist[ftype].iatoms.data();
2285 for (int i = 0; i < ilist[ftype].size();)
2287 if (ftype == F_VSITEN)
2289 /* The 3 below is from 1+NRAL(ftype)=3 */
2290 inc = ip[iat[i]].vsiten.n * 3;
2293 if (iat[1 + i] < tData->rangeStart || iat[1 + i] >= tData->rangeEnd)
2295 /* This vsite belongs to a different thread */
2300 /* We would like to assign this vsite to task thread,
2301 * but it might depend on atoms outside the atom range of thread
2302 * or on another vsite not assigned to task thread.
2305 if (ftype != F_VSITEN)
2307 for (int j = i + 2; j < i + nral1; j++)
2309 /* Do a range check to avoid a harmless race on taskIndex */
2310 if (iat[j] < tData->rangeStart || iat[j] >= tData->rangeEnd || taskIndex[iat[j]] != thread)
2312 if (!tData->useInterdependentTask || ptype[iat[j]] == eptVSite)
2314 /* At least one constructing atom is a vsite
2315 * that is not assigned to the same thread.
2316 * Put this vsite into a separate task.
2322 /* There are constructing atoms outside our range,
2323 * put this vsite into a second task to be executed
2324 * on the same thread. During construction no barrier
2325 * is needed between the two tasks on the same thread.
2326 * During spreading we need to run this task with
2327 * an additional thread-local intermediate force buffer
2328 * (or atomic reduction) and a barrier between the two
2331 task = nthread + thread;
2337 for (int j = i + 2; j < i + inc; j += 3)
2339 /* Do a range check to avoid a harmless race on taskIndex */
2340 if (iat[j] < tData->rangeStart || iat[j] >= tData->rangeEnd || taskIndex[iat[j]] != thread)
2342 GMX_ASSERT(ptype[iat[j]] != eptVSite,
2343 "A vsite to be assigned in assignVsitesToThread has a vsite as "
2344 "a constructing atom that does not belong to our task, such "
2345 "vsites should be assigned to the single 'master' task");
2347 task = nthread + thread;
2352 /* Update this vsite's thread index entry */
2353 taskIndex[iat[1 + i]] = task;
2355 if (task == thread || task == nthread + thread)
2357 /* Copy this vsite to the thread data struct of thread */
2358 InteractionList* il_task;
2361 il_task = &tData->ilist[ftype];
2365 il_task = &tData->idTask.ilist[ftype];
2367 /* Copy the vsite data to the thread-task local array */
2368 il_task->push_back(iat[i], nral1 - 1, iat + i + 1);
2369 if (task == nthread + thread)
2371 /* This vsite write outside our own task force block.
2372 * Put it into the interdependent task list and flag
2373 * the atoms involved for reduction.
2375 tData->idTask.vsite.push_back(iat[i + 1]);
2376 if (ftype != F_VSITEN)
2378 for (int j = i + 2; j < i + nral1; j++)
2380 flagAtom(&tData->idTask, iat[j], nthread, natperthread);
2385 for (int j = i + 2; j < i + inc; j += 3)
2387 flagAtom(&tData->idTask, iat[j], nthread, natperthread);
2398 /*! \brief Assign all vsites with taskIndex[]==task to task tData */
2399 static void assignVsitesToSingleTask(VsiteThread* tData,
2401 gmx::ArrayRef<const int> taskIndex,
2402 ArrayRef<const InteractionList> ilist,
2403 ArrayRef<const t_iparams> ip)
2405 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2407 tData->ilist[ftype].clear();
2408 tData->idTask.ilist[ftype].clear();
2410 int nral1 = 1 + NRAL(ftype);
2412 const int* iat = ilist[ftype].iatoms.data();
2413 InteractionList* il_task = &tData->ilist[ftype];
2415 for (int i = 0; i < ilist[ftype].size();)
2417 if (ftype == F_VSITEN)
2419 /* The 3 below is from 1+NRAL(ftype)=3 */
2420 inc = ip[iat[i]].vsiten.n * 3;
2422 /* Check if the vsite is assigned to our task */
2423 if (taskIndex[iat[1 + i]] == task)
2425 /* Copy the vsite data to the thread-task local array */
2426 il_task->push_back(iat[i], inc - 1, iat + i + 1);
2434 void ThreadingInfo::setVirtualSites(ArrayRef<const InteractionList> ilists,
2435 ArrayRef<const t_iparams> iparams,
2436 const t_mdatoms& mdatoms,
2437 const bool useDomdec)
2439 if (numThreads_ <= 1)
2445 /* The current way of distributing the vsites over threads in primitive.
2446 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
2447 * without taking into account how the vsites are distributed.
2448 * Without domain decomposition we at least tighten the upper bound
2449 * of the range (useful for common systems such as a vsite-protein
2451 * With domain decomposition, as long as the vsites are distributed
2452 * uniformly in each domain along the major dimension, usually x,
2453 * it will also perform well.
2455 int vsite_atom_range;
2459 vsite_atom_range = -1;
2460 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2463 if (ftype != F_VSITEN)
2465 int nral1 = 1 + NRAL(ftype);
2466 ArrayRef<const int> iat = ilists[ftype].iatoms;
2467 for (int i = 0; i < ilists[ftype].size(); i += nral1)
2469 for (int j = i + 1; j < i + nral1; j++)
2471 vsite_atom_range = std::max(vsite_atom_range, iat[j]);
2479 ArrayRef<const int> iat = ilists[ftype].iatoms;
2482 while (i < ilists[ftype].size())
2484 /* The 3 below is from 1+NRAL(ftype)=3 */
2485 vs_ind_end = i + iparams[iat[i]].vsiten.n * 3;
2487 vsite_atom_range = std::max(vsite_atom_range, iat[i + 1]);
2488 while (i < vs_ind_end)
2490 vsite_atom_range = std::max(vsite_atom_range, iat[i + 2]);
2498 natperthread = (vsite_atom_range + numThreads_ - 1) / numThreads_;
2502 /* Any local or not local atom could be involved in virtual sites.
2503 * But since we usually have very few non-local virtual sites
2504 * (only non-local vsites that depend on local vsites),
2505 * we distribute the local atom range equally over the threads.
2506 * When assigning vsites to threads, we should take care that the last
2507 * threads also covers the non-local range.
2509 vsite_atom_range = mdatoms.nr;
2510 natperthread = (mdatoms.homenr + numThreads_ - 1) / numThreads_;
2515 fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n",
2516 mdatoms.nr, vsite_atom_range, natperthread);
2519 /* To simplify the vsite assignment, we make an index which tells us
2520 * to which task particles, both non-vsites and vsites, are assigned.
2522 taskIndex_.resize(mdatoms.nr);
2524 /* Initialize the task index array. Here we assign the non-vsite
2525 * particles to task=thread, so we easily figure out if vsites
2526 * depend on local and/or non-local particles in assignVsitesToThread.
2530 for (int i = 0; i < mdatoms.nr; i++)
2532 if (mdatoms.ptype[i] == eptVSite)
2534 /* vsites are not assigned to a task yet */
2539 /* assign non-vsite particles to task thread */
2540 taskIndex_[i] = thread;
2542 if (i == (thread + 1) * natperthread && thread < numThreads_)
2549 #pragma omp parallel num_threads(numThreads_)
2553 int thread = gmx_omp_get_thread_num();
2554 VsiteThread& tData = *tData_[thread];
2556 /* Clear the buffer use flags that were set before */
2557 if (tData.useInterdependentTask)
2559 InterdependentTask& idTask = tData.idTask;
2561 /* To avoid an extra OpenMP barrier in spread_vsite_f,
2562 * we clear the force buffer at the next step,
2563 * so we need to do it here as well.
2565 clearTaskForceBufferUsedElements(&idTask);
2567 idTask.vsite.resize(0);
2568 for (int t = 0; t < numThreads_; t++)
2570 AtomIndex& atomIndex = idTask.atomIndex[t];
2571 int natom = atomIndex.atom.size();
2572 for (int i = 0; i < natom; i++)
2574 idTask.use[atomIndex.atom[i]] = false;
2576 atomIndex.atom.resize(0);
2581 /* To avoid large f_buf allocations of #threads*vsite_atom_range
2582 * we don't use task2 with more than 200000 atoms. This doesn't
2583 * affect performance, since with such a large range relatively few
2584 * vsites will end up in the separate task.
2585 * Note that useTask2 should be the same for all threads.
2587 tData.useInterdependentTask = (vsite_atom_range <= 200000);
2588 if (tData.useInterdependentTask)
2590 size_t natoms_use_in_vsites = vsite_atom_range;
2591 InterdependentTask& idTask = tData.idTask;
2592 /* To avoid resizing and re-clearing every nstlist steps,
2593 * we never down size the force buffer.
2595 if (natoms_use_in_vsites > idTask.force.size() || natoms_use_in_vsites > idTask.use.size())
2597 idTask.force.resize(natoms_use_in_vsites, { 0, 0, 0 });
2598 idTask.use.resize(natoms_use_in_vsites, false);
2602 /* Assign all vsites that can execute independently on threads */
2603 tData.rangeStart = thread * natperthread;
2604 if (thread < numThreads_ - 1)
2606 tData.rangeEnd = (thread + 1) * natperthread;
2610 /* The last thread should cover up to the end of the range */
2611 tData.rangeEnd = mdatoms.nr;
2613 assignVsitesToThread(&tData, thread, numThreads_, natperthread, taskIndex_, ilists,
2614 iparams, mdatoms.ptype);
2616 if (tData.useInterdependentTask)
2618 /* In the worst case, all tasks write to force ranges of
2619 * all other tasks, leading to #tasks^2 scaling (this is only
2620 * the overhead, the actual flops remain constant).
2621 * But in most cases there is far less coupling. To improve
2622 * scaling at high thread counts we therefore construct
2623 * an index to only loop over the actually affected tasks.
2625 InterdependentTask& idTask = tData.idTask;
2627 /* Ensure assignVsitesToThread finished on other threads */
2630 idTask.spreadTask.resize(0);
2631 idTask.reduceTask.resize(0);
2632 for (int t = 0; t < numThreads_; t++)
2634 /* Do we write to the force buffer of task t? */
2635 if (!idTask.atomIndex[t].atom.empty())
2637 idTask.spreadTask.push_back(t);
2639 /* Does task t write to our force buffer? */
2640 if (!tData_[t]->idTask.atomIndex[thread].atom.empty())
2642 idTask.reduceTask.push_back(t);
2647 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2649 /* Assign all remaining vsites, that will have taskIndex[]=2*vsite->nthreads,
2650 * to a single task that will not run in parallel with other tasks.
2652 assignVsitesToSingleTask(tData_[numThreads_].get(), 2 * numThreads_, taskIndex_, ilists, iparams);
2654 if (debug && numThreads_ > 1)
2656 fprintf(debug, "virtual site useInterdependentTask %d, nuse:\n",
2657 static_cast<int>(tData_[0]->useInterdependentTask));
2658 for (int th = 0; th < numThreads_ + 1; th++)
2660 fprintf(debug, " %4d", tData_[th]->idTask.nuse);
2662 fprintf(debug, "\n");
2664 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2666 if (!ilists[ftype].empty())
2668 fprintf(debug, "%-20s thread dist:", interaction_function[ftype].longname);
2669 for (int th = 0; th < numThreads_ + 1; th++)
2671 fprintf(debug, " %4d %4d ", tData_[th]->ilist[ftype].size(),
2672 tData_[th]->idTask.ilist[ftype].size());
2674 fprintf(debug, "\n");
2680 int nrOrig = vsiteIlistNrCount(ilists);
2682 for (int th = 0; th < numThreads_ + 1; th++)
2684 nrThreaded += vsiteIlistNrCount(tData_[th]->ilist) + vsiteIlistNrCount(tData_[th]->idTask.ilist);
2686 GMX_ASSERT(nrThreaded == nrOrig,
2687 "The number of virtual sites assigned to all thread task has to match the total "
2688 "number of virtual sites");
2692 void VirtualSitesHandler::Impl::setVirtualSites(ArrayRef<const InteractionList> ilists,
2693 const t_mdatoms& mdatoms)
2697 threadingInfo_.setVirtualSites(ilists, iparams_, mdatoms, domainInfo_.useDomdec());
2700 void VirtualSitesHandler::setVirtualSites(ArrayRef<const InteractionList> ilists, const t_mdatoms& mdatoms)
2702 impl_->setVirtualSites(ilists, mdatoms);