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_vsite2(const rvec xi, const rvec xj, rvec x, real a, const t_pbc* pbc)
345 pbc_dx_aiuc(pbc, xj, xi, dx);
346 x[XX] = xi[XX] + a * dx[XX];
347 x[YY] = xi[YY] + a * dx[YY];
348 x[ZZ] = xi[ZZ] + a * dx[ZZ];
352 x[XX] = b * xi[XX] + a * xj[XX];
353 x[YY] = b * xi[YY] + a * xj[YY];
354 x[ZZ] = b * xi[ZZ] + a * xj[ZZ];
358 /* TOTAL: 10 flops */
361 static void constr_vsite2FD(const rvec xi, const rvec xj, rvec x, real a, const t_pbc* pbc)
364 pbc_rvec_sub(pbc, xj, xi, xij);
367 const real b = a * inverseNorm(xij);
370 x[XX] = xi[XX] + b * xij[XX];
371 x[YY] = xi[YY] + b * xij[YY];
372 x[ZZ] = xi[ZZ] + b * xij[ZZ];
375 /* TOTAL: 25 flops */
378 static void constr_vsite3(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, const t_pbc* pbc)
387 pbc_dx_aiuc(pbc, xj, xi, dxj);
388 pbc_dx_aiuc(pbc, xk, xi, dxk);
389 x[XX] = xi[XX] + a * dxj[XX] + b * dxk[XX];
390 x[YY] = xi[YY] + a * dxj[YY] + b * dxk[YY];
391 x[ZZ] = xi[ZZ] + a * dxj[ZZ] + b * dxk[ZZ];
395 x[XX] = c * xi[XX] + a * xj[XX] + b * xk[XX];
396 x[YY] = c * xi[YY] + a * xj[YY] + b * xk[YY];
397 x[ZZ] = c * xi[ZZ] + a * xj[ZZ] + b * xk[ZZ];
401 /* TOTAL: 17 flops */
404 static void constr_vsite3FD(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, const t_pbc* pbc)
409 pbc_rvec_sub(pbc, xj, xi, xij);
410 pbc_rvec_sub(pbc, xk, xj, xjk);
413 /* temp goes from i to a point on the line jk */
414 temp[XX] = xij[XX] + a * xjk[XX];
415 temp[YY] = xij[YY] + a * xjk[YY];
416 temp[ZZ] = xij[ZZ] + a * xjk[ZZ];
419 c = b * inverseNorm(temp);
422 x[XX] = xi[XX] + c * temp[XX];
423 x[YY] = xi[YY] + c * temp[YY];
424 x[ZZ] = xi[ZZ] + c * temp[ZZ];
427 /* TOTAL: 34 flops */
430 static void constr_vsite3FAD(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, const t_pbc* pbc)
433 real a1, b1, c1, invdij;
435 pbc_rvec_sub(pbc, xj, xi, xij);
436 pbc_rvec_sub(pbc, xk, xj, xjk);
439 invdij = inverseNorm(xij);
440 c1 = invdij * invdij * iprod(xij, xjk);
441 xp[XX] = xjk[XX] - c1 * xij[XX];
442 xp[YY] = xjk[YY] - c1 * xij[YY];
443 xp[ZZ] = xjk[ZZ] - c1 * xij[ZZ];
445 b1 = b * inverseNorm(xp);
448 x[XX] = xi[XX] + a1 * xij[XX] + b1 * xp[XX];
449 x[YY] = xi[YY] + a1 * xij[YY] + b1 * xp[YY];
450 x[ZZ] = xi[ZZ] + a1 * xij[ZZ] + b1 * xp[ZZ];
453 /* TOTAL: 63 flops */
457 constr_vsite3OUT(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, real c, const t_pbc* pbc)
461 pbc_rvec_sub(pbc, xj, xi, xij);
462 pbc_rvec_sub(pbc, xk, xi, xik);
463 cprod(xij, xik, temp);
466 x[XX] = xi[XX] + a * xij[XX] + b * xik[XX] + c * temp[XX];
467 x[YY] = xi[YY] + a * xij[YY] + b * xik[YY] + c * temp[YY];
468 x[ZZ] = xi[ZZ] + a * xij[ZZ] + b * xik[ZZ] + c * temp[ZZ];
471 /* TOTAL: 33 flops */
474 static void constr_vsite4FD(const rvec xi,
484 rvec xij, xjk, xjl, temp;
487 pbc_rvec_sub(pbc, xj, xi, xij);
488 pbc_rvec_sub(pbc, xk, xj, xjk);
489 pbc_rvec_sub(pbc, xl, xj, xjl);
492 /* temp goes from i to a point on the plane jkl */
493 temp[XX] = xij[XX] + a * xjk[XX] + b * xjl[XX];
494 temp[YY] = xij[YY] + a * xjk[YY] + b * xjl[YY];
495 temp[ZZ] = xij[ZZ] + a * xjk[ZZ] + b * xjl[ZZ];
498 d = c * inverseNorm(temp);
501 x[XX] = xi[XX] + d * temp[XX];
502 x[YY] = xi[YY] + d * temp[YY];
503 x[ZZ] = xi[ZZ] + d * temp[ZZ];
506 /* TOTAL: 43 flops */
509 static void constr_vsite4FDN(const rvec xi,
519 rvec xij, xik, xil, ra, rb, rja, rjb, rm;
522 pbc_rvec_sub(pbc, xj, xi, xij);
523 pbc_rvec_sub(pbc, xk, xi, xik);
524 pbc_rvec_sub(pbc, xl, xi, xil);
527 ra[XX] = a * xik[XX];
528 ra[YY] = a * xik[YY];
529 ra[ZZ] = a * xik[ZZ];
531 rb[XX] = b * xil[XX];
532 rb[YY] = b * xil[YY];
533 rb[ZZ] = b * xil[ZZ];
537 rvec_sub(ra, xij, rja);
538 rvec_sub(rb, xij, rjb);
544 d = c * inverseNorm(rm);
547 x[XX] = xi[XX] + d * rm[XX];
548 x[YY] = xi[YY] + d * rm[YY];
549 x[ZZ] = xi[ZZ] + d * rm[ZZ];
552 /* TOTAL: 47 flops */
555 static int constr_vsiten(const t_iatom* ia, ArrayRef<const t_iparams> ip, ArrayRef<RVec> x, const t_pbc* pbc)
562 n3 = 3 * ip[ia[0]].vsiten.n;
565 copy_rvec(x[ai], x1);
567 for (int i = 3; i < n3; i += 3)
570 a = ip[ia[i]].vsiten.a;
573 pbc_dx_aiuc(pbc, x[ai], x1, dx);
577 rvec_sub(x[ai], x1, dx);
579 dsum[XX] += a * dx[XX];
580 dsum[YY] += a * dx[YY];
581 dsum[ZZ] += a * dx[ZZ];
585 x[av][XX] = x1[XX] + dsum[XX];
586 x[av][YY] = x1[YY] + dsum[YY];
587 x[av][ZZ] = x1[ZZ] + dsum[ZZ];
594 //! PBC modes for vsite construction and spreading
597 all, //!< Apply normal, simple PBC for all vsites
598 none //!< No PBC treatment needed
601 /*! \brief Returns the PBC mode based on the system PBC and vsite properties
603 * \param[in] pbcPtr A pointer to a PBC struct or nullptr when no PBC treatment is required
605 static PbcMode getPbcMode(const t_pbc* pbcPtr)
607 if (pbcPtr == nullptr)
609 return PbcMode::none;
617 /*! \brief Executes the vsite construction task for a single thread
619 * \param[in,out] x Coordinates to construct vsites for
620 * \param[in] dt Time step, needed when v is not empty
621 * \param[in,out] v When not empty, velocities are generated for virtual sites
622 * \param[in] ip Interaction parameters for all interaction, only vsite parameters are used
623 * \param[in] ilist The interaction lists, only vsites are usesd
624 * \param[in] pbc_null PBC struct, used for PBC distance calculations when !=nullptr
626 static void construct_vsites_thread(ArrayRef<RVec> x,
629 ArrayRef<const t_iparams> ip,
630 ArrayRef<const InteractionList> ilist,
631 const t_pbc* pbc_null)
643 const PbcMode pbcMode = getPbcMode(pbc_null);
644 /* We need another pbc pointer, as with charge groups we switch per vsite */
645 const t_pbc* pbc_null2 = pbc_null;
647 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
649 if (ilist[ftype].empty())
655 int nra = interaction_function[ftype].nratoms;
657 int nr = ilist[ftype].size();
659 const t_iatom* ia = ilist[ftype].iatoms.data();
661 for (int i = 0; i < nr;)
664 /* The vsite and constructing atoms */
667 /* Constants for constructing vsites */
668 real a1 = ip[tp].vsite.a;
669 /* Copy the old position */
671 copy_rvec(x[avsite], xv);
673 /* Construct the vsite depending on type */
680 constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
684 constr_vsite2FD(x[ai], x[aj], x[avsite], a1, pbc_null2);
690 constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
696 constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
702 constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
709 constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
717 constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1, pbc_null2);
725 constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1, pbc_null2);
727 case F_VSITEN: inc = constr_vsiten(ia, ip, x, pbc_null2); break;
729 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d", ftype, __FILE__, __LINE__);
732 if (pbcMode == PbcMode::all)
734 /* Keep the vsite in the same periodic image as before */
736 int ishift = pbc_dx_aiuc(pbc_null, x[avsite], xv, dx);
737 if (ishift != CENTRAL)
739 rvec_add(xv, dx, x[avsite]);
744 /* Calculate velocity of vsite... */
746 rvec_sub(x[avsite], xv, vv);
747 svmul(inv_dt, vv, v[avsite]);
750 /* Increment loop variables */
758 /*! \brief Dispatch the vsite construction tasks for all threads
760 * \param[in] threadingInfo Used to divide work over threads when != nullptr
761 * \param[in,out] x Coordinates to construct vsites for
762 * \param[in] dt Time step, needed when v is not empty
763 * \param[in,out] v When not empty, velocities are generated for virtual sites
764 * \param[in] ip Interaction parameters for all interaction, only vsite parameters are used
765 * \param[in] ilist The interaction lists, only vsites are usesd
766 * \param[in] domainInfo Information about PBC and DD
767 * \param[in] box Used for PBC when PBC is set in domainInfo
769 static void construct_vsites(const ThreadingInfo* threadingInfo,
773 ArrayRef<const t_iparams> ip,
774 ArrayRef<const InteractionList> ilist,
775 const DomainInfo& domainInfo,
778 const bool useDomdec = domainInfo.useDomdec();
780 t_pbc pbc, *pbc_null;
782 /* We only need to do pbc when we have inter update-group vsites.
783 * Note that with domain decomposition we do not need to apply PBC here
784 * when we have at least 3 domains along each dimension. Currently we
785 * do not optimize this case.
787 if (domainInfo.pbcType_ != PbcType::No && domainInfo.useMolPbc_)
789 /* This is wasting some CPU time as we now do this multiple times
793 clear_ivec(null_ivec);
794 pbc_null = set_pbc_dd(&pbc, domainInfo.pbcType_,
795 useDomdec ? domainInfo.domdec_->numCells : null_ivec, FALSE, box);
804 dd_move_x_vsites(*domainInfo.domdec_, box, as_rvec_array(x.data()));
807 if (threadingInfo == nullptr || threadingInfo->numThreads() == 1)
809 construct_vsites_thread(x, dt, v, ip, ilist, pbc_null);
813 #pragma omp parallel num_threads(threadingInfo->numThreads())
817 const int th = gmx_omp_get_thread_num();
818 const VsiteThread& tData = threadingInfo->threadData(th);
819 GMX_ASSERT(tData.rangeStart >= 0,
820 "The thread data should be initialized before calling construct_vsites");
822 construct_vsites_thread(x, dt, v, ip, tData.ilist, pbc_null);
823 if (tData.useInterdependentTask)
825 /* Here we don't need a barrier (unlike the spreading),
826 * since both tasks only construct vsites from particles,
827 * or local vsites, not from non-local vsites.
829 construct_vsites_thread(x, dt, v, ip, tData.idTask.ilist, pbc_null);
832 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
834 /* Now we can construct the vsites that might depend on other vsites */
835 construct_vsites_thread(x, dt, v, ip, threadingInfo->threadDataNonLocalDependent().ilist, pbc_null);
839 void VirtualSitesHandler::Impl::construct(ArrayRef<RVec> x, real dt, ArrayRef<RVec> v, const matrix box) const
841 construct_vsites(&threadingInfo_, x, dt, v, iparams_, ilists_, domainInfo_, box);
844 void VirtualSitesHandler::construct(ArrayRef<RVec> x, real dt, ArrayRef<RVec> v, const matrix box) const
846 impl_->construct(x, dt, v, box);
849 void constructVirtualSites(ArrayRef<RVec> x, ArrayRef<const t_iparams> ip, ArrayRef<const InteractionList> ilist)
853 const DomainInfo domainInfo;
854 construct_vsites(nullptr, x, 0, {}, ip, ilist, domainInfo, nullptr);
858 /* Force spreading routines */
860 template<VirialHandling virialHandling>
861 static void spread_vsite2(const t_iatom ia[],
863 ArrayRef<const RVec> x,
865 ArrayRef<RVec> fshift,
875 svmul(1 - a, f[av], fi);
883 if (virialHandling == VirialHandling::Pbc)
889 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
890 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
898 if (siv != CENTRAL || sij != CENTRAL)
900 rvec_inc(fshift[siv], f[av]);
901 rvec_dec(fshift[CENTRAL], fi);
902 rvec_dec(fshift[sij], fj);
906 /* TOTAL: 13 flops */
909 void constructVirtualSitesGlobal(const gmx_mtop_t& mtop, gmx::ArrayRef<gmx::RVec> x)
911 GMX_ASSERT(x.ssize() >= mtop.natoms, "x should contain the whole system");
912 GMX_ASSERT(!mtop.moleculeBlockIndices.empty(),
913 "molblock indices are needed in constructVsitesGlobal");
915 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
917 const gmx_molblock_t& molb = mtop.molblock[mb];
918 const gmx_moltype_t& molt = mtop.moltype[molb.type];
919 if (vsiteIlistNrCount(molt.ilist) > 0)
921 int atomOffset = mtop.moleculeBlockIndices[mb].globalAtomStart;
922 for (int mol = 0; mol < molb.nmol; mol++)
924 constructVirtualSites(x.subArray(atomOffset, molt.atoms.nr), mtop.ffparams.iparams,
926 atomOffset += molt.atoms.nr;
932 template<VirialHandling virialHandling>
933 static void spread_vsite2FD(const t_iatom ia[],
935 ArrayRef<const RVec> x,
937 ArrayRef<RVec> fshift,
941 const int av = ia[1];
942 const int ai = ia[2];
943 const int aj = ia[3];
945 copy_rvec(f[av], fv);
948 int sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
951 const real invDistance = inverseNorm(xij);
952 const real b = a * invDistance;
955 const real fproj = iprod(xij, fv) * invDistance * invDistance;
958 fj[XX] = b * (fv[XX] - fproj * xij[XX]);
959 fj[YY] = b * (fv[YY] - fproj * xij[YY]);
960 fj[ZZ] = b * (fv[ZZ] - fproj * xij[ZZ]);
963 /* b is already calculated in constr_vsite2FD
964 storing b somewhere will save flops. */
966 f[ai][XX] += fv[XX] - fj[XX];
967 f[ai][YY] += fv[YY] - fj[YY];
968 f[ai][ZZ] += fv[ZZ] - fj[ZZ];
974 if (virialHandling == VirialHandling::Pbc)
980 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
987 if (svi != CENTRAL || sji != CENTRAL)
989 rvec_dec(fshift[svi], fv);
990 fshift[CENTRAL][XX] += fv[XX] - fj[XX];
991 fshift[CENTRAL][YY] += fv[YY] - fj[YY];
992 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ];
993 fshift[sji][XX] += fj[XX];
994 fshift[sji][YY] += fj[YY];
995 fshift[sji][ZZ] += fj[ZZ];
999 if (virialHandling == VirialHandling::NonLinear)
1001 /* Under this condition, the virial for the current forces is not
1002 * calculated from the redistributed forces. This means that
1003 * the effect of non-linear virtual site constructions on the virial
1004 * needs to be added separately. This contribution can be calculated
1005 * in many ways, but the simplest and cheapest way is to use
1006 * the first constructing atom ai as a reference position in space:
1007 * subtract (xv-xi)*fv and add (xj-xi)*fj.
1011 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1013 for (int i = 0; i < DIM; i++)
1015 for (int j = 0; j < DIM; j++)
1017 /* As xix is a linear combination of j and k, use that here */
1018 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * fj[j];
1023 /* TOTAL: 38 flops */
1026 template<VirialHandling virialHandling>
1027 static void spread_vsite3(const t_iatom ia[],
1030 ArrayRef<const RVec> x,
1032 ArrayRef<RVec> fshift,
1035 rvec fi, fj, fk, dx;
1043 svmul(1 - a - b, f[av], fi);
1044 svmul(a, f[av], fj);
1045 svmul(b, f[av], fk);
1048 rvec_inc(f[ai], fi);
1049 rvec_inc(f[aj], fj);
1050 rvec_inc(f[ak], fk);
1053 if (virialHandling == VirialHandling::Pbc)
1060 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
1061 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
1062 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
1071 if (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL)
1073 rvec_inc(fshift[siv], f[av]);
1074 rvec_dec(fshift[CENTRAL], fi);
1075 rvec_dec(fshift[sij], fj);
1076 rvec_dec(fshift[sik], fk);
1080 /* TOTAL: 20 flops */
1083 template<VirialHandling virialHandling>
1084 static void spread_vsite3FD(const t_iatom ia[],
1087 ArrayRef<const RVec> x,
1089 ArrayRef<RVec> fshift,
1094 rvec xvi, xij, xjk, xix, fv, temp;
1095 t_iatom av, ai, aj, ak;
1102 copy_rvec(f[av], fv);
1104 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1105 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1108 /* xix goes from i to point x on the line jk */
1109 xix[XX] = xij[XX] + a * xjk[XX];
1110 xix[YY] = xij[YY] + a * xjk[YY];
1111 xix[ZZ] = xij[ZZ] + a * xjk[ZZ];
1114 const real invDistance = inverseNorm(xix);
1115 const real c = b * invDistance;
1116 /* 4 + ?10? flops */
1118 fproj = iprod(xix, fv) * invDistance * invDistance; /* = (xix . f)/(xix . xix) */
1120 temp[XX] = c * (fv[XX] - fproj * xix[XX]);
1121 temp[YY] = c * (fv[YY] - fproj * xix[YY]);
1122 temp[ZZ] = c * (fv[ZZ] - fproj * xix[ZZ]);
1125 /* c is already calculated in constr_vsite3FD
1126 storing c somewhere will save 26 flops! */
1129 f[ai][XX] += fv[XX] - temp[XX];
1130 f[ai][YY] += fv[YY] - temp[YY];
1131 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
1132 f[aj][XX] += a1 * temp[XX];
1133 f[aj][YY] += a1 * temp[YY];
1134 f[aj][ZZ] += a1 * temp[ZZ];
1135 f[ak][XX] += a * temp[XX];
1136 f[ak][YY] += a * temp[YY];
1137 f[ak][ZZ] += a * temp[ZZ];
1140 if (virialHandling == VirialHandling::Pbc)
1145 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1152 if (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL)
1154 rvec_dec(fshift[svi], fv);
1155 fshift[CENTRAL][XX] += fv[XX] - (1 + a) * temp[XX];
1156 fshift[CENTRAL][YY] += fv[YY] - (1 + a) * temp[YY];
1157 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a) * temp[ZZ];
1158 fshift[sji][XX] += temp[XX];
1159 fshift[sji][YY] += temp[YY];
1160 fshift[sji][ZZ] += temp[ZZ];
1161 fshift[skj][XX] += a * temp[XX];
1162 fshift[skj][YY] += a * temp[YY];
1163 fshift[skj][ZZ] += a * temp[ZZ];
1167 if (virialHandling == VirialHandling::NonLinear)
1169 /* Under this condition, the virial for the current forces is not
1170 * calculated from the redistributed forces. This means that
1171 * the effect of non-linear virtual site constructions on the virial
1172 * needs to be added separately. This contribution can be calculated
1173 * in many ways, but the simplest and cheapest way is to use
1174 * the first constructing atom ai as a reference position in space:
1175 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
1179 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1181 for (int i = 0; i < DIM; i++)
1183 for (int j = 0; j < DIM; j++)
1185 /* As xix is a linear combination of j and k, use that here */
1186 dxdf[i][j] += -xiv[i] * fv[j] + xix[i] * temp[j];
1191 /* TOTAL: 61 flops */
1194 template<VirialHandling virialHandling>
1195 static void spread_vsite3FAD(const t_iatom ia[],
1198 ArrayRef<const RVec> x,
1200 ArrayRef<RVec> fshift,
1204 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
1205 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
1206 t_iatom av, ai, aj, ak;
1213 copy_rvec(f[ia[1]], fv);
1215 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1216 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1219 invdij = inverseNorm(xij);
1220 invdij2 = invdij * invdij;
1221 c1 = iprod(xij, xjk) * invdij2;
1222 xperp[XX] = xjk[XX] - c1 * xij[XX];
1223 xperp[YY] = xjk[YY] - c1 * xij[YY];
1224 xperp[ZZ] = xjk[ZZ] - c1 * xij[ZZ];
1225 /* xperp in plane ijk, perp. to ij */
1226 invdp = inverseNorm(xperp);
1231 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
1232 storing them somewhere will save 45 flops! */
1234 fproj = iprod(xij, fv) * invdij2;
1235 svmul(fproj, xij, Fpij); /* proj. f on xij */
1236 svmul(iprod(xperp, fv) * invdp * invdp, xperp, Fppp); /* proj. f on xperp */
1237 svmul(b1 * fproj, xperp, f3);
1240 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
1241 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
1242 for (int d = 0; d < DIM; d++)
1250 f[ai][XX] += fv[XX] - f1[XX] + c1 * f2[XX] + f3[XX];
1251 f[ai][YY] += fv[YY] - f1[YY] + c1 * f2[YY] + f3[YY];
1252 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1 * f2[ZZ] + f3[ZZ];
1253 f[aj][XX] += f1[XX] - c2 * f2[XX] - f3[XX];
1254 f[aj][YY] += f1[YY] - c2 * f2[YY] - f3[YY];
1255 f[aj][ZZ] += f1[ZZ] - c2 * f2[ZZ] - f3[ZZ];
1256 f[ak][XX] += f2[XX];
1257 f[ak][YY] += f2[YY];
1258 f[ak][ZZ] += f2[ZZ];
1261 if (virialHandling == VirialHandling::Pbc)
1267 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1274 if (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL)
1276 rvec_dec(fshift[svi], fv);
1277 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1 - c1) * f2[XX] + f3[XX];
1278 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1 - c1) * f2[YY] + f3[YY];
1279 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1 - c1) * f2[ZZ] + f3[ZZ];
1280 fshift[sji][XX] += f1[XX] - c1 * f2[XX] - f3[XX];
1281 fshift[sji][YY] += f1[YY] - c1 * f2[YY] - f3[YY];
1282 fshift[sji][ZZ] += f1[ZZ] - c1 * f2[ZZ] - f3[ZZ];
1283 fshift[skj][XX] += f2[XX];
1284 fshift[skj][YY] += f2[YY];
1285 fshift[skj][ZZ] += f2[ZZ];
1289 if (virialHandling == VirialHandling::NonLinear)
1292 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1294 for (int i = 0; i < DIM; i++)
1296 for (int j = 0; j < DIM; j++)
1298 /* Note that xik=xij+xjk, so we have to add xij*f2 */
1299 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * (f1[j] + (1 - c2) * f2[j] - f3[j])
1305 /* TOTAL: 113 flops */
1308 template<VirialHandling virialHandling>
1309 static void spread_vsite3OUT(const t_iatom ia[],
1313 ArrayRef<const RVec> x,
1315 ArrayRef<RVec> fshift,
1319 rvec xvi, xij, xik, fv, fj, fk;
1329 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1330 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1333 copy_rvec(f[av], fv);
1340 fj[XX] = a * fv[XX] - xik[ZZ] * cfy + xik[YY] * cfz;
1341 fj[YY] = xik[ZZ] * cfx + a * fv[YY] - xik[XX] * cfz;
1342 fj[ZZ] = -xik[YY] * cfx + xik[XX] * cfy + a * fv[ZZ];
1344 fk[XX] = b * fv[XX] + xij[ZZ] * cfy - xij[YY] * cfz;
1345 fk[YY] = -xij[ZZ] * cfx + b * fv[YY] + xij[XX] * cfz;
1346 fk[ZZ] = xij[YY] * cfx - xij[XX] * cfy + b * fv[ZZ];
1349 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
1350 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
1351 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1352 rvec_inc(f[aj], fj);
1353 rvec_inc(f[ak], fk);
1356 if (virialHandling == VirialHandling::Pbc)
1361 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1368 if (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL)
1370 rvec_dec(fshift[svi], fv);
1371 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
1372 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
1373 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1374 rvec_inc(fshift[sji], fj);
1375 rvec_inc(fshift[ski], fk);
1379 if (virialHandling == VirialHandling::NonLinear)
1383 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1385 for (int i = 0; i < DIM; i++)
1387 for (int j = 0; j < DIM; j++)
1389 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * fj[j] + xik[i] * fk[j];
1394 /* TOTAL: 54 flops */
1397 template<VirialHandling virialHandling>
1398 static void spread_vsite4FD(const t_iatom ia[],
1402 ArrayRef<const RVec> x,
1404 ArrayRef<RVec> fshift,
1409 rvec xvi, xij, xjk, xjl, xix, fv, temp;
1410 int av, ai, aj, ak, al;
1411 int sji, skj, slj, m;
1419 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1420 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1421 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1424 /* xix goes from i to point x on the plane jkl */
1425 for (m = 0; m < DIM; m++)
1427 xix[m] = xij[m] + a * xjk[m] + b * xjl[m];
1431 const real invDistance = inverseNorm(xix);
1432 const real d = c * invDistance;
1433 /* 4 + ?10? flops */
1435 copy_rvec(f[av], fv);
1437 fproj = iprod(xix, fv) * invDistance * invDistance; /* = (xix . f)/(xix . xix) */
1439 for (m = 0; m < DIM; m++)
1441 temp[m] = d * (fv[m] - fproj * xix[m]);
1445 /* c is already calculated in constr_vsite3FD
1446 storing c somewhere will save 35 flops! */
1449 for (m = 0; m < DIM; m++)
1451 f[ai][m] += fv[m] - temp[m];
1452 f[aj][m] += a1 * temp[m];
1453 f[ak][m] += a * temp[m];
1454 f[al][m] += b * temp[m];
1458 if (virialHandling == VirialHandling::Pbc)
1463 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1470 if (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL)
1472 rvec_dec(fshift[svi], fv);
1473 for (m = 0; m < DIM; m++)
1475 fshift[CENTRAL][m] += fv[m] - (1 + a + b) * temp[m];
1476 fshift[sji][m] += temp[m];
1477 fshift[skj][m] += a * temp[m];
1478 fshift[slj][m] += b * temp[m];
1483 if (virialHandling == VirialHandling::NonLinear)
1488 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1490 for (i = 0; i < DIM; i++)
1492 for (j = 0; j < DIM; j++)
1494 dxdf[i][j] += -xiv[i] * fv[j] + xix[i] * temp[j];
1499 /* TOTAL: 77 flops */
1502 template<VirialHandling virialHandling>
1503 static void spread_vsite4FDN(const t_iatom ia[],
1507 ArrayRef<const RVec> x,
1509 ArrayRef<RVec> fshift,
1513 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1514 rvec fv, fj, fk, fl;
1517 int av, ai, aj, ak, al;
1520 /* DEBUG: check atom indices */
1527 copy_rvec(f[av], fv);
1529 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1530 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1531 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1534 ra[XX] = a * xik[XX];
1535 ra[YY] = a * xik[YY];
1536 ra[ZZ] = a * xik[ZZ];
1538 rb[XX] = b * xil[XX];
1539 rb[YY] = b * xil[YY];
1540 rb[ZZ] = b * xil[ZZ];
1544 rvec_sub(ra, xij, rja);
1545 rvec_sub(rb, xij, rjb);
1546 rvec_sub(rb, ra, rab);
1549 cprod(rja, rjb, rm);
1552 invrm = inverseNorm(rm);
1553 denom = invrm * invrm;
1556 cfx = c * invrm * fv[XX];
1557 cfy = c * invrm * fv[YY];
1558 cfz = c * invrm * fv[ZZ];
1569 fj[XX] = (-rm[XX] * rt[XX]) * cfx + (rab[ZZ] - rm[YY] * rt[XX]) * cfy
1570 + (-rab[YY] - rm[ZZ] * rt[XX]) * cfz;
1571 fj[YY] = (-rab[ZZ] - rm[XX] * rt[YY]) * cfx + (-rm[YY] * rt[YY]) * cfy
1572 + (rab[XX] - rm[ZZ] * rt[YY]) * cfz;
1573 fj[ZZ] = (rab[YY] - rm[XX] * rt[ZZ]) * cfx + (-rab[XX] - rm[YY] * rt[ZZ]) * cfy
1574 + (-rm[ZZ] * rt[ZZ]) * cfz;
1580 rt[XX] *= denom * a;
1581 rt[YY] *= denom * a;
1582 rt[ZZ] *= denom * a;
1585 fk[XX] = (-rm[XX] * rt[XX]) * cfx + (-a * rjb[ZZ] - rm[YY] * rt[XX]) * cfy
1586 + (a * rjb[YY] - rm[ZZ] * rt[XX]) * cfz;
1587 fk[YY] = (a * rjb[ZZ] - rm[XX] * rt[YY]) * cfx + (-rm[YY] * rt[YY]) * cfy
1588 + (-a * rjb[XX] - rm[ZZ] * rt[YY]) * cfz;
1589 fk[ZZ] = (-a * rjb[YY] - rm[XX] * rt[ZZ]) * cfx + (a * rjb[XX] - rm[YY] * rt[ZZ]) * cfy
1590 + (-rm[ZZ] * rt[ZZ]) * cfz;
1596 rt[XX] *= denom * b;
1597 rt[YY] *= denom * b;
1598 rt[ZZ] *= denom * b;
1601 fl[XX] = (-rm[XX] * rt[XX]) * cfx + (b * rja[ZZ] - rm[YY] * rt[XX]) * cfy
1602 + (-b * rja[YY] - rm[ZZ] * rt[XX]) * cfz;
1603 fl[YY] = (-b * rja[ZZ] - rm[XX] * rt[YY]) * cfx + (-rm[YY] * rt[YY]) * cfy
1604 + (b * rja[XX] - rm[ZZ] * rt[YY]) * cfz;
1605 fl[ZZ] = (b * rja[YY] - rm[XX] * rt[ZZ]) * cfx + (-b * rja[XX] - rm[YY] * rt[ZZ]) * cfy
1606 + (-rm[ZZ] * rt[ZZ]) * cfz;
1609 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1610 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1611 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1612 rvec_inc(f[aj], fj);
1613 rvec_inc(f[ak], fk);
1614 rvec_inc(f[al], fl);
1617 if (virialHandling == VirialHandling::Pbc)
1622 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1629 if (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL)
1631 rvec_dec(fshift[svi], fv);
1632 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1633 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1634 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1635 rvec_inc(fshift[sij], fj);
1636 rvec_inc(fshift[sik], fk);
1637 rvec_inc(fshift[sil], fl);
1641 if (virialHandling == VirialHandling::NonLinear)
1646 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1648 for (i = 0; i < DIM; i++)
1650 for (j = 0; j < DIM; j++)
1652 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * fj[j] + xik[i] * fk[j] + xil[i] * fl[j];
1657 /* Total: 207 flops (Yuck!) */
1660 template<VirialHandling virialHandling>
1661 static int spread_vsiten(const t_iatom ia[],
1662 ArrayRef<const t_iparams> ip,
1663 ArrayRef<const RVec> x,
1665 ArrayRef<RVec> fshift,
1673 n3 = 3 * ip[ia[0]].vsiten.n;
1675 copy_rvec(x[av], xv);
1677 for (i = 0; i < n3; i += 3)
1682 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1688 a = ip[ia[i]].vsiten.a;
1689 svmul(a, f[av], fi);
1690 rvec_inc(f[ai], fi);
1692 if (virialHandling == VirialHandling::Pbc && siv != CENTRAL)
1694 rvec_inc(fshift[siv], fi);
1695 rvec_dec(fshift[CENTRAL], fi);
1705 //! Returns the number of virtual sites in the interaction list, for VSITEN the number of atoms
1706 static int vsite_count(ArrayRef<const InteractionList> ilist, int ftype)
1708 if (ftype == F_VSITEN)
1710 return ilist[ftype].size() / 3;
1714 return ilist[ftype].size() / (1 + interaction_function[ftype].nratoms);
1718 //! Executes the force spreading task for a single thread
1719 template<VirialHandling virialHandling>
1720 static void spreadForceForThread(ArrayRef<const RVec> x,
1722 ArrayRef<RVec> fshift,
1724 ArrayRef<const t_iparams> ip,
1725 ArrayRef<const InteractionList> ilist,
1726 const t_pbc* pbc_null)
1728 const PbcMode pbcMode = getPbcMode(pbc_null);
1729 /* We need another pbc pointer, as with charge groups we switch per vsite */
1730 const t_pbc* pbc_null2 = pbc_null;
1731 gmx::ArrayRef<const int> vsite_pbc;
1733 /* this loop goes backwards to be able to build *
1734 * higher type vsites from lower types */
1735 for (int ftype = c_ftypeVsiteEnd - 1; ftype >= c_ftypeVsiteStart; ftype--)
1737 if (ilist[ftype].empty())
1743 int nra = interaction_function[ftype].nratoms;
1745 int nr = ilist[ftype].size();
1747 const t_iatom* ia = ilist[ftype].iatoms.data();
1749 if (pbcMode == PbcMode::all)
1751 pbc_null2 = pbc_null;
1754 for (int i = 0; i < nr;)
1758 /* Constants for constructing */
1760 a1 = ip[tp].vsite.a;
1761 /* Construct the vsite depending on type */
1765 spread_vsite2<virialHandling>(ia, a1, x, f, fshift, pbc_null2);
1768 spread_vsite2FD<virialHandling>(ia, a1, x, f, fshift, dxdf, pbc_null2);
1771 b1 = ip[tp].vsite.b;
1772 spread_vsite3<virialHandling>(ia, a1, b1, x, f, fshift, pbc_null2);
1775 b1 = ip[tp].vsite.b;
1776 spread_vsite3FD<virialHandling>(ia, a1, b1, x, f, fshift, dxdf, pbc_null2);
1779 b1 = ip[tp].vsite.b;
1780 spread_vsite3FAD<virialHandling>(ia, a1, b1, x, f, fshift, dxdf, pbc_null2);
1783 b1 = ip[tp].vsite.b;
1784 c1 = ip[tp].vsite.c;
1785 spread_vsite3OUT<virialHandling>(ia, a1, b1, c1, x, f, fshift, dxdf, pbc_null2);
1788 b1 = ip[tp].vsite.b;
1789 c1 = ip[tp].vsite.c;
1790 spread_vsite4FD<virialHandling>(ia, a1, b1, c1, x, f, fshift, dxdf, pbc_null2);
1793 b1 = ip[tp].vsite.b;
1794 c1 = ip[tp].vsite.c;
1795 spread_vsite4FDN<virialHandling>(ia, a1, b1, c1, x, f, fshift, dxdf, pbc_null2);
1798 inc = spread_vsiten<virialHandling>(ia, ip, x, f, fshift, pbc_null2);
1801 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d", ftype, __FILE__, __LINE__);
1803 clear_rvec(f[ia[1]]);
1805 /* Increment loop variables */
1813 //! Wrapper function for calling the templated thread-local spread function
1814 static void spreadForceWrapper(ArrayRef<const RVec> x,
1816 const VirialHandling virialHandling,
1817 ArrayRef<RVec> fshift,
1819 const bool clearDxdf,
1820 ArrayRef<const t_iparams> ip,
1821 ArrayRef<const InteractionList> ilist,
1822 const t_pbc* pbc_null)
1824 if (virialHandling == VirialHandling::NonLinear && clearDxdf)
1829 switch (virialHandling)
1831 case VirialHandling::None:
1832 spreadForceForThread<VirialHandling::None>(x, f, fshift, dxdf, ip, ilist, pbc_null);
1834 case VirialHandling::Pbc:
1835 spreadForceForThread<VirialHandling::Pbc>(x, f, fshift, dxdf, ip, ilist, pbc_null);
1837 case VirialHandling::NonLinear:
1838 spreadForceForThread<VirialHandling::NonLinear>(x, f, fshift, dxdf, ip, ilist, pbc_null);
1843 //! Clears the task force buffer elements that are written by task idTask
1844 static void clearTaskForceBufferUsedElements(InterdependentTask* idTask)
1846 int ntask = idTask->spreadTask.size();
1847 for (int ti = 0; ti < ntask; ti++)
1849 const AtomIndex* atomList = &idTask->atomIndex[idTask->spreadTask[ti]];
1850 int natom = atomList->atom.size();
1851 RVec* force = idTask->force.data();
1852 for (int i = 0; i < natom; i++)
1854 clear_rvec(force[atomList->atom[i]]);
1859 void VirtualSitesHandler::Impl::spreadForces(ArrayRef<const RVec> x,
1861 const VirialHandling virialHandling,
1862 ArrayRef<RVec> fshift,
1866 gmx_wallcycle* wcycle)
1868 wallcycle_start(wcycle, ewcVSITESPREAD);
1870 const bool useDomdec = domainInfo_.useDomdec();
1872 t_pbc pbc, *pbc_null;
1874 if (domainInfo_.useMolPbc_)
1876 /* This is wasting some CPU time as we now do this multiple times
1879 pbc_null = set_pbc_dd(&pbc, domainInfo_.pbcType_,
1880 useDomdec ? domainInfo_.domdec_->numCells : nullptr, FALSE, box);
1889 dd_clear_f_vsites(*domainInfo_.domdec_, f);
1892 const int numThreads = threadingInfo_.numThreads();
1894 if (numThreads == 1)
1897 spreadForceWrapper(x, f, virialHandling, fshift, dxdf, true, iparams_, ilists_, pbc_null);
1899 if (virialHandling == VirialHandling::NonLinear)
1901 for (int i = 0; i < DIM; i++)
1903 for (int j = 0; j < DIM; j++)
1905 virial[i][j] += -0.5 * dxdf[i][j];
1912 /* First spread the vsites that might depend on non-local vsites */
1913 auto& nlDependentVSites = threadingInfo_.threadDataNonLocalDependent();
1914 spreadForceWrapper(x, f, virialHandling, fshift, nlDependentVSites.dxdf, true, iparams_,
1915 nlDependentVSites.ilist, pbc_null);
1917 #pragma omp parallel num_threads(numThreads)
1921 int thread = gmx_omp_get_thread_num();
1922 VsiteThread& tData = threadingInfo_.threadData(thread);
1924 ArrayRef<RVec> fshift_t;
1925 if (virialHandling == VirialHandling::Pbc)
1933 fshift_t = tData.fshift;
1935 for (int i = 0; i < SHIFTS; i++)
1937 clear_rvec(fshift_t[i]);
1942 if (tData.useInterdependentTask)
1944 /* Spread the vsites that spread outside our local range.
1945 * This is done using a thread-local force buffer force.
1946 * First we need to copy the input vsite forces to force.
1948 InterdependentTask* idTask = &tData.idTask;
1950 /* Clear the buffer elements set by our task during
1951 * the last call to spread_vsite_f.
1953 clearTaskForceBufferUsedElements(idTask);
1955 int nvsite = idTask->vsite.size();
1956 for (int i = 0; i < nvsite; i++)
1958 copy_rvec(f[idTask->vsite[i]], idTask->force[idTask->vsite[i]]);
1960 spreadForceWrapper(x, idTask->force, virialHandling, fshift_t, tData.dxdf, true,
1961 iparams_, tData.idTask.ilist, pbc_null);
1963 /* We need a barrier before reducing forces below
1964 * that have been produced by a different thread above.
1968 /* Loop over all thread task and reduce forces they
1969 * produced on atoms that fall in our range.
1970 * Note that atomic reduction would be a simpler solution,
1971 * but that might not have good support on all platforms.
1973 int ntask = idTask->reduceTask.size();
1974 for (int ti = 0; ti < ntask; ti++)
1976 const InterdependentTask& idt_foreign =
1977 threadingInfo_.threadData(idTask->reduceTask[ti]).idTask;
1978 const AtomIndex& atomList = idt_foreign.atomIndex[thread];
1979 const RVec* f_foreign = idt_foreign.force.data();
1981 for (int ind : atomList.atom)
1983 rvec_inc(f[ind], f_foreign[ind]);
1984 /* Clearing of f_foreign is done at the next step */
1987 /* Clear the vsite forces, both in f and force */
1988 for (int i = 0; i < nvsite; i++)
1990 int ind = tData.idTask.vsite[i];
1992 clear_rvec(tData.idTask.force[ind]);
1996 /* Spread the vsites that spread locally only */
1997 spreadForceWrapper(x, f, virialHandling, fshift_t, tData.dxdf, false, iparams_,
1998 tData.ilist, pbc_null);
2000 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2003 if (virialHandling == VirialHandling::Pbc)
2005 for (int th = 1; th < numThreads; th++)
2007 for (int i = 0; i < SHIFTS; i++)
2009 rvec_inc(fshift[i], threadingInfo_.threadData(th).fshift[i]);
2014 if (virialHandling == VirialHandling::NonLinear)
2016 for (int th = 0; th < numThreads + 1; th++)
2018 /* MSVC doesn't like matrix references, so we use a pointer */
2019 const matrix& dxdf = threadingInfo_.threadData(th).dxdf;
2021 for (int i = 0; i < DIM; i++)
2023 for (int j = 0; j < DIM; j++)
2025 virial[i][j] += -0.5 * dxdf[i][j];
2034 dd_move_f_vsites(*domainInfo_.domdec_, f, fshift);
2037 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(ilists_, F_VSITE2));
2038 inc_nrnb(nrnb, eNR_VSITE2FD, vsite_count(ilists_, F_VSITE2FD));
2039 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(ilists_, F_VSITE3));
2040 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(ilists_, F_VSITE3FD));
2041 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(ilists_, F_VSITE3FAD));
2042 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(ilists_, F_VSITE3OUT));
2043 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(ilists_, F_VSITE4FD));
2044 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(ilists_, F_VSITE4FDN));
2045 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(ilists_, F_VSITEN));
2047 wallcycle_stop(wcycle, ewcVSITESPREAD);
2050 /*! \brief Returns the an array with group indices for each atom
2052 * \param[in] grouping The paritioning of the atom range into atom groups
2054 static std::vector<int> makeAtomToGroupMapping(const gmx::RangePartitioning& grouping)
2056 std::vector<int> atomToGroup(grouping.fullRange().end(), 0);
2058 for (int group = 0; group < grouping.numBlocks(); group++)
2060 auto block = grouping.block(group);
2061 std::fill(atomToGroup.begin() + block.begin(), atomToGroup.begin() + block.end(), group);
2067 int countNonlinearVsites(const gmx_mtop_t& mtop)
2069 int numNonlinearVsites = 0;
2070 for (const gmx_molblock_t& molb : mtop.molblock)
2072 const gmx_moltype_t& molt = mtop.moltype[molb.type];
2074 for (const auto& ilist : extractILists(molt.ilist, IF_VSITE))
2076 if (ilist.functionType != F_VSITE2 && ilist.functionType != F_VSITE3
2077 && ilist.functionType != F_VSITEN)
2079 numNonlinearVsites += molb.nmol * ilist.iatoms.size() / (1 + NRAL(ilist.functionType));
2084 return numNonlinearVsites;
2087 void VirtualSitesHandler::spreadForces(ArrayRef<const RVec> x,
2089 const VirialHandling virialHandling,
2090 ArrayRef<RVec> fshift,
2094 gmx_wallcycle* wcycle)
2096 impl_->spreadForces(x, f, virialHandling, fshift, virial, nrnb, box, wcycle);
2099 int countInterUpdategroupVsites(const gmx_mtop_t& mtop,
2100 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingPerMoleculetype)
2102 int n_intercg_vsite = 0;
2103 for (const gmx_molblock_t& molb : mtop.molblock)
2105 const gmx_moltype_t& molt = mtop.moltype[molb.type];
2107 std::vector<int> atomToGroup;
2108 if (!updateGroupingPerMoleculetype.empty())
2110 atomToGroup = makeAtomToGroupMapping(updateGroupingPerMoleculetype[molb.type]);
2112 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2114 const int nral = NRAL(ftype);
2115 const InteractionList& il = molt.ilist[ftype];
2116 for (int i = 0; i < il.size(); i += 1 + nral)
2118 bool isInterGroup = atomToGroup.empty();
2121 const int group = atomToGroup[il.iatoms[1 + i]];
2122 for (int a = 1; a < nral; a++)
2124 if (atomToGroup[il.iatoms[1 + a]] != group)
2126 isInterGroup = true;
2133 n_intercg_vsite += molb.nmol;
2139 return n_intercg_vsite;
2142 std::unique_ptr<VirtualSitesHandler> makeVirtualSitesHandler(const gmx_mtop_t& mtop,
2143 const t_commrec* cr,
2146 GMX_RELEASE_ASSERT(cr != nullptr, "We need a valid commrec");
2148 std::unique_ptr<VirtualSitesHandler> vsite;
2150 /* check if there are vsites */
2152 for (int ftype = 0; ftype < F_NRE; ftype++)
2154 if (interaction_function[ftype].flags & IF_VSITE)
2156 GMX_ASSERT(ftype >= c_ftypeVsiteStart && ftype < c_ftypeVsiteEnd,
2157 "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
2159 nvsite += gmx_mtop_ftype_count(&mtop, ftype);
2163 GMX_ASSERT(ftype < c_ftypeVsiteStart || ftype >= c_ftypeVsiteEnd,
2164 "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
2173 return std::make_unique<VirtualSitesHandler>(mtop, cr->dd, pbcType);
2176 ThreadingInfo::ThreadingInfo() : numThreads_(gmx_omp_nthreads_get(emntVSITE))
2178 if (numThreads_ > 1)
2180 /* We need one extra thread data structure for the overlap vsites */
2181 tData_.resize(numThreads_ + 1);
2182 #pragma omp parallel for num_threads(numThreads_) schedule(static)
2183 for (int thread = 0; thread < numThreads_; thread++)
2187 tData_[thread] = std::make_unique<VsiteThread>();
2189 InterdependentTask& idTask = tData_[thread]->idTask;
2191 idTask.atomIndex.resize(numThreads_);
2193 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2195 if (numThreads_ > 1)
2197 tData_[numThreads_] = std::make_unique<VsiteThread>();
2202 //! Returns the number of inter update-group vsites
2203 static int getNumInterUpdategroupVsites(const gmx_mtop_t& mtop, const gmx_domdec_t* domdec)
2205 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingPerMoleculetype;
2208 updateGroupingPerMoleculetype = getUpdateGroupingPerMoleculetype(*domdec);
2211 return countInterUpdategroupVsites(mtop, updateGroupingPerMoleculetype);
2214 VirtualSitesHandler::Impl::Impl(const gmx_mtop_t& mtop, gmx_domdec_t* domdec, const PbcType pbcType) :
2215 numInterUpdategroupVirtualSites_(getNumInterUpdategroupVsites(mtop, domdec)),
2216 domainInfo_({ pbcType, pbcType != PbcType::No && numInterUpdategroupVirtualSites_ > 0, domdec }),
2217 iparams_(mtop.ffparams.iparams)
2221 VirtualSitesHandler::VirtualSitesHandler(const gmx_mtop_t& mtop, gmx_domdec_t* domdec, const PbcType pbcType) :
2222 impl_(new Impl(mtop, domdec, pbcType))
2226 //! Flag that atom \p atom which is home in another task, if it has not already been added before
2227 static inline void flagAtom(InterdependentTask* idTask, const int atom, const int numThreads, const int numAtomsPerThread)
2229 if (!idTask->use[atom])
2231 idTask->use[atom] = true;
2232 int thread = atom / numAtomsPerThread;
2233 /* Assign all non-local atom force writes to thread 0 */
2234 if (thread >= numThreads)
2238 idTask->atomIndex[thread].atom.push_back(atom);
2242 /*! \brief Here we try to assign all vsites that are in our local range.
2244 * Our task local atom range is tData->rangeStart - tData->rangeEnd.
2245 * Vsites that depend only on local atoms, as indicated by taskIndex[]==thread,
2246 * are assigned to task tData->ilist. Vsites that depend on non-local atoms
2247 * but not on other vsites are assigned to task tData->id_task.ilist.
2248 * taskIndex[] is set for all vsites in our range, either to our local tasks
2249 * or to the single last task as taskIndex[]=2*nthreads.
2251 static void assignVsitesToThread(VsiteThread* tData,
2255 gmx::ArrayRef<int> taskIndex,
2256 ArrayRef<const InteractionList> ilist,
2257 ArrayRef<const t_iparams> ip,
2258 const unsigned short* ptype)
2260 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2262 tData->ilist[ftype].clear();
2263 tData->idTask.ilist[ftype].clear();
2265 int nral1 = 1 + NRAL(ftype);
2267 const int* iat = ilist[ftype].iatoms.data();
2268 for (int i = 0; i < ilist[ftype].size();)
2270 if (ftype == F_VSITEN)
2272 /* The 3 below is from 1+NRAL(ftype)=3 */
2273 inc = ip[iat[i]].vsiten.n * 3;
2276 if (iat[1 + i] < tData->rangeStart || iat[1 + i] >= tData->rangeEnd)
2278 /* This vsite belongs to a different thread */
2283 /* We would like to assign this vsite to task thread,
2284 * but it might depend on atoms outside the atom range of thread
2285 * or on another vsite not assigned to task thread.
2288 if (ftype != F_VSITEN)
2290 for (int j = i + 2; j < i + nral1; j++)
2292 /* Do a range check to avoid a harmless race on taskIndex */
2293 if (iat[j] < tData->rangeStart || iat[j] >= tData->rangeEnd || taskIndex[iat[j]] != thread)
2295 if (!tData->useInterdependentTask || ptype[iat[j]] == eptVSite)
2297 /* At least one constructing atom is a vsite
2298 * that is not assigned to the same thread.
2299 * Put this vsite into a separate task.
2305 /* There are constructing atoms outside our range,
2306 * put this vsite into a second task to be executed
2307 * on the same thread. During construction no barrier
2308 * is needed between the two tasks on the same thread.
2309 * During spreading we need to run this task with
2310 * an additional thread-local intermediate force buffer
2311 * (or atomic reduction) and a barrier between the two
2314 task = nthread + thread;
2320 for (int j = i + 2; j < i + inc; j += 3)
2322 /* Do a range check to avoid a harmless race on taskIndex */
2323 if (iat[j] < tData->rangeStart || iat[j] >= tData->rangeEnd || taskIndex[iat[j]] != thread)
2325 GMX_ASSERT(ptype[iat[j]] != eptVSite,
2326 "A vsite to be assigned in assignVsitesToThread has a vsite as "
2327 "a constructing atom that does not belong to our task, such "
2328 "vsites should be assigned to the single 'master' task");
2330 task = nthread + thread;
2335 /* Update this vsite's thread index entry */
2336 taskIndex[iat[1 + i]] = task;
2338 if (task == thread || task == nthread + thread)
2340 /* Copy this vsite to the thread data struct of thread */
2341 InteractionList* il_task;
2344 il_task = &tData->ilist[ftype];
2348 il_task = &tData->idTask.ilist[ftype];
2350 /* Copy the vsite data to the thread-task local array */
2351 il_task->push_back(iat[i], nral1 - 1, iat + i + 1);
2352 if (task == nthread + thread)
2354 /* This vsite write outside our own task force block.
2355 * Put it into the interdependent task list and flag
2356 * the atoms involved for reduction.
2358 tData->idTask.vsite.push_back(iat[i + 1]);
2359 if (ftype != F_VSITEN)
2361 for (int j = i + 2; j < i + nral1; j++)
2363 flagAtom(&tData->idTask, iat[j], nthread, natperthread);
2368 for (int j = i + 2; j < i + inc; j += 3)
2370 flagAtom(&tData->idTask, iat[j], nthread, natperthread);
2381 /*! \brief Assign all vsites with taskIndex[]==task to task tData */
2382 static void assignVsitesToSingleTask(VsiteThread* tData,
2384 gmx::ArrayRef<const int> taskIndex,
2385 ArrayRef<const InteractionList> ilist,
2386 ArrayRef<const t_iparams> ip)
2388 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2390 tData->ilist[ftype].clear();
2391 tData->idTask.ilist[ftype].clear();
2393 int nral1 = 1 + NRAL(ftype);
2395 const int* iat = ilist[ftype].iatoms.data();
2396 InteractionList* il_task = &tData->ilist[ftype];
2398 for (int i = 0; i < ilist[ftype].size();)
2400 if (ftype == F_VSITEN)
2402 /* The 3 below is from 1+NRAL(ftype)=3 */
2403 inc = ip[iat[i]].vsiten.n * 3;
2405 /* Check if the vsite is assigned to our task */
2406 if (taskIndex[iat[1 + i]] == task)
2408 /* Copy the vsite data to the thread-task local array */
2409 il_task->push_back(iat[i], inc - 1, iat + i + 1);
2417 void ThreadingInfo::setVirtualSites(ArrayRef<const InteractionList> ilists,
2418 ArrayRef<const t_iparams> iparams,
2419 const t_mdatoms& mdatoms,
2420 const bool useDomdec)
2422 if (numThreads_ <= 1)
2428 /* The current way of distributing the vsites over threads in primitive.
2429 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
2430 * without taking into account how the vsites are distributed.
2431 * Without domain decomposition we at least tighten the upper bound
2432 * of the range (useful for common systems such as a vsite-protein
2434 * With domain decomposition, as long as the vsites are distributed
2435 * uniformly in each domain along the major dimension, usually x,
2436 * it will also perform well.
2438 int vsite_atom_range;
2442 vsite_atom_range = -1;
2443 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2446 if (ftype != F_VSITEN)
2448 int nral1 = 1 + NRAL(ftype);
2449 ArrayRef<const int> iat = ilists[ftype].iatoms;
2450 for (int i = 0; i < ilists[ftype].size(); i += nral1)
2452 for (int j = i + 1; j < i + nral1; j++)
2454 vsite_atom_range = std::max(vsite_atom_range, iat[j]);
2462 ArrayRef<const int> iat = ilists[ftype].iatoms;
2465 while (i < ilists[ftype].size())
2467 /* The 3 below is from 1+NRAL(ftype)=3 */
2468 vs_ind_end = i + iparams[iat[i]].vsiten.n * 3;
2470 vsite_atom_range = std::max(vsite_atom_range, iat[i + 1]);
2471 while (i < vs_ind_end)
2473 vsite_atom_range = std::max(vsite_atom_range, iat[i + 2]);
2481 natperthread = (vsite_atom_range + numThreads_ - 1) / numThreads_;
2485 /* Any local or not local atom could be involved in virtual sites.
2486 * But since we usually have very few non-local virtual sites
2487 * (only non-local vsites that depend on local vsites),
2488 * we distribute the local atom range equally over the threads.
2489 * When assigning vsites to threads, we should take care that the last
2490 * threads also covers the non-local range.
2492 vsite_atom_range = mdatoms.nr;
2493 natperthread = (mdatoms.homenr + numThreads_ - 1) / numThreads_;
2498 fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n",
2499 mdatoms.nr, vsite_atom_range, natperthread);
2502 /* To simplify the vsite assignment, we make an index which tells us
2503 * to which task particles, both non-vsites and vsites, are assigned.
2505 taskIndex_.resize(mdatoms.nr);
2507 /* Initialize the task index array. Here we assign the non-vsite
2508 * particles to task=thread, so we easily figure out if vsites
2509 * depend on local and/or non-local particles in assignVsitesToThread.
2513 for (int i = 0; i < mdatoms.nr; i++)
2515 if (mdatoms.ptype[i] == eptVSite)
2517 /* vsites are not assigned to a task yet */
2522 /* assign non-vsite particles to task thread */
2523 taskIndex_[i] = thread;
2525 if (i == (thread + 1) * natperthread && thread < numThreads_)
2532 #pragma omp parallel num_threads(numThreads_)
2536 int thread = gmx_omp_get_thread_num();
2537 VsiteThread& tData = *tData_[thread];
2539 /* Clear the buffer use flags that were set before */
2540 if (tData.useInterdependentTask)
2542 InterdependentTask& idTask = tData.idTask;
2544 /* To avoid an extra OpenMP barrier in spread_vsite_f,
2545 * we clear the force buffer at the next step,
2546 * so we need to do it here as well.
2548 clearTaskForceBufferUsedElements(&idTask);
2550 idTask.vsite.resize(0);
2551 for (int t = 0; t < numThreads_; t++)
2553 AtomIndex& atomIndex = idTask.atomIndex[t];
2554 int natom = atomIndex.atom.size();
2555 for (int i = 0; i < natom; i++)
2557 idTask.use[atomIndex.atom[i]] = false;
2559 atomIndex.atom.resize(0);
2564 /* To avoid large f_buf allocations of #threads*vsite_atom_range
2565 * we don't use task2 with more than 200000 atoms. This doesn't
2566 * affect performance, since with such a large range relatively few
2567 * vsites will end up in the separate task.
2568 * Note that useTask2 should be the same for all threads.
2570 tData.useInterdependentTask = (vsite_atom_range <= 200000);
2571 if (tData.useInterdependentTask)
2573 size_t natoms_use_in_vsites = vsite_atom_range;
2574 InterdependentTask& idTask = tData.idTask;
2575 /* To avoid resizing and re-clearing every nstlist steps,
2576 * we never down size the force buffer.
2578 if (natoms_use_in_vsites > idTask.force.size() || natoms_use_in_vsites > idTask.use.size())
2580 idTask.force.resize(natoms_use_in_vsites, { 0, 0, 0 });
2581 idTask.use.resize(natoms_use_in_vsites, false);
2585 /* Assign all vsites that can execute independently on threads */
2586 tData.rangeStart = thread * natperthread;
2587 if (thread < numThreads_ - 1)
2589 tData.rangeEnd = (thread + 1) * natperthread;
2593 /* The last thread should cover up to the end of the range */
2594 tData.rangeEnd = mdatoms.nr;
2596 assignVsitesToThread(&tData, thread, numThreads_, natperthread, taskIndex_, ilists,
2597 iparams, mdatoms.ptype);
2599 if (tData.useInterdependentTask)
2601 /* In the worst case, all tasks write to force ranges of
2602 * all other tasks, leading to #tasks^2 scaling (this is only
2603 * the overhead, the actual flops remain constant).
2604 * But in most cases there is far less coupling. To improve
2605 * scaling at high thread counts we therefore construct
2606 * an index to only loop over the actually affected tasks.
2608 InterdependentTask& idTask = tData.idTask;
2610 /* Ensure assignVsitesToThread finished on other threads */
2613 idTask.spreadTask.resize(0);
2614 idTask.reduceTask.resize(0);
2615 for (int t = 0; t < numThreads_; t++)
2617 /* Do we write to the force buffer of task t? */
2618 if (!idTask.atomIndex[t].atom.empty())
2620 idTask.spreadTask.push_back(t);
2622 /* Does task t write to our force buffer? */
2623 if (!tData_[t]->idTask.atomIndex[thread].atom.empty())
2625 idTask.reduceTask.push_back(t);
2630 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2632 /* Assign all remaining vsites, that will have taskIndex[]=2*vsite->nthreads,
2633 * to a single task that will not run in parallel with other tasks.
2635 assignVsitesToSingleTask(tData_[numThreads_].get(), 2 * numThreads_, taskIndex_, ilists, iparams);
2637 if (debug && numThreads_ > 1)
2639 fprintf(debug, "virtual site useInterdependentTask %d, nuse:\n",
2640 static_cast<int>(tData_[0]->useInterdependentTask));
2641 for (int th = 0; th < numThreads_ + 1; th++)
2643 fprintf(debug, " %4d", tData_[th]->idTask.nuse);
2645 fprintf(debug, "\n");
2647 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2649 if (!ilists[ftype].empty())
2651 fprintf(debug, "%-20s thread dist:", interaction_function[ftype].longname);
2652 for (int th = 0; th < numThreads_ + 1; th++)
2654 fprintf(debug, " %4d %4d ", tData_[th]->ilist[ftype].size(),
2655 tData_[th]->idTask.ilist[ftype].size());
2657 fprintf(debug, "\n");
2663 int nrOrig = vsiteIlistNrCount(ilists);
2665 for (int th = 0; th < numThreads_ + 1; th++)
2667 nrThreaded += vsiteIlistNrCount(tData_[th]->ilist) + vsiteIlistNrCount(tData_[th]->idTask.ilist);
2669 GMX_ASSERT(nrThreaded == nrOrig,
2670 "The number of virtual sites assigned to all thread task has to match the total "
2671 "number of virtual sites");
2675 void VirtualSitesHandler::Impl::setVirtualSites(ArrayRef<const InteractionList> ilists,
2676 const t_mdatoms& mdatoms)
2680 threadingInfo_.setVirtualSites(ilists, iparams_, mdatoms, domainInfo_.useDomdec());
2683 void VirtualSitesHandler::setVirtualSites(ArrayRef<const InteractionList> ilists, const t_mdatoms& mdatoms)
2685 impl_->setVirtualSites(ilists, mdatoms);