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.
38 /*! \libinternal \file
39 * \brief Implements the VirtualSitesHandler class and vsite standalone functions
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_mdlib
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_struct.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/gmxlib/nrnb.h"
60 #include "gromacs/math/functions.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/mdlib/gmx_omp_nthreads.h"
63 #include "gromacs/mdtypes/commrec.h"
64 #include "gromacs/mdtypes/mdatom.h"
65 #include "gromacs/pbcutil/ishift.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/timing/wallcycle.h"
68 #include "gromacs/topology/ifunc.h"
69 #include "gromacs/topology/mtop_util.h"
70 #include "gromacs/topology/topology.h"
71 #include "gromacs/utility/exceptions.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/gmxassert.h"
74 #include "gromacs/utility/gmxomp.h"
76 /* The strategy used here for assigning virtual sites to (thread-)tasks
79 * We divide the atom range that vsites operate on (natoms_local with DD,
80 * 0 - last atom involved in vsites without DD) equally over all threads.
82 * Vsites in the local range constructed from atoms in the local range
83 * and/or other vsites that are fully local are assigned to a simple,
86 * Vsites that are not assigned after using the above criterion get assigned
87 * to a so called "interdependent" thread task when none of the constructing
88 * atoms is a vsite. These tasks are called interdependent, because one task
89 * accesses atoms assigned to a different task/thread.
90 * Note that this option is turned off with large (local) atom counts
91 * to avoid high memory usage.
93 * Any remaining vsites are assigned to a separate master thread task.
99 //! VirialHandling is often used outside VirtualSitesHandler class members
100 using VirialHandling = VirtualSitesHandler::VirialHandling;
103 * \brief Information on PBC and domain decomposition for virtual sites
108 //! Constructs without PBC and DD
109 DomainInfo() = default;
111 //! Constructs with PBC and DD, if !=nullptr
112 DomainInfo(PbcType pbcType, bool haveInterUpdateGroupVirtualSites, gmx_domdec_t* domdec) :
114 useMolPbc_(pbcType != PbcType::No && haveInterUpdateGroupVirtualSites),
119 //! Returns whether we are using domain decomposition with more than 1 DD rank
120 bool useDomdec() const { return (domdec_ != nullptr); }
123 const PbcType pbcType_ = PbcType::No;
124 //! Whether molecules are broken over PBC
125 const bool useMolPbc_ = false;
126 //! Pointer to the domain decomposition struct, nullptr without PP DD
127 const gmx_domdec_t* domdec_ = nullptr;
131 * \brief List of atom indices belonging to a task
135 //! List of atom indices
136 std::vector<int> atom;
140 * \brief Data structure for thread tasks that use constructing atoms outside their own atom range
142 struct InterdependentTask
144 //! The interaction lists, only vsite entries are used
145 InteractionLists ilist;
146 //! Thread/task-local force buffer
147 std::vector<RVec> force;
148 //! The atom indices of the vsites of our task
149 std::vector<int> vsite;
150 //! Flags if elements in force are spread to or not
151 std::vector<bool> use;
152 //! The number of entries set to true in use
154 //! Array of atoms indices, size nthreads, covering all nuse set elements in use
155 std::vector<AtomIndex> atomIndex;
156 //! List of tasks (force blocks) this task spread forces to
157 std::vector<int> spreadTask;
158 //! List of tasks that write to this tasks force block range
159 std::vector<int> reduceTask;
163 * \brief Vsite thread task data structure
167 //! Start of atom range of this task
169 //! End of atom range of this task
171 //! The interaction lists, only vsite entries are used
172 std::array<InteractionList, F_NRE> ilist;
173 //! Local fshift accumulation buffer
174 std::array<RVec, SHIFTS> fshift;
175 //! Local virial dx*df accumulation buffer
177 //! 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
178 bool useInterdependentTask;
179 //! Data for vsites that involve constructing atoms in the atom range of other threads/tasks
180 InterdependentTask idTask;
182 /*! \brief Constructor */
187 for (auto& elem : fshift)
189 elem = { 0.0_real, 0.0_real, 0.0_real };
192 useInterdependentTask = false;
198 * \brief Information on how the virtual site work is divided over thread tasks
203 //! Constructor, retrieves the number of threads to use from gmx_omp_nthreads.h
206 //! Returns the number of threads to use for vsite operations
207 int numThreads() const { return numThreads_; }
209 //! Returns the thread data for the given thread
210 const VsiteThread& threadData(int threadIndex) const { return *tData_[threadIndex]; }
212 //! Returns the thread data for the given thread
213 VsiteThread& threadData(int threadIndex) { return *tData_[threadIndex]; }
215 //! Returns the thread data for vsites that depend on non-local vsites
216 const VsiteThread& threadDataNonLocalDependent() const { return *tData_[numThreads_]; }
218 //! Returns the thread data for vsites that depend on non-local vsites
219 VsiteThread& threadDataNonLocalDependent() { return *tData_[numThreads_]; }
221 //! Set VSites and distribute VSite work over threads, should be called after DD partitioning
222 void setVirtualSites(ArrayRef<const InteractionList> ilist,
223 ArrayRef<const t_iparams> iparams,
224 const t_mdatoms& mdatoms,
228 //! Number of threads used for vsite operations
229 const int numThreads_;
230 //! Thread local vsites and work structs
231 std::vector<std::unique_ptr<VsiteThread>> tData_;
232 //! Work array for dividing vsites over threads
233 std::vector<int> taskIndex_;
237 * \brief Impl class for VirtualSitesHandler
239 class VirtualSitesHandler::Impl
242 //! Constructor, domdec should be nullptr without DD
243 Impl(const gmx_mtop_t& mtop, gmx_domdec_t* domdec, PbcType pbcType);
245 //! Returns the number of virtual sites acting over multiple update groups
246 int numInterUpdategroupVirtualSites() const { return numInterUpdategroupVirtualSites_; }
248 //! Set VSites and distribute VSite work over threads, should be called after DD partitioning
249 void setVirtualSites(ArrayRef<const InteractionList> ilist, const t_mdatoms& mdatoms);
251 /*! \brief Create positions of vsite atoms based for the local system
253 * \param[in,out] x The coordinates
254 * \param[in] dt The time step
255 * \param[in,out] v When != nullptr, velocities for vsites are set as displacement/dt
256 * \param[in] box The box
258 void construct(ArrayRef<RVec> x, real dt, ArrayRef<RVec> v, const matrix box) const;
260 /*! \brief Spread the force operating on the vsite atoms on the surrounding atoms.
262 * vsite should point to a valid object.
263 * The virialHandling parameter determines how virial contributions are handled.
264 * If this is set to Linear, shift forces are accumulated into fshift.
265 * If this is set to NonLinear, non-linear contributions are added to virial.
266 * This non-linear correction is required when the virial is not calculated
267 * afterwards from the particle position and forces, but in a different way,
268 * as for instance for the PME mesh contribution.
270 void spreadForces(ArrayRef<const RVec> x,
272 VirialHandling virialHandling,
273 ArrayRef<RVec> fshift,
277 gmx_wallcycle* wcycle);
280 //! The number of vsites that cross update groups, when =0 no PBC treatment is needed
281 const int numInterUpdategroupVirtualSites_;
282 //! PBC and DD information
283 const DomainInfo domainInfo_;
284 //! The interaction parameters
285 const ArrayRef<const t_iparams> iparams_;
286 //! The interaction lists
287 ArrayRef<const InteractionList> ilists_;
288 //! Information for handling vsite threading
289 ThreadingInfo threadingInfo_;
292 VirtualSitesHandler::~VirtualSitesHandler() = default;
294 int VirtualSitesHandler::numInterUpdategroupVirtualSites() const
296 return impl_->numInterUpdategroupVirtualSites();
299 /*! \brief Returns the sum of the vsite ilist sizes over all vsite types
301 * \param[in] ilist The interaction list
303 static int vsiteIlistNrCount(ArrayRef<const InteractionList> ilist)
306 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
308 nr += ilist[ftype].size();
314 //! Computes the distance between xi and xj, pbc is used when pbc!=nullptr
315 static int pbc_rvec_sub(const t_pbc* pbc, const rvec xi, const rvec xj, rvec dx)
319 return pbc_dx_aiuc(pbc, xi, xj, dx);
323 rvec_sub(xi, xj, dx);
328 //! Returns the 1/norm(x)
329 static inline real inverseNorm(const rvec x)
331 return gmx::invsqrt(iprod(x, x));
335 /* Vsite construction routines */
337 static void constr_vsite1(const rvec xi, rvec x)
344 static void constr_vsite2(const rvec xi, const rvec xj, rvec x, real a, const t_pbc* pbc)
352 pbc_dx_aiuc(pbc, xj, xi, dx);
353 x[XX] = xi[XX] + a * dx[XX];
354 x[YY] = xi[YY] + a * dx[YY];
355 x[ZZ] = xi[ZZ] + a * dx[ZZ];
359 x[XX] = b * xi[XX] + a * xj[XX];
360 x[YY] = b * xi[YY] + a * xj[YY];
361 x[ZZ] = b * xi[ZZ] + a * xj[ZZ];
365 /* TOTAL: 10 flops */
368 static void constr_vsite2FD(const rvec xi, const rvec xj, rvec x, real a, const t_pbc* pbc)
371 pbc_rvec_sub(pbc, xj, xi, xij);
374 const real b = a * inverseNorm(xij);
377 x[XX] = xi[XX] + b * xij[XX];
378 x[YY] = xi[YY] + b * xij[YY];
379 x[ZZ] = xi[ZZ] + b * xij[ZZ];
382 /* TOTAL: 25 flops */
385 static void constr_vsite3(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, const t_pbc* pbc)
394 pbc_dx_aiuc(pbc, xj, xi, dxj);
395 pbc_dx_aiuc(pbc, xk, xi, dxk);
396 x[XX] = xi[XX] + a * dxj[XX] + b * dxk[XX];
397 x[YY] = xi[YY] + a * dxj[YY] + b * dxk[YY];
398 x[ZZ] = xi[ZZ] + a * dxj[ZZ] + b * dxk[ZZ];
402 x[XX] = c * xi[XX] + a * xj[XX] + b * xk[XX];
403 x[YY] = c * xi[YY] + a * xj[YY] + b * xk[YY];
404 x[ZZ] = c * xi[ZZ] + a * xj[ZZ] + b * xk[ZZ];
408 /* TOTAL: 17 flops */
411 static void constr_vsite3FD(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, const t_pbc* pbc)
416 pbc_rvec_sub(pbc, xj, xi, xij);
417 pbc_rvec_sub(pbc, xk, xj, xjk);
420 /* temp goes from i to a point on the line jk */
421 temp[XX] = xij[XX] + a * xjk[XX];
422 temp[YY] = xij[YY] + a * xjk[YY];
423 temp[ZZ] = xij[ZZ] + a * xjk[ZZ];
426 c = b * inverseNorm(temp);
429 x[XX] = xi[XX] + c * temp[XX];
430 x[YY] = xi[YY] + c * temp[YY];
431 x[ZZ] = xi[ZZ] + c * temp[ZZ];
434 /* TOTAL: 34 flops */
437 static void constr_vsite3FAD(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, const t_pbc* pbc)
440 real a1, b1, c1, invdij;
442 pbc_rvec_sub(pbc, xj, xi, xij);
443 pbc_rvec_sub(pbc, xk, xj, xjk);
446 invdij = inverseNorm(xij);
447 c1 = invdij * invdij * iprod(xij, xjk);
448 xp[XX] = xjk[XX] - c1 * xij[XX];
449 xp[YY] = xjk[YY] - c1 * xij[YY];
450 xp[ZZ] = xjk[ZZ] - c1 * xij[ZZ];
452 b1 = b * inverseNorm(xp);
455 x[XX] = xi[XX] + a1 * xij[XX] + b1 * xp[XX];
456 x[YY] = xi[YY] + a1 * xij[YY] + b1 * xp[YY];
457 x[ZZ] = xi[ZZ] + a1 * xij[ZZ] + b1 * xp[ZZ];
460 /* TOTAL: 63 flops */
464 constr_vsite3OUT(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, real c, const t_pbc* pbc)
468 pbc_rvec_sub(pbc, xj, xi, xij);
469 pbc_rvec_sub(pbc, xk, xi, xik);
470 cprod(xij, xik, temp);
473 x[XX] = xi[XX] + a * xij[XX] + b * xik[XX] + c * temp[XX];
474 x[YY] = xi[YY] + a * xij[YY] + b * xik[YY] + c * temp[YY];
475 x[ZZ] = xi[ZZ] + a * xij[ZZ] + b * xik[ZZ] + c * temp[ZZ];
478 /* TOTAL: 33 flops */
481 static void constr_vsite4FD(const rvec xi,
491 rvec xij, xjk, xjl, temp;
494 pbc_rvec_sub(pbc, xj, xi, xij);
495 pbc_rvec_sub(pbc, xk, xj, xjk);
496 pbc_rvec_sub(pbc, xl, xj, xjl);
499 /* temp goes from i to a point on the plane jkl */
500 temp[XX] = xij[XX] + a * xjk[XX] + b * xjl[XX];
501 temp[YY] = xij[YY] + a * xjk[YY] + b * xjl[YY];
502 temp[ZZ] = xij[ZZ] + a * xjk[ZZ] + b * xjl[ZZ];
505 d = c * inverseNorm(temp);
508 x[XX] = xi[XX] + d * temp[XX];
509 x[YY] = xi[YY] + d * temp[YY];
510 x[ZZ] = xi[ZZ] + d * temp[ZZ];
513 /* TOTAL: 43 flops */
516 static void constr_vsite4FDN(const rvec xi,
526 rvec xij, xik, xil, ra, rb, rja, rjb, rm;
529 pbc_rvec_sub(pbc, xj, xi, xij);
530 pbc_rvec_sub(pbc, xk, xi, xik);
531 pbc_rvec_sub(pbc, xl, xi, xil);
534 ra[XX] = a * xik[XX];
535 ra[YY] = a * xik[YY];
536 ra[ZZ] = a * xik[ZZ];
538 rb[XX] = b * xil[XX];
539 rb[YY] = b * xil[YY];
540 rb[ZZ] = b * xil[ZZ];
544 rvec_sub(ra, xij, rja);
545 rvec_sub(rb, xij, rjb);
551 d = c * inverseNorm(rm);
554 x[XX] = xi[XX] + d * rm[XX];
555 x[YY] = xi[YY] + d * rm[YY];
556 x[ZZ] = xi[ZZ] + d * rm[ZZ];
559 /* TOTAL: 47 flops */
562 static int constr_vsiten(const t_iatom* ia, ArrayRef<const t_iparams> ip, ArrayRef<RVec> x, const t_pbc* pbc)
569 n3 = 3 * ip[ia[0]].vsiten.n;
572 copy_rvec(x[ai], x1);
574 for (int i = 3; i < n3; i += 3)
577 a = ip[ia[i]].vsiten.a;
580 pbc_dx_aiuc(pbc, x[ai], x1, dx);
584 rvec_sub(x[ai], x1, dx);
586 dsum[XX] += a * dx[XX];
587 dsum[YY] += a * dx[YY];
588 dsum[ZZ] += a * dx[ZZ];
592 x[av][XX] = x1[XX] + dsum[XX];
593 x[av][YY] = x1[YY] + dsum[YY];
594 x[av][ZZ] = x1[ZZ] + dsum[ZZ];
601 //! PBC modes for vsite construction and spreading
604 all, //!< Apply normal, simple PBC for all vsites
605 none //!< No PBC treatment needed
608 /*! \brief Returns the PBC mode based on the system PBC and vsite properties
610 * \param[in] pbcPtr A pointer to a PBC struct or nullptr when no PBC treatment is required
612 static PbcMode getPbcMode(const t_pbc* pbcPtr)
614 if (pbcPtr == nullptr)
616 return PbcMode::none;
624 /*! \brief Executes the vsite construction task for a single thread
626 * \param[in,out] x Coordinates to construct vsites for
627 * \param[in] dt Time step, needed when v is not empty
628 * \param[in,out] v When not empty, velocities are generated for virtual sites
629 * \param[in] ip Interaction parameters for all interaction, only vsite parameters are used
630 * \param[in] ilist The interaction lists, only vsites are usesd
631 * \param[in] pbc_null PBC struct, used for PBC distance calculations when !=nullptr
633 static void construct_vsites_thread(ArrayRef<RVec> x,
636 ArrayRef<const t_iparams> ip,
637 ArrayRef<const InteractionList> ilist,
638 const t_pbc* pbc_null)
650 const PbcMode pbcMode = getPbcMode(pbc_null);
651 /* We need another pbc pointer, as with charge groups we switch per vsite */
652 const t_pbc* pbc_null2 = pbc_null;
654 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
656 if (ilist[ftype].empty())
662 int nra = interaction_function[ftype].nratoms;
664 int nr = ilist[ftype].size();
666 const t_iatom* ia = ilist[ftype].iatoms.data();
668 for (int i = 0; i < nr;)
671 /* The vsite and constructing atoms */
674 /* Constants for constructing vsites */
675 real a1 = ip[tp].vsite.a;
676 /* Copy the old position */
678 copy_rvec(x[avsite], xv);
680 /* Construct the vsite depending on type */
685 case F_VSITE1: constr_vsite1(x[ai], x[avsite]); break;
688 constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
692 constr_vsite2FD(x[ai], x[aj], x[avsite], a1, pbc_null2);
698 constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
704 constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
710 constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
717 constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
725 constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1, pbc_null2);
733 constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1, pbc_null2);
735 case F_VSITEN: inc = constr_vsiten(ia, ip, x, pbc_null2); break;
737 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d", ftype, __FILE__, __LINE__);
740 if (pbcMode == PbcMode::all)
742 /* Keep the vsite in the same periodic image as before */
744 int ishift = pbc_dx_aiuc(pbc_null, x[avsite], xv, dx);
745 if (ishift != CENTRAL)
747 rvec_add(xv, dx, x[avsite]);
752 /* Calculate velocity of vsite... */
754 rvec_sub(x[avsite], xv, vv);
755 svmul(inv_dt, vv, v[avsite]);
758 /* Increment loop variables */
766 /*! \brief Dispatch the vsite construction tasks for all threads
768 * \param[in] threadingInfo Used to divide work over threads when != nullptr
769 * \param[in,out] x Coordinates to construct vsites for
770 * \param[in] dt Time step, needed when v is not empty
771 * \param[in,out] v When not empty, velocities are generated for virtual sites
772 * \param[in] ip Interaction parameters for all interaction, only vsite parameters are used
773 * \param[in] ilist The interaction lists, only vsites are usesd
774 * \param[in] domainInfo Information about PBC and DD
775 * \param[in] box Used for PBC when PBC is set in domainInfo
777 static void construct_vsites(const ThreadingInfo* threadingInfo,
781 ArrayRef<const t_iparams> ip,
782 ArrayRef<const InteractionList> ilist,
783 const DomainInfo& domainInfo,
786 const bool useDomdec = domainInfo.useDomdec();
788 t_pbc pbc, *pbc_null;
790 /* We only need to do pbc when we have inter update-group vsites.
791 * Note that with domain decomposition we do not need to apply PBC here
792 * when we have at least 3 domains along each dimension. Currently we
793 * do not optimize this case.
795 if (domainInfo.pbcType_ != PbcType::No && domainInfo.useMolPbc_)
797 /* This is wasting some CPU time as we now do this multiple times
801 clear_ivec(null_ivec);
802 pbc_null = set_pbc_dd(&pbc, domainInfo.pbcType_,
803 useDomdec ? domainInfo.domdec_->numCells : null_ivec, FALSE, box);
812 dd_move_x_vsites(*domainInfo.domdec_, box, as_rvec_array(x.data()));
815 if (threadingInfo == nullptr || threadingInfo->numThreads() == 1)
817 construct_vsites_thread(x, dt, v, ip, ilist, pbc_null);
821 #pragma omp parallel num_threads(threadingInfo->numThreads())
825 const int th = gmx_omp_get_thread_num();
826 const VsiteThread& tData = threadingInfo->threadData(th);
827 GMX_ASSERT(tData.rangeStart >= 0,
828 "The thread data should be initialized before calling construct_vsites");
830 construct_vsites_thread(x, dt, v, ip, tData.ilist, pbc_null);
831 if (tData.useInterdependentTask)
833 /* Here we don't need a barrier (unlike the spreading),
834 * since both tasks only construct vsites from particles,
835 * or local vsites, not from non-local vsites.
837 construct_vsites_thread(x, dt, v, ip, tData.idTask.ilist, pbc_null);
840 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
842 /* Now we can construct the vsites that might depend on other vsites */
843 construct_vsites_thread(x, dt, v, ip, threadingInfo->threadDataNonLocalDependent().ilist, pbc_null);
847 void VirtualSitesHandler::Impl::construct(ArrayRef<RVec> x, real dt, ArrayRef<RVec> v, const matrix box) const
849 construct_vsites(&threadingInfo_, x, dt, v, iparams_, ilists_, domainInfo_, box);
852 void VirtualSitesHandler::construct(ArrayRef<RVec> x, real dt, ArrayRef<RVec> v, const matrix box) const
854 impl_->construct(x, dt, v, box);
857 void constructVirtualSites(ArrayRef<RVec> x, ArrayRef<const t_iparams> ip, ArrayRef<const InteractionList> ilist)
861 const DomainInfo domainInfo;
862 construct_vsites(nullptr, x, 0, {}, ip, ilist, domainInfo, nullptr);
866 /* Force spreading routines */
868 static void spread_vsite1(const t_iatom ia[], ArrayRef<RVec> f)
870 const int av = ia[1];
871 const int ai = ia[2];
876 template<VirialHandling virialHandling>
877 static void spread_vsite2(const t_iatom ia[],
879 ArrayRef<const RVec> x,
881 ArrayRef<RVec> fshift,
891 svmul(1 - a, f[av], fi);
899 if (virialHandling == VirialHandling::Pbc)
905 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
906 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
914 if (siv != CENTRAL || sij != CENTRAL)
916 rvec_inc(fshift[siv], f[av]);
917 rvec_dec(fshift[CENTRAL], fi);
918 rvec_dec(fshift[sij], fj);
922 /* TOTAL: 13 flops */
925 void constructVirtualSitesGlobal(const gmx_mtop_t& mtop, gmx::ArrayRef<gmx::RVec> x)
927 GMX_ASSERT(x.ssize() >= mtop.natoms, "x should contain the whole system");
928 GMX_ASSERT(!mtop.moleculeBlockIndices.empty(),
929 "molblock indices are needed in constructVsitesGlobal");
931 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
933 const gmx_molblock_t& molb = mtop.molblock[mb];
934 const gmx_moltype_t& molt = mtop.moltype[molb.type];
935 if (vsiteIlistNrCount(molt.ilist) > 0)
937 int atomOffset = mtop.moleculeBlockIndices[mb].globalAtomStart;
938 for (int mol = 0; mol < molb.nmol; mol++)
940 constructVirtualSites(x.subArray(atomOffset, molt.atoms.nr), mtop.ffparams.iparams,
942 atomOffset += molt.atoms.nr;
948 template<VirialHandling virialHandling>
949 static void spread_vsite2FD(const t_iatom ia[],
951 ArrayRef<const RVec> x,
953 ArrayRef<RVec> fshift,
957 const int av = ia[1];
958 const int ai = ia[2];
959 const int aj = ia[3];
961 copy_rvec(f[av], fv);
964 int sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
967 const real invDistance = inverseNorm(xij);
968 const real b = a * invDistance;
971 const real fproj = iprod(xij, fv) * invDistance * invDistance;
974 fj[XX] = b * (fv[XX] - fproj * xij[XX]);
975 fj[YY] = b * (fv[YY] - fproj * xij[YY]);
976 fj[ZZ] = b * (fv[ZZ] - fproj * xij[ZZ]);
979 /* b is already calculated in constr_vsite2FD
980 storing b somewhere will save flops. */
982 f[ai][XX] += fv[XX] - fj[XX];
983 f[ai][YY] += fv[YY] - fj[YY];
984 f[ai][ZZ] += fv[ZZ] - fj[ZZ];
990 if (virialHandling == VirialHandling::Pbc)
996 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1003 if (svi != CENTRAL || sji != CENTRAL)
1005 rvec_dec(fshift[svi], fv);
1006 fshift[CENTRAL][XX] += fv[XX] - fj[XX];
1007 fshift[CENTRAL][YY] += fv[YY] - fj[YY];
1008 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ];
1009 fshift[sji][XX] += fj[XX];
1010 fshift[sji][YY] += fj[YY];
1011 fshift[sji][ZZ] += fj[ZZ];
1015 if (virialHandling == VirialHandling::NonLinear)
1017 /* Under this condition, the virial for the current forces is not
1018 * calculated from the redistributed forces. This means that
1019 * the effect of non-linear virtual site constructions on the virial
1020 * needs to be added separately. This contribution can be calculated
1021 * in many ways, but the simplest and cheapest way is to use
1022 * the first constructing atom ai as a reference position in space:
1023 * subtract (xv-xi)*fv and add (xj-xi)*fj.
1027 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1029 for (int i = 0; i < DIM; i++)
1031 for (int j = 0; j < DIM; j++)
1033 /* As xix is a linear combination of j and k, use that here */
1034 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * fj[j];
1039 /* TOTAL: 38 flops */
1042 template<VirialHandling virialHandling>
1043 static void spread_vsite3(const t_iatom ia[],
1046 ArrayRef<const RVec> x,
1048 ArrayRef<RVec> fshift,
1051 rvec fi, fj, fk, dx;
1059 svmul(1 - a - b, f[av], fi);
1060 svmul(a, f[av], fj);
1061 svmul(b, f[av], fk);
1064 rvec_inc(f[ai], fi);
1065 rvec_inc(f[aj], fj);
1066 rvec_inc(f[ak], fk);
1069 if (virialHandling == VirialHandling::Pbc)
1076 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
1077 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
1078 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
1087 if (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL)
1089 rvec_inc(fshift[siv], f[av]);
1090 rvec_dec(fshift[CENTRAL], fi);
1091 rvec_dec(fshift[sij], fj);
1092 rvec_dec(fshift[sik], fk);
1096 /* TOTAL: 20 flops */
1099 template<VirialHandling virialHandling>
1100 static void spread_vsite3FD(const t_iatom ia[],
1103 ArrayRef<const RVec> x,
1105 ArrayRef<RVec> fshift,
1110 rvec xvi, xij, xjk, xix, fv, temp;
1111 t_iatom av, ai, aj, ak;
1118 copy_rvec(f[av], fv);
1120 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1121 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1124 /* xix goes from i to point x on the line jk */
1125 xix[XX] = xij[XX] + a * xjk[XX];
1126 xix[YY] = xij[YY] + a * xjk[YY];
1127 xix[ZZ] = xij[ZZ] + a * xjk[ZZ];
1130 const real invDistance = inverseNorm(xix);
1131 const real c = b * invDistance;
1132 /* 4 + ?10? flops */
1134 fproj = iprod(xix, fv) * invDistance * invDistance; /* = (xix . f)/(xix . xix) */
1136 temp[XX] = c * (fv[XX] - fproj * xix[XX]);
1137 temp[YY] = c * (fv[YY] - fproj * xix[YY]);
1138 temp[ZZ] = c * (fv[ZZ] - fproj * xix[ZZ]);
1141 /* c is already calculated in constr_vsite3FD
1142 storing c somewhere will save 26 flops! */
1145 f[ai][XX] += fv[XX] - temp[XX];
1146 f[ai][YY] += fv[YY] - temp[YY];
1147 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
1148 f[aj][XX] += a1 * temp[XX];
1149 f[aj][YY] += a1 * temp[YY];
1150 f[aj][ZZ] += a1 * temp[ZZ];
1151 f[ak][XX] += a * temp[XX];
1152 f[ak][YY] += a * temp[YY];
1153 f[ak][ZZ] += a * temp[ZZ];
1156 if (virialHandling == VirialHandling::Pbc)
1161 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1168 if (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL)
1170 rvec_dec(fshift[svi], fv);
1171 fshift[CENTRAL][XX] += fv[XX] - (1 + a) * temp[XX];
1172 fshift[CENTRAL][YY] += fv[YY] - (1 + a) * temp[YY];
1173 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a) * temp[ZZ];
1174 fshift[sji][XX] += temp[XX];
1175 fshift[sji][YY] += temp[YY];
1176 fshift[sji][ZZ] += temp[ZZ];
1177 fshift[skj][XX] += a * temp[XX];
1178 fshift[skj][YY] += a * temp[YY];
1179 fshift[skj][ZZ] += a * temp[ZZ];
1183 if (virialHandling == VirialHandling::NonLinear)
1185 /* Under this condition, the virial for the current forces is not
1186 * calculated from the redistributed forces. This means that
1187 * the effect of non-linear virtual site constructions on the virial
1188 * needs to be added separately. This contribution can be calculated
1189 * in many ways, but the simplest and cheapest way is to use
1190 * the first constructing atom ai as a reference position in space:
1191 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
1195 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1197 for (int i = 0; i < DIM; i++)
1199 for (int j = 0; j < DIM; j++)
1201 /* As xix is a linear combination of j and k, use that here */
1202 dxdf[i][j] += -xiv[i] * fv[j] + xix[i] * temp[j];
1207 /* TOTAL: 61 flops */
1210 template<VirialHandling virialHandling>
1211 static void spread_vsite3FAD(const t_iatom ia[],
1214 ArrayRef<const RVec> x,
1216 ArrayRef<RVec> fshift,
1220 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
1221 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
1222 t_iatom av, ai, aj, ak;
1229 copy_rvec(f[ia[1]], fv);
1231 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1232 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1235 invdij = inverseNorm(xij);
1236 invdij2 = invdij * invdij;
1237 c1 = iprod(xij, xjk) * invdij2;
1238 xperp[XX] = xjk[XX] - c1 * xij[XX];
1239 xperp[YY] = xjk[YY] - c1 * xij[YY];
1240 xperp[ZZ] = xjk[ZZ] - c1 * xij[ZZ];
1241 /* xperp in plane ijk, perp. to ij */
1242 invdp = inverseNorm(xperp);
1247 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
1248 storing them somewhere will save 45 flops! */
1250 fproj = iprod(xij, fv) * invdij2;
1251 svmul(fproj, xij, Fpij); /* proj. f on xij */
1252 svmul(iprod(xperp, fv) * invdp * invdp, xperp, Fppp); /* proj. f on xperp */
1253 svmul(b1 * fproj, xperp, f3);
1256 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
1257 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
1258 for (int d = 0; d < DIM; d++)
1266 f[ai][XX] += fv[XX] - f1[XX] + c1 * f2[XX] + f3[XX];
1267 f[ai][YY] += fv[YY] - f1[YY] + c1 * f2[YY] + f3[YY];
1268 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1 * f2[ZZ] + f3[ZZ];
1269 f[aj][XX] += f1[XX] - c2 * f2[XX] - f3[XX];
1270 f[aj][YY] += f1[YY] - c2 * f2[YY] - f3[YY];
1271 f[aj][ZZ] += f1[ZZ] - c2 * f2[ZZ] - f3[ZZ];
1272 f[ak][XX] += f2[XX];
1273 f[ak][YY] += f2[YY];
1274 f[ak][ZZ] += f2[ZZ];
1277 if (virialHandling == VirialHandling::Pbc)
1283 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1290 if (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL)
1292 rvec_dec(fshift[svi], fv);
1293 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1 - c1) * f2[XX] + f3[XX];
1294 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1 - c1) * f2[YY] + f3[YY];
1295 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1 - c1) * f2[ZZ] + f3[ZZ];
1296 fshift[sji][XX] += f1[XX] - c1 * f2[XX] - f3[XX];
1297 fshift[sji][YY] += f1[YY] - c1 * f2[YY] - f3[YY];
1298 fshift[sji][ZZ] += f1[ZZ] - c1 * f2[ZZ] - f3[ZZ];
1299 fshift[skj][XX] += f2[XX];
1300 fshift[skj][YY] += f2[YY];
1301 fshift[skj][ZZ] += f2[ZZ];
1305 if (virialHandling == VirialHandling::NonLinear)
1308 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1310 for (int i = 0; i < DIM; i++)
1312 for (int j = 0; j < DIM; j++)
1314 /* Note that xik=xij+xjk, so we have to add xij*f2 */
1315 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * (f1[j] + (1 - c2) * f2[j] - f3[j])
1321 /* TOTAL: 113 flops */
1324 template<VirialHandling virialHandling>
1325 static void spread_vsite3OUT(const t_iatom ia[],
1329 ArrayRef<const RVec> x,
1331 ArrayRef<RVec> fshift,
1335 rvec xvi, xij, xik, fv, fj, fk;
1345 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1346 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1349 copy_rvec(f[av], fv);
1356 fj[XX] = a * fv[XX] - xik[ZZ] * cfy + xik[YY] * cfz;
1357 fj[YY] = xik[ZZ] * cfx + a * fv[YY] - xik[XX] * cfz;
1358 fj[ZZ] = -xik[YY] * cfx + xik[XX] * cfy + a * fv[ZZ];
1360 fk[XX] = b * fv[XX] + xij[ZZ] * cfy - xij[YY] * cfz;
1361 fk[YY] = -xij[ZZ] * cfx + b * fv[YY] + xij[XX] * cfz;
1362 fk[ZZ] = xij[YY] * cfx - xij[XX] * cfy + b * fv[ZZ];
1365 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
1366 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
1367 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1368 rvec_inc(f[aj], fj);
1369 rvec_inc(f[ak], fk);
1372 if (virialHandling == VirialHandling::Pbc)
1377 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1384 if (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL)
1386 rvec_dec(fshift[svi], fv);
1387 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
1388 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
1389 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1390 rvec_inc(fshift[sji], fj);
1391 rvec_inc(fshift[ski], fk);
1395 if (virialHandling == VirialHandling::NonLinear)
1399 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1401 for (int i = 0; i < DIM; i++)
1403 for (int j = 0; j < DIM; j++)
1405 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * fj[j] + xik[i] * fk[j];
1410 /* TOTAL: 54 flops */
1413 template<VirialHandling virialHandling>
1414 static void spread_vsite4FD(const t_iatom ia[],
1418 ArrayRef<const RVec> x,
1420 ArrayRef<RVec> fshift,
1425 rvec xvi, xij, xjk, xjl, xix, fv, temp;
1426 int av, ai, aj, ak, al;
1427 int sji, skj, slj, m;
1435 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1436 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1437 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1440 /* xix goes from i to point x on the plane jkl */
1441 for (m = 0; m < DIM; m++)
1443 xix[m] = xij[m] + a * xjk[m] + b * xjl[m];
1447 const real invDistance = inverseNorm(xix);
1448 const real d = c * invDistance;
1449 /* 4 + ?10? flops */
1451 copy_rvec(f[av], fv);
1453 fproj = iprod(xix, fv) * invDistance * invDistance; /* = (xix . f)/(xix . xix) */
1455 for (m = 0; m < DIM; m++)
1457 temp[m] = d * (fv[m] - fproj * xix[m]);
1461 /* c is already calculated in constr_vsite3FD
1462 storing c somewhere will save 35 flops! */
1465 for (m = 0; m < DIM; m++)
1467 f[ai][m] += fv[m] - temp[m];
1468 f[aj][m] += a1 * temp[m];
1469 f[ak][m] += a * temp[m];
1470 f[al][m] += b * temp[m];
1474 if (virialHandling == VirialHandling::Pbc)
1479 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1486 if (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL)
1488 rvec_dec(fshift[svi], fv);
1489 for (m = 0; m < DIM; m++)
1491 fshift[CENTRAL][m] += fv[m] - (1 + a + b) * temp[m];
1492 fshift[sji][m] += temp[m];
1493 fshift[skj][m] += a * temp[m];
1494 fshift[slj][m] += b * temp[m];
1499 if (virialHandling == VirialHandling::NonLinear)
1504 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1506 for (i = 0; i < DIM; i++)
1508 for (j = 0; j < DIM; j++)
1510 dxdf[i][j] += -xiv[i] * fv[j] + xix[i] * temp[j];
1515 /* TOTAL: 77 flops */
1518 template<VirialHandling virialHandling>
1519 static void spread_vsite4FDN(const t_iatom ia[],
1523 ArrayRef<const RVec> x,
1525 ArrayRef<RVec> fshift,
1529 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1530 rvec fv, fj, fk, fl;
1533 int av, ai, aj, ak, al;
1536 /* DEBUG: check atom indices */
1543 copy_rvec(f[av], fv);
1545 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1546 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1547 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1550 ra[XX] = a * xik[XX];
1551 ra[YY] = a * xik[YY];
1552 ra[ZZ] = a * xik[ZZ];
1554 rb[XX] = b * xil[XX];
1555 rb[YY] = b * xil[YY];
1556 rb[ZZ] = b * xil[ZZ];
1560 rvec_sub(ra, xij, rja);
1561 rvec_sub(rb, xij, rjb);
1562 rvec_sub(rb, ra, rab);
1565 cprod(rja, rjb, rm);
1568 invrm = inverseNorm(rm);
1569 denom = invrm * invrm;
1572 cfx = c * invrm * fv[XX];
1573 cfy = c * invrm * fv[YY];
1574 cfz = c * invrm * fv[ZZ];
1585 fj[XX] = (-rm[XX] * rt[XX]) * cfx + (rab[ZZ] - rm[YY] * rt[XX]) * cfy
1586 + (-rab[YY] - rm[ZZ] * rt[XX]) * cfz;
1587 fj[YY] = (-rab[ZZ] - rm[XX] * rt[YY]) * cfx + (-rm[YY] * rt[YY]) * cfy
1588 + (rab[XX] - rm[ZZ] * rt[YY]) * cfz;
1589 fj[ZZ] = (rab[YY] - rm[XX] * rt[ZZ]) * cfx + (-rab[XX] - rm[YY] * rt[ZZ]) * cfy
1590 + (-rm[ZZ] * rt[ZZ]) * cfz;
1596 rt[XX] *= denom * a;
1597 rt[YY] *= denom * a;
1598 rt[ZZ] *= denom * a;
1601 fk[XX] = (-rm[XX] * rt[XX]) * cfx + (-a * rjb[ZZ] - rm[YY] * rt[XX]) * cfy
1602 + (a * rjb[YY] - rm[ZZ] * rt[XX]) * cfz;
1603 fk[YY] = (a * rjb[ZZ] - rm[XX] * rt[YY]) * cfx + (-rm[YY] * rt[YY]) * cfy
1604 + (-a * rjb[XX] - rm[ZZ] * rt[YY]) * cfz;
1605 fk[ZZ] = (-a * rjb[YY] - rm[XX] * rt[ZZ]) * cfx + (a * rjb[XX] - rm[YY] * rt[ZZ]) * cfy
1606 + (-rm[ZZ] * rt[ZZ]) * cfz;
1612 rt[XX] *= denom * b;
1613 rt[YY] *= denom * b;
1614 rt[ZZ] *= denom * b;
1617 fl[XX] = (-rm[XX] * rt[XX]) * cfx + (b * rja[ZZ] - rm[YY] * rt[XX]) * cfy
1618 + (-b * rja[YY] - rm[ZZ] * rt[XX]) * cfz;
1619 fl[YY] = (-b * rja[ZZ] - rm[XX] * rt[YY]) * cfx + (-rm[YY] * rt[YY]) * cfy
1620 + (b * rja[XX] - rm[ZZ] * rt[YY]) * cfz;
1621 fl[ZZ] = (b * rja[YY] - rm[XX] * rt[ZZ]) * cfx + (-b * rja[XX] - rm[YY] * rt[ZZ]) * cfy
1622 + (-rm[ZZ] * rt[ZZ]) * cfz;
1625 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1626 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1627 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1628 rvec_inc(f[aj], fj);
1629 rvec_inc(f[ak], fk);
1630 rvec_inc(f[al], fl);
1633 if (virialHandling == VirialHandling::Pbc)
1638 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1645 if (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL)
1647 rvec_dec(fshift[svi], fv);
1648 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1649 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1650 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1651 rvec_inc(fshift[sij], fj);
1652 rvec_inc(fshift[sik], fk);
1653 rvec_inc(fshift[sil], fl);
1657 if (virialHandling == VirialHandling::NonLinear)
1662 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1664 for (i = 0; i < DIM; i++)
1666 for (j = 0; j < DIM; j++)
1668 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * fj[j] + xik[i] * fk[j] + xil[i] * fl[j];
1673 /* Total: 207 flops (Yuck!) */
1676 template<VirialHandling virialHandling>
1677 static int spread_vsiten(const t_iatom ia[],
1678 ArrayRef<const t_iparams> ip,
1679 ArrayRef<const RVec> x,
1681 ArrayRef<RVec> fshift,
1689 n3 = 3 * ip[ia[0]].vsiten.n;
1691 copy_rvec(x[av], xv);
1693 for (i = 0; i < n3; i += 3)
1698 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1704 a = ip[ia[i]].vsiten.a;
1705 svmul(a, f[av], fi);
1706 rvec_inc(f[ai], fi);
1708 if (virialHandling == VirialHandling::Pbc && siv != CENTRAL)
1710 rvec_inc(fshift[siv], fi);
1711 rvec_dec(fshift[CENTRAL], fi);
1721 //! Returns the number of virtual sites in the interaction list, for VSITEN the number of atoms
1722 static int vsite_count(ArrayRef<const InteractionList> ilist, int ftype)
1724 if (ftype == F_VSITEN)
1726 return ilist[ftype].size() / 3;
1730 return ilist[ftype].size() / (1 + interaction_function[ftype].nratoms);
1734 //! Executes the force spreading task for a single thread
1735 template<VirialHandling virialHandling>
1736 static void spreadForceForThread(ArrayRef<const RVec> x,
1738 ArrayRef<RVec> fshift,
1740 ArrayRef<const t_iparams> ip,
1741 ArrayRef<const InteractionList> ilist,
1742 const t_pbc* pbc_null)
1744 const PbcMode pbcMode = getPbcMode(pbc_null);
1745 /* We need another pbc pointer, as with charge groups we switch per vsite */
1746 const t_pbc* pbc_null2 = pbc_null;
1747 gmx::ArrayRef<const int> vsite_pbc;
1749 /* this loop goes backwards to be able to build *
1750 * higher type vsites from lower types */
1751 for (int ftype = c_ftypeVsiteEnd - 1; ftype >= c_ftypeVsiteStart; ftype--)
1753 if (ilist[ftype].empty())
1759 int nra = interaction_function[ftype].nratoms;
1761 int nr = ilist[ftype].size();
1763 const t_iatom* ia = ilist[ftype].iatoms.data();
1765 if (pbcMode == PbcMode::all)
1767 pbc_null2 = pbc_null;
1770 for (int i = 0; i < nr;)
1774 /* Constants for constructing */
1776 a1 = ip[tp].vsite.a;
1777 /* Construct the vsite depending on type */
1780 case F_VSITE1: spread_vsite1(ia, f); break;
1782 spread_vsite2<virialHandling>(ia, a1, x, f, fshift, pbc_null2);
1785 spread_vsite2FD<virialHandling>(ia, a1, x, f, fshift, dxdf, pbc_null2);
1788 b1 = ip[tp].vsite.b;
1789 spread_vsite3<virialHandling>(ia, a1, b1, x, f, fshift, pbc_null2);
1792 b1 = ip[tp].vsite.b;
1793 spread_vsite3FD<virialHandling>(ia, a1, b1, x, f, fshift, dxdf, pbc_null2);
1796 b1 = ip[tp].vsite.b;
1797 spread_vsite3FAD<virialHandling>(ia, a1, b1, x, f, fshift, dxdf, pbc_null2);
1800 b1 = ip[tp].vsite.b;
1801 c1 = ip[tp].vsite.c;
1802 spread_vsite3OUT<virialHandling>(ia, a1, b1, c1, x, f, fshift, dxdf, pbc_null2);
1805 b1 = ip[tp].vsite.b;
1806 c1 = ip[tp].vsite.c;
1807 spread_vsite4FD<virialHandling>(ia, a1, b1, c1, x, f, fshift, dxdf, pbc_null2);
1810 b1 = ip[tp].vsite.b;
1811 c1 = ip[tp].vsite.c;
1812 spread_vsite4FDN<virialHandling>(ia, a1, b1, c1, x, f, fshift, dxdf, pbc_null2);
1815 inc = spread_vsiten<virialHandling>(ia, ip, x, f, fshift, pbc_null2);
1818 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d", ftype, __FILE__, __LINE__);
1820 clear_rvec(f[ia[1]]);
1822 /* Increment loop variables */
1830 //! Wrapper function for calling the templated thread-local spread function
1831 static void spreadForceWrapper(ArrayRef<const RVec> x,
1833 const VirialHandling virialHandling,
1834 ArrayRef<RVec> fshift,
1836 const bool clearDxdf,
1837 ArrayRef<const t_iparams> ip,
1838 ArrayRef<const InteractionList> ilist,
1839 const t_pbc* pbc_null)
1841 if (virialHandling == VirialHandling::NonLinear && clearDxdf)
1846 switch (virialHandling)
1848 case VirialHandling::None:
1849 spreadForceForThread<VirialHandling::None>(x, f, fshift, dxdf, ip, ilist, pbc_null);
1851 case VirialHandling::Pbc:
1852 spreadForceForThread<VirialHandling::Pbc>(x, f, fshift, dxdf, ip, ilist, pbc_null);
1854 case VirialHandling::NonLinear:
1855 spreadForceForThread<VirialHandling::NonLinear>(x, f, fshift, dxdf, ip, ilist, pbc_null);
1860 //! Clears the task force buffer elements that are written by task idTask
1861 static void clearTaskForceBufferUsedElements(InterdependentTask* idTask)
1863 int ntask = idTask->spreadTask.size();
1864 for (int ti = 0; ti < ntask; ti++)
1866 const AtomIndex* atomList = &idTask->atomIndex[idTask->spreadTask[ti]];
1867 int natom = atomList->atom.size();
1868 RVec* force = idTask->force.data();
1869 for (int i = 0; i < natom; i++)
1871 clear_rvec(force[atomList->atom[i]]);
1876 void VirtualSitesHandler::Impl::spreadForces(ArrayRef<const RVec> x,
1878 const VirialHandling virialHandling,
1879 ArrayRef<RVec> fshift,
1883 gmx_wallcycle* wcycle)
1885 wallcycle_start(wcycle, ewcVSITESPREAD);
1887 const bool useDomdec = domainInfo_.useDomdec();
1889 t_pbc pbc, *pbc_null;
1891 if (domainInfo_.useMolPbc_)
1893 /* This is wasting some CPU time as we now do this multiple times
1896 pbc_null = set_pbc_dd(&pbc, domainInfo_.pbcType_,
1897 useDomdec ? domainInfo_.domdec_->numCells : nullptr, FALSE, box);
1906 dd_clear_f_vsites(*domainInfo_.domdec_, f);
1909 const int numThreads = threadingInfo_.numThreads();
1911 if (numThreads == 1)
1914 spreadForceWrapper(x, f, virialHandling, fshift, dxdf, true, iparams_, ilists_, pbc_null);
1916 if (virialHandling == VirialHandling::NonLinear)
1918 for (int i = 0; i < DIM; i++)
1920 for (int j = 0; j < DIM; j++)
1922 virial[i][j] += -0.5 * dxdf[i][j];
1929 /* First spread the vsites that might depend on non-local vsites */
1930 auto& nlDependentVSites = threadingInfo_.threadDataNonLocalDependent();
1931 spreadForceWrapper(x, f, virialHandling, fshift, nlDependentVSites.dxdf, true, iparams_,
1932 nlDependentVSites.ilist, pbc_null);
1934 #pragma omp parallel num_threads(numThreads)
1938 int thread = gmx_omp_get_thread_num();
1939 VsiteThread& tData = threadingInfo_.threadData(thread);
1941 ArrayRef<RVec> fshift_t;
1942 if (virialHandling == VirialHandling::Pbc)
1950 fshift_t = tData.fshift;
1952 for (int i = 0; i < SHIFTS; i++)
1954 clear_rvec(fshift_t[i]);
1959 if (tData.useInterdependentTask)
1961 /* Spread the vsites that spread outside our local range.
1962 * This is done using a thread-local force buffer force.
1963 * First we need to copy the input vsite forces to force.
1965 InterdependentTask* idTask = &tData.idTask;
1967 /* Clear the buffer elements set by our task during
1968 * the last call to spread_vsite_f.
1970 clearTaskForceBufferUsedElements(idTask);
1972 int nvsite = idTask->vsite.size();
1973 for (int i = 0; i < nvsite; i++)
1975 copy_rvec(f[idTask->vsite[i]], idTask->force[idTask->vsite[i]]);
1977 spreadForceWrapper(x, idTask->force, virialHandling, fshift_t, tData.dxdf, true,
1978 iparams_, tData.idTask.ilist, pbc_null);
1980 /* We need a barrier before reducing forces below
1981 * that have been produced by a different thread above.
1985 /* Loop over all thread task and reduce forces they
1986 * produced on atoms that fall in our range.
1987 * Note that atomic reduction would be a simpler solution,
1988 * but that might not have good support on all platforms.
1990 int ntask = idTask->reduceTask.size();
1991 for (int ti = 0; ti < ntask; ti++)
1993 const InterdependentTask& idt_foreign =
1994 threadingInfo_.threadData(idTask->reduceTask[ti]).idTask;
1995 const AtomIndex& atomList = idt_foreign.atomIndex[thread];
1996 const RVec* f_foreign = idt_foreign.force.data();
1998 for (int ind : atomList.atom)
2000 rvec_inc(f[ind], f_foreign[ind]);
2001 /* Clearing of f_foreign is done at the next step */
2004 /* Clear the vsite forces, both in f and force */
2005 for (int i = 0; i < nvsite; i++)
2007 int ind = tData.idTask.vsite[i];
2009 clear_rvec(tData.idTask.force[ind]);
2013 /* Spread the vsites that spread locally only */
2014 spreadForceWrapper(x, f, virialHandling, fshift_t, tData.dxdf, false, iparams_,
2015 tData.ilist, pbc_null);
2017 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2020 if (virialHandling == VirialHandling::Pbc)
2022 for (int th = 1; th < numThreads; th++)
2024 for (int i = 0; i < SHIFTS; i++)
2026 rvec_inc(fshift[i], threadingInfo_.threadData(th).fshift[i]);
2031 if (virialHandling == VirialHandling::NonLinear)
2033 for (int th = 0; th < numThreads + 1; th++)
2035 /* MSVC doesn't like matrix references, so we use a pointer */
2036 const matrix& dxdf = threadingInfo_.threadData(th).dxdf;
2038 for (int i = 0; i < DIM; i++)
2040 for (int j = 0; j < DIM; j++)
2042 virial[i][j] += -0.5 * dxdf[i][j];
2051 dd_move_f_vsites(*domainInfo_.domdec_, f, fshift);
2054 inc_nrnb(nrnb, eNR_VSITE1, vsite_count(ilists_, F_VSITE1));
2055 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(ilists_, F_VSITE2));
2056 inc_nrnb(nrnb, eNR_VSITE2FD, vsite_count(ilists_, F_VSITE2FD));
2057 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(ilists_, F_VSITE3));
2058 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(ilists_, F_VSITE3FD));
2059 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(ilists_, F_VSITE3FAD));
2060 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(ilists_, F_VSITE3OUT));
2061 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(ilists_, F_VSITE4FD));
2062 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(ilists_, F_VSITE4FDN));
2063 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(ilists_, F_VSITEN));
2065 wallcycle_stop(wcycle, ewcVSITESPREAD);
2068 /*! \brief Returns the an array with group indices for each atom
2070 * \param[in] grouping The paritioning of the atom range into atom groups
2072 static std::vector<int> makeAtomToGroupMapping(const gmx::RangePartitioning& grouping)
2074 std::vector<int> atomToGroup(grouping.fullRange().end(), 0);
2076 for (int group = 0; group < grouping.numBlocks(); group++)
2078 auto block = grouping.block(group);
2079 std::fill(atomToGroup.begin() + block.begin(), atomToGroup.begin() + block.end(), group);
2085 int countNonlinearVsites(const gmx_mtop_t& mtop)
2087 int numNonlinearVsites = 0;
2088 for (const gmx_molblock_t& molb : mtop.molblock)
2090 const gmx_moltype_t& molt = mtop.moltype[molb.type];
2092 for (const auto& ilist : extractILists(molt.ilist, IF_VSITE))
2094 if (ilist.functionType != F_VSITE2 && ilist.functionType != F_VSITE3
2095 && ilist.functionType != F_VSITEN)
2097 numNonlinearVsites += molb.nmol * ilist.iatoms.size() / (1 + NRAL(ilist.functionType));
2102 return numNonlinearVsites;
2105 void VirtualSitesHandler::spreadForces(ArrayRef<const RVec> x,
2107 const VirialHandling virialHandling,
2108 ArrayRef<RVec> fshift,
2112 gmx_wallcycle* wcycle)
2114 impl_->spreadForces(x, f, virialHandling, fshift, virial, nrnb, box, wcycle);
2117 int countInterUpdategroupVsites(const gmx_mtop_t& mtop,
2118 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingPerMoleculetype)
2120 int n_intercg_vsite = 0;
2121 for (const gmx_molblock_t& molb : mtop.molblock)
2123 const gmx_moltype_t& molt = mtop.moltype[molb.type];
2125 std::vector<int> atomToGroup;
2126 if (!updateGroupingPerMoleculetype.empty())
2128 atomToGroup = makeAtomToGroupMapping(updateGroupingPerMoleculetype[molb.type]);
2130 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2132 const int nral = NRAL(ftype);
2133 const InteractionList& il = molt.ilist[ftype];
2134 for (int i = 0; i < il.size(); i += 1 + nral)
2136 bool isInterGroup = atomToGroup.empty();
2139 const int group = atomToGroup[il.iatoms[1 + i]];
2140 for (int a = 1; a < nral; a++)
2142 if (atomToGroup[il.iatoms[1 + a]] != group)
2144 isInterGroup = true;
2151 n_intercg_vsite += molb.nmol;
2157 return n_intercg_vsite;
2160 std::unique_ptr<VirtualSitesHandler> makeVirtualSitesHandler(const gmx_mtop_t& mtop,
2161 const t_commrec* cr,
2164 GMX_RELEASE_ASSERT(cr != nullptr, "We need a valid commrec");
2166 std::unique_ptr<VirtualSitesHandler> vsite;
2168 /* check if there are vsites */
2170 for (int ftype = 0; ftype < F_NRE; ftype++)
2172 if (interaction_function[ftype].flags & IF_VSITE)
2174 GMX_ASSERT(ftype >= c_ftypeVsiteStart && ftype < c_ftypeVsiteEnd,
2175 "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
2177 nvsite += gmx_mtop_ftype_count(&mtop, ftype);
2181 GMX_ASSERT(ftype < c_ftypeVsiteStart || ftype >= c_ftypeVsiteEnd,
2182 "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
2191 return std::make_unique<VirtualSitesHandler>(mtop, cr->dd, pbcType);
2194 ThreadingInfo::ThreadingInfo() : numThreads_(gmx_omp_nthreads_get(emntVSITE))
2196 if (numThreads_ > 1)
2198 /* We need one extra thread data structure for the overlap vsites */
2199 tData_.resize(numThreads_ + 1);
2200 #pragma omp parallel for num_threads(numThreads_) schedule(static)
2201 for (int thread = 0; thread < numThreads_; thread++)
2205 tData_[thread] = std::make_unique<VsiteThread>();
2207 InterdependentTask& idTask = tData_[thread]->idTask;
2209 idTask.atomIndex.resize(numThreads_);
2211 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2213 if (numThreads_ > 1)
2215 tData_[numThreads_] = std::make_unique<VsiteThread>();
2220 //! Returns the number of inter update-group vsites
2221 static int getNumInterUpdategroupVsites(const gmx_mtop_t& mtop, const gmx_domdec_t* domdec)
2223 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingPerMoleculetype;
2226 updateGroupingPerMoleculetype = getUpdateGroupingPerMoleculetype(*domdec);
2229 return countInterUpdategroupVsites(mtop, updateGroupingPerMoleculetype);
2232 VirtualSitesHandler::Impl::Impl(const gmx_mtop_t& mtop, gmx_domdec_t* domdec, const PbcType pbcType) :
2233 numInterUpdategroupVirtualSites_(getNumInterUpdategroupVsites(mtop, domdec)),
2234 domainInfo_({ pbcType, pbcType != PbcType::No && numInterUpdategroupVirtualSites_ > 0, domdec }),
2235 iparams_(mtop.ffparams.iparams)
2239 VirtualSitesHandler::VirtualSitesHandler(const gmx_mtop_t& mtop, gmx_domdec_t* domdec, const PbcType pbcType) :
2240 impl_(new Impl(mtop, domdec, pbcType))
2244 //! Flag that atom \p atom which is home in another task, if it has not already been added before
2245 static inline void flagAtom(InterdependentTask* idTask, const int atom, const int numThreads, const int numAtomsPerThread)
2247 if (!idTask->use[atom])
2249 idTask->use[atom] = true;
2250 int thread = atom / numAtomsPerThread;
2251 /* Assign all non-local atom force writes to thread 0 */
2252 if (thread >= numThreads)
2256 idTask->atomIndex[thread].atom.push_back(atom);
2260 /*! \brief Here we try to assign all vsites that are in our local range.
2262 * Our task local atom range is tData->rangeStart - tData->rangeEnd.
2263 * Vsites that depend only on local atoms, as indicated by taskIndex[]==thread,
2264 * are assigned to task tData->ilist. Vsites that depend on non-local atoms
2265 * but not on other vsites are assigned to task tData->id_task.ilist.
2266 * taskIndex[] is set for all vsites in our range, either to our local tasks
2267 * or to the single last task as taskIndex[]=2*nthreads.
2269 static void assignVsitesToThread(VsiteThread* tData,
2273 gmx::ArrayRef<int> taskIndex,
2274 ArrayRef<const InteractionList> ilist,
2275 ArrayRef<const t_iparams> ip,
2276 const unsigned short* ptype)
2278 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2280 tData->ilist[ftype].clear();
2281 tData->idTask.ilist[ftype].clear();
2283 int nral1 = 1 + NRAL(ftype);
2285 const int* iat = ilist[ftype].iatoms.data();
2286 for (int i = 0; i < ilist[ftype].size();)
2288 if (ftype == F_VSITEN)
2290 /* The 3 below is from 1+NRAL(ftype)=3 */
2291 inc = ip[iat[i]].vsiten.n * 3;
2294 if (iat[1 + i] < tData->rangeStart || iat[1 + i] >= tData->rangeEnd)
2296 /* This vsite belongs to a different thread */
2301 /* We would like to assign this vsite to task thread,
2302 * but it might depend on atoms outside the atom range of thread
2303 * or on another vsite not assigned to task thread.
2306 if (ftype != F_VSITEN)
2308 for (int j = i + 2; j < i + nral1; j++)
2310 /* Do a range check to avoid a harmless race on taskIndex */
2311 if (iat[j] < tData->rangeStart || iat[j] >= tData->rangeEnd || taskIndex[iat[j]] != thread)
2313 if (!tData->useInterdependentTask || ptype[iat[j]] == eptVSite)
2315 /* At least one constructing atom is a vsite
2316 * that is not assigned to the same thread.
2317 * Put this vsite into a separate task.
2323 /* There are constructing atoms outside our range,
2324 * put this vsite into a second task to be executed
2325 * on the same thread. During construction no barrier
2326 * is needed between the two tasks on the same thread.
2327 * During spreading we need to run this task with
2328 * an additional thread-local intermediate force buffer
2329 * (or atomic reduction) and a barrier between the two
2332 task = nthread + thread;
2338 for (int j = i + 2; j < i + inc; j += 3)
2340 /* Do a range check to avoid a harmless race on taskIndex */
2341 if (iat[j] < tData->rangeStart || iat[j] >= tData->rangeEnd || taskIndex[iat[j]] != thread)
2343 GMX_ASSERT(ptype[iat[j]] != eptVSite,
2344 "A vsite to be assigned in assignVsitesToThread has a vsite as "
2345 "a constructing atom that does not belong to our task, such "
2346 "vsites should be assigned to the single 'master' task");
2348 task = nthread + thread;
2353 /* Update this vsite's thread index entry */
2354 taskIndex[iat[1 + i]] = task;
2356 if (task == thread || task == nthread + thread)
2358 /* Copy this vsite to the thread data struct of thread */
2359 InteractionList* il_task;
2362 il_task = &tData->ilist[ftype];
2366 il_task = &tData->idTask.ilist[ftype];
2368 /* Copy the vsite data to the thread-task local array */
2369 il_task->push_back(iat[i], nral1 - 1, iat + i + 1);
2370 if (task == nthread + thread)
2372 /* This vsite write outside our own task force block.
2373 * Put it into the interdependent task list and flag
2374 * the atoms involved for reduction.
2376 tData->idTask.vsite.push_back(iat[i + 1]);
2377 if (ftype != F_VSITEN)
2379 for (int j = i + 2; j < i + nral1; j++)
2381 flagAtom(&tData->idTask, iat[j], nthread, natperthread);
2386 for (int j = i + 2; j < i + inc; j += 3)
2388 flagAtom(&tData->idTask, iat[j], nthread, natperthread);
2399 /*! \brief Assign all vsites with taskIndex[]==task to task tData */
2400 static void assignVsitesToSingleTask(VsiteThread* tData,
2402 gmx::ArrayRef<const int> taskIndex,
2403 ArrayRef<const InteractionList> ilist,
2404 ArrayRef<const t_iparams> ip)
2406 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2408 tData->ilist[ftype].clear();
2409 tData->idTask.ilist[ftype].clear();
2411 int nral1 = 1 + NRAL(ftype);
2413 const int* iat = ilist[ftype].iatoms.data();
2414 InteractionList* il_task = &tData->ilist[ftype];
2416 for (int i = 0; i < ilist[ftype].size();)
2418 if (ftype == F_VSITEN)
2420 /* The 3 below is from 1+NRAL(ftype)=3 */
2421 inc = ip[iat[i]].vsiten.n * 3;
2423 /* Check if the vsite is assigned to our task */
2424 if (taskIndex[iat[1 + i]] == task)
2426 /* Copy the vsite data to the thread-task local array */
2427 il_task->push_back(iat[i], inc - 1, iat + i + 1);
2435 void ThreadingInfo::setVirtualSites(ArrayRef<const InteractionList> ilists,
2436 ArrayRef<const t_iparams> iparams,
2437 const t_mdatoms& mdatoms,
2438 const bool useDomdec)
2440 if (numThreads_ <= 1)
2446 /* The current way of distributing the vsites over threads in primitive.
2447 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
2448 * without taking into account how the vsites are distributed.
2449 * Without domain decomposition we at least tighten the upper bound
2450 * of the range (useful for common systems such as a vsite-protein
2452 * With domain decomposition, as long as the vsites are distributed
2453 * uniformly in each domain along the major dimension, usually x,
2454 * it will also perform well.
2456 int vsite_atom_range;
2460 vsite_atom_range = -1;
2461 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2464 if (ftype != F_VSITEN)
2466 int nral1 = 1 + NRAL(ftype);
2467 ArrayRef<const int> iat = ilists[ftype].iatoms;
2468 for (int i = 0; i < ilists[ftype].size(); i += nral1)
2470 for (int j = i + 1; j < i + nral1; j++)
2472 vsite_atom_range = std::max(vsite_atom_range, iat[j]);
2480 ArrayRef<const int> iat = ilists[ftype].iatoms;
2483 while (i < ilists[ftype].size())
2485 /* The 3 below is from 1+NRAL(ftype)=3 */
2486 vs_ind_end = i + iparams[iat[i]].vsiten.n * 3;
2488 vsite_atom_range = std::max(vsite_atom_range, iat[i + 1]);
2489 while (i < vs_ind_end)
2491 vsite_atom_range = std::max(vsite_atom_range, iat[i + 2]);
2499 natperthread = (vsite_atom_range + numThreads_ - 1) / numThreads_;
2503 /* Any local or not local atom could be involved in virtual sites.
2504 * But since we usually have very few non-local virtual sites
2505 * (only non-local vsites that depend on local vsites),
2506 * we distribute the local atom range equally over the threads.
2507 * When assigning vsites to threads, we should take care that the last
2508 * threads also covers the non-local range.
2510 vsite_atom_range = mdatoms.nr;
2511 natperthread = (mdatoms.homenr + numThreads_ - 1) / numThreads_;
2516 fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n",
2517 mdatoms.nr, vsite_atom_range, natperthread);
2520 /* To simplify the vsite assignment, we make an index which tells us
2521 * to which task particles, both non-vsites and vsites, are assigned.
2523 taskIndex_.resize(mdatoms.nr);
2525 /* Initialize the task index array. Here we assign the non-vsite
2526 * particles to task=thread, so we easily figure out if vsites
2527 * depend on local and/or non-local particles in assignVsitesToThread.
2531 for (int i = 0; i < mdatoms.nr; i++)
2533 if (mdatoms.ptype[i] == eptVSite)
2535 /* vsites are not assigned to a task yet */
2540 /* assign non-vsite particles to task thread */
2541 taskIndex_[i] = thread;
2543 if (i == (thread + 1) * natperthread && thread < numThreads_)
2550 #pragma omp parallel num_threads(numThreads_)
2554 int thread = gmx_omp_get_thread_num();
2555 VsiteThread& tData = *tData_[thread];
2557 /* Clear the buffer use flags that were set before */
2558 if (tData.useInterdependentTask)
2560 InterdependentTask& idTask = tData.idTask;
2562 /* To avoid an extra OpenMP barrier in spread_vsite_f,
2563 * we clear the force buffer at the next step,
2564 * so we need to do it here as well.
2566 clearTaskForceBufferUsedElements(&idTask);
2568 idTask.vsite.resize(0);
2569 for (int t = 0; t < numThreads_; t++)
2571 AtomIndex& atomIndex = idTask.atomIndex[t];
2572 int natom = atomIndex.atom.size();
2573 for (int i = 0; i < natom; i++)
2575 idTask.use[atomIndex.atom[i]] = false;
2577 atomIndex.atom.resize(0);
2582 /* To avoid large f_buf allocations of #threads*vsite_atom_range
2583 * we don't use task2 with more than 200000 atoms. This doesn't
2584 * affect performance, since with such a large range relatively few
2585 * vsites will end up in the separate task.
2586 * Note that useTask2 should be the same for all threads.
2588 tData.useInterdependentTask = (vsite_atom_range <= 200000);
2589 if (tData.useInterdependentTask)
2591 size_t natoms_use_in_vsites = vsite_atom_range;
2592 InterdependentTask& idTask = tData.idTask;
2593 /* To avoid resizing and re-clearing every nstlist steps,
2594 * we never down size the force buffer.
2596 if (natoms_use_in_vsites > idTask.force.size() || natoms_use_in_vsites > idTask.use.size())
2598 idTask.force.resize(natoms_use_in_vsites, { 0, 0, 0 });
2599 idTask.use.resize(natoms_use_in_vsites, false);
2603 /* Assign all vsites that can execute independently on threads */
2604 tData.rangeStart = thread * natperthread;
2605 if (thread < numThreads_ - 1)
2607 tData.rangeEnd = (thread + 1) * natperthread;
2611 /* The last thread should cover up to the end of the range */
2612 tData.rangeEnd = mdatoms.nr;
2614 assignVsitesToThread(&tData, thread, numThreads_, natperthread, taskIndex_, ilists,
2615 iparams, mdatoms.ptype);
2617 if (tData.useInterdependentTask)
2619 /* In the worst case, all tasks write to force ranges of
2620 * all other tasks, leading to #tasks^2 scaling (this is only
2621 * the overhead, the actual flops remain constant).
2622 * But in most cases there is far less coupling. To improve
2623 * scaling at high thread counts we therefore construct
2624 * an index to only loop over the actually affected tasks.
2626 InterdependentTask& idTask = tData.idTask;
2628 /* Ensure assignVsitesToThread finished on other threads */
2631 idTask.spreadTask.resize(0);
2632 idTask.reduceTask.resize(0);
2633 for (int t = 0; t < numThreads_; t++)
2635 /* Do we write to the force buffer of task t? */
2636 if (!idTask.atomIndex[t].atom.empty())
2638 idTask.spreadTask.push_back(t);
2640 /* Does task t write to our force buffer? */
2641 if (!tData_[t]->idTask.atomIndex[thread].atom.empty())
2643 idTask.reduceTask.push_back(t);
2648 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2650 /* Assign all remaining vsites, that will have taskIndex[]=2*vsite->nthreads,
2651 * to a single task that will not run in parallel with other tasks.
2653 assignVsitesToSingleTask(tData_[numThreads_].get(), 2 * numThreads_, taskIndex_, ilists, iparams);
2655 if (debug && numThreads_ > 1)
2657 fprintf(debug, "virtual site useInterdependentTask %d, nuse:\n",
2658 static_cast<int>(tData_[0]->useInterdependentTask));
2659 for (int th = 0; th < numThreads_ + 1; th++)
2661 fprintf(debug, " %4d", tData_[th]->idTask.nuse);
2663 fprintf(debug, "\n");
2665 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2667 if (!ilists[ftype].empty())
2669 fprintf(debug, "%-20s thread dist:", interaction_function[ftype].longname);
2670 for (int th = 0; th < numThreads_ + 1; th++)
2672 fprintf(debug, " %4d %4d ", tData_[th]->ilist[ftype].size(),
2673 tData_[th]->idTask.ilist[ftype].size());
2675 fprintf(debug, "\n");
2681 int nrOrig = vsiteIlistNrCount(ilists);
2683 for (int th = 0; th < numThreads_ + 1; th++)
2685 nrThreaded += vsiteIlistNrCount(tData_[th]->ilist) + vsiteIlistNrCount(tData_[th]->idTask.ilist);
2687 GMX_ASSERT(nrThreaded == nrOrig,
2688 "The number of virtual sites assigned to all thread task has to match the total "
2689 "number of virtual sites");
2693 void VirtualSitesHandler::Impl::setVirtualSites(ArrayRef<const InteractionList> ilists,
2694 const t_mdatoms& mdatoms)
2698 threadingInfo_.setVirtualSites(ilists, iparams_, mdatoms, domainInfo_.useDomdec());
2701 void VirtualSitesHandler::setVirtualSites(ArrayRef<const InteractionList> ilists, const t_mdatoms& mdatoms)
2703 impl_->setVirtualSites(ilists, mdatoms);