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,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Implements the VirtualSitesHandler class and vsite standalone functions
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_mdlib
55 #include "gromacs/domdec/domdec.h"
56 #include "gromacs/domdec/domdec_struct.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/gmxlib/nrnb.h"
59 #include "gromacs/math/functions.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/mdlib/gmx_omp_nthreads.h"
62 #include "gromacs/mdtypes/commrec.h"
63 #include "gromacs/mdtypes/mdatom.h"
64 #include "gromacs/pbcutil/ishift.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/timing/wallcycle.h"
67 #include "gromacs/topology/ifunc.h"
68 #include "gromacs/topology/mtop_util.h"
69 #include "gromacs/topology/topology.h"
70 #include "gromacs/utility/exceptions.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/gmxassert.h"
73 #include "gromacs/utility/gmxomp.h"
75 /* The strategy used here for assigning virtual sites to (thread-)tasks
78 * We divide the atom range that vsites operate on (natoms_local with DD,
79 * 0 - last atom involved in vsites without DD) equally over all threads.
81 * Vsites in the local range constructed from atoms in the local range
82 * and/or other vsites that are fully local are assigned to a simple,
85 * Vsites that are not assigned after using the above criterion get assigned
86 * to a so called "interdependent" thread task when none of the constructing
87 * atoms is a vsite. These tasks are called interdependent, because one task
88 * accesses atoms assigned to a different task/thread.
89 * Note that this option is turned off with large (local) atom counts
90 * to avoid high memory usage.
92 * Any remaining vsites are assigned to a separate master thread task.
97 //! VirialHandling is often used outside VirtualSitesHandler class members
98 using VirialHandling = VirtualSitesHandler::VirialHandling;
100 /*! \brief Information on PBC and domain decomposition for virtual sites
105 //! Constructs without PBC and DD
106 DomainInfo() = default;
108 //! Constructs with PBC and DD, if !=nullptr
109 DomainInfo(PbcType pbcType, bool haveInterUpdateGroupVirtualSites, gmx_domdec_t* domdec) :
111 useMolPbc_(pbcType != PbcType::No && haveInterUpdateGroupVirtualSites),
116 //! Returns whether we are using domain decomposition with more than 1 DD rank
117 bool useDomdec() const { return (domdec_ != nullptr); }
120 const PbcType pbcType_ = PbcType::No;
121 //! Whether molecules are broken over PBC
122 const bool useMolPbc_ = false;
123 //! Pointer to the domain decomposition struct, nullptr without PP DD
124 const gmx_domdec_t* domdec_ = nullptr;
127 /*! \brief List of atom indices belonging to a task
131 //! List of atom indices
132 std::vector<int> atom;
135 /*! \brief Data structure for thread tasks that use constructing atoms outside their own atom range
137 struct InterdependentTask
139 //! The interaction lists, only vsite entries are used
140 InteractionLists ilist;
141 //! Thread/task-local force buffer
142 std::vector<RVec> force;
143 //! The atom indices of the vsites of our task
144 std::vector<int> vsite;
145 //! Flags if elements in force are spread to or not
146 std::vector<bool> use;
147 //! The number of entries set to true in use
149 //! Array of atoms indices, size nthreads, covering all nuse set elements in use
150 std::vector<AtomIndex> atomIndex;
151 //! List of tasks (force blocks) this task spread forces to
152 std::vector<int> spreadTask;
153 //! List of tasks that write to this tasks force block range
154 std::vector<int> reduceTask;
157 /*! \brief Vsite thread task data structure
161 //! Start of atom range of this task
163 //! End of atom range of this task
165 //! The interaction lists, only vsite entries are used
166 std::array<InteractionList, F_NRE> ilist;
167 //! Local fshift accumulation buffer
168 std::array<RVec, SHIFTS> fshift;
169 //! Local virial dx*df accumulation buffer
171 //! 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
172 bool useInterdependentTask;
173 //! Data for vsites that involve constructing atoms in the atom range of other threads/tasks
174 InterdependentTask idTask;
176 /*! \brief Constructor */
181 for (auto& elem : fshift)
183 elem = { 0.0_real, 0.0_real, 0.0_real };
186 useInterdependentTask = false;
191 /*! \brief Information on how the virtual site work is divided over thread tasks
196 //! Constructor, retrieves the number of threads to use from gmx_omp_nthreads.h
199 //! Returns the number of threads to use for vsite operations
200 int numThreads() const { return numThreads_; }
202 //! Returns the thread data for the given thread
203 const VsiteThread& threadData(int threadIndex) const { return *tData_[threadIndex]; }
205 //! Returns the thread data for the given thread
206 VsiteThread& threadData(int threadIndex) { return *tData_[threadIndex]; }
208 //! Returns the thread data for vsites that depend on non-local vsites
209 const VsiteThread& threadDataNonLocalDependent() const { return *tData_[numThreads_]; }
211 //! Returns the thread data for vsites that depend on non-local vsites
212 VsiteThread& threadDataNonLocalDependent() { return *tData_[numThreads_]; }
214 //! Set VSites and distribute VSite work over threads, should be called after DD partitioning
215 void setVirtualSites(ArrayRef<const InteractionList> ilist,
216 ArrayRef<const t_iparams> iparams,
217 const t_mdatoms& mdatoms,
221 //! Number of threads used for vsite operations
222 const int numThreads_;
223 //! Thread local vsites and work structs
224 std::vector<std::unique_ptr<VsiteThread>> tData_;
225 //! Work array for dividing vsites over threads
226 std::vector<int> taskIndex_;
229 /*! \brief Impl class for VirtualSitesHandler
231 class VirtualSitesHandler::Impl
234 //! Constructor, domdec should be nullptr without DD
235 Impl(const gmx_mtop_t& mtop, gmx_domdec_t* domdec, PbcType pbcType);
237 //! Returns the number of virtual sites acting over multiple update groups
238 int numInterUpdategroupVirtualSites() const { return numInterUpdategroupVirtualSites_; }
240 //! Set VSites and distribute VSite work over threads, should be called after DD partitioning
241 void setVirtualSites(ArrayRef<const InteractionList> ilist, const t_mdatoms& mdatoms);
243 /*! \brief Create positions of vsite atoms based for the local system
245 * \param[in,out] x The coordinates
246 * \param[in,out] v The velocities, needed if operation requires it
247 * \param[in] box The box
248 * \param[in] operation Whether we calculate positions, velocities, or both
250 void construct(ArrayRef<RVec> x, ArrayRef<RVec> v, const matrix box, VSiteOperation operation) const;
252 /*! \brief Spread the force operating on the vsite atoms on the surrounding atoms.
254 * vsite should point to a valid object.
255 * The virialHandling parameter determines how virial contributions are handled.
256 * If this is set to Linear, shift forces are accumulated into fshift.
257 * If this is set to NonLinear, non-linear contributions are added to virial.
258 * This non-linear correction is required when the virial is not calculated
259 * afterwards from the particle position and forces, but in a different way,
260 * as for instance for the PME mesh contribution.
262 void spreadForces(ArrayRef<const RVec> x,
264 VirialHandling virialHandling,
265 ArrayRef<RVec> fshift,
269 gmx_wallcycle* wcycle);
272 //! The number of vsites that cross update groups, when =0 no PBC treatment is needed
273 const int numInterUpdategroupVirtualSites_;
274 //! PBC and DD information
275 const DomainInfo domainInfo_;
276 //! The interaction parameters
277 const ArrayRef<const t_iparams> iparams_;
278 //! The interaction lists
279 ArrayRef<const InteractionList> ilists_;
280 //! Information for handling vsite threading
281 ThreadingInfo threadingInfo_;
284 VirtualSitesHandler::~VirtualSitesHandler() = default;
286 int VirtualSitesHandler::numInterUpdategroupVirtualSites() const
288 return impl_->numInterUpdategroupVirtualSites();
291 /*! \brief Returns the sum of the vsite ilist sizes over all vsite types
293 * \param[in] ilist The interaction list
295 static int vsiteIlistNrCount(ArrayRef<const InteractionList> ilist)
298 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
300 nr += ilist[ftype].size();
306 //! Computes the distance between xi and xj, pbc is used when pbc!=nullptr
307 static int pbc_rvec_sub(const t_pbc* pbc, const rvec xi, const rvec xj, rvec dx)
311 return pbc_dx_aiuc(pbc, xi, xj, dx);
315 rvec_sub(xi, xj, dx);
320 //! Returns the 1/norm(x)
321 static inline real inverseNorm(const rvec x)
323 return gmx::invsqrt(iprod(x, x));
326 //! Whether we're calculating the virtual site position
327 enum class VSiteCalculatePosition
332 //! Whether we're calculating the virtual site velocity
333 enum class VSiteCalculateVelocity
340 /* Vsite construction routines */
342 // GCC 8 falsely flags unused variables if constexpr prunes a code path, fixed in GCC 9
343 // https://gcc.gnu.org/bugzilla/show_bug.cgi?id=85827
345 GCC_DIAGNOSTIC_IGNORE(-Wunused-but-set-parameter)
348 template<VSiteCalculatePosition calculatePosition, VSiteCalculateVelocity calculateVelocity>
349 static void constr_vsite1(const rvec xi, rvec x, const rvec vi, rvec v)
351 if (calculatePosition == VSiteCalculatePosition::Yes)
356 if (calculateVelocity == VSiteCalculateVelocity::Yes)
362 template<VSiteCalculatePosition calculatePosition, VSiteCalculateVelocity calculateVelocity>
364 constr_vsite2(const rvec xi, const rvec xj, rvec x, real a, const t_pbc* pbc, const rvec vi, const rvec vj, rvec v)
366 const real b = 1 - a;
369 if (calculatePosition == VSiteCalculatePosition::Yes)
374 pbc_dx_aiuc(pbc, xj, xi, dx);
375 x[XX] = xi[XX] + a * dx[XX];
376 x[YY] = xi[YY] + a * dx[YY];
377 x[ZZ] = xi[ZZ] + a * dx[ZZ];
381 x[XX] = b * xi[XX] + a * xj[XX];
382 x[YY] = b * xi[YY] + a * xj[YY];
383 x[ZZ] = b * xi[ZZ] + a * xj[ZZ];
386 /* TOTAL: 10 flops */
388 if (calculateVelocity == VSiteCalculateVelocity::Yes)
390 v[XX] = b * vi[XX] + a * vj[XX];
391 v[YY] = b * vi[YY] + a * vj[YY];
392 v[ZZ] = b * vi[ZZ] + a * vj[ZZ];
396 template<VSiteCalculatePosition calculatePosition, VSiteCalculateVelocity calculateVelocity>
398 constr_vsite2FD(const rvec xi, const rvec xj, rvec x, real a, const t_pbc* pbc, const rvec vi, const rvec vj, rvec v)
401 pbc_rvec_sub(pbc, xj, xi, xij);
404 const real invNormXij = inverseNorm(xij);
405 const real b = a * invNormXij;
408 if (calculatePosition == VSiteCalculatePosition::Yes)
410 x[XX] = xi[XX] + b * xij[XX];
411 x[YY] = xi[YY] + b * xij[YY];
412 x[ZZ] = xi[ZZ] + b * xij[ZZ];
414 /* TOTAL: 25 flops */
416 if (calculateVelocity == VSiteCalculateVelocity::Yes)
419 rvec_sub(vj, vi, vij);
420 const real vijDotXij = iprod(vij, xij);
422 v[XX] = vi[XX] + b * (vij[XX] - xij[XX] * vijDotXij * invNormXij * invNormXij);
423 v[YY] = vi[YY] + b * (vij[YY] - xij[YY] * vijDotXij * invNormXij * invNormXij);
424 v[ZZ] = vi[ZZ] + b * (vij[ZZ] - xij[ZZ] * vijDotXij * invNormXij * invNormXij);
428 template<VSiteCalculatePosition calculatePosition, VSiteCalculateVelocity calculateVelocity>
429 static void constr_vsite3(const rvec xi,
441 const real c = 1 - a - b;
444 if (calculatePosition == VSiteCalculatePosition::Yes)
450 pbc_dx_aiuc(pbc, xj, xi, dxj);
451 pbc_dx_aiuc(pbc, xk, xi, dxk);
452 x[XX] = xi[XX] + a * dxj[XX] + b * dxk[XX];
453 x[YY] = xi[YY] + a * dxj[YY] + b * dxk[YY];
454 x[ZZ] = xi[ZZ] + a * dxj[ZZ] + b * dxk[ZZ];
458 x[XX] = c * xi[XX] + a * xj[XX] + b * xk[XX];
459 x[YY] = c * xi[YY] + a * xj[YY] + b * xk[YY];
460 x[ZZ] = c * xi[ZZ] + a * xj[ZZ] + b * xk[ZZ];
463 /* TOTAL: 17 flops */
465 if (calculateVelocity == VSiteCalculateVelocity::Yes)
467 v[XX] = c * vi[XX] + a * vj[XX] + b * vk[XX];
468 v[YY] = c * vi[YY] + a * vj[YY] + b * vk[YY];
469 v[ZZ] = c * vi[ZZ] + a * vj[ZZ] + b * vk[ZZ];
473 template<VSiteCalculatePosition calculatePosition, VSiteCalculateVelocity calculateVelocity>
474 static void constr_vsite3FD(const rvec xi,
488 pbc_rvec_sub(pbc, xj, xi, xij);
489 pbc_rvec_sub(pbc, xk, xj, xjk);
492 /* temp goes from i to a point on the line jk */
493 temp[XX] = xij[XX] + a * xjk[XX];
494 temp[YY] = xij[YY] + a * xjk[YY];
495 temp[ZZ] = xij[ZZ] + a * xjk[ZZ];
498 const real invNormTemp = inverseNorm(temp);
499 const real c = b * invNormTemp;
502 if (calculatePosition == VSiteCalculatePosition::Yes)
504 x[XX] = xi[XX] + c * temp[XX];
505 x[YY] = xi[YY] + c * temp[YY];
506 x[ZZ] = xi[ZZ] + c * temp[ZZ];
508 /* TOTAL: 34 flops */
510 if (calculateVelocity == VSiteCalculateVelocity::Yes)
514 rvec_sub(vj, vi, vij);
515 rvec_sub(vk, vj, vjk);
516 const rvec tempV = { vij[XX] + a * vjk[XX], vij[YY] + a * vjk[YY], vij[ZZ] + a * vjk[ZZ] };
517 const real tempDotTempV = iprod(temp, tempV);
519 v[XX] = vi[XX] + c * (tempV[XX] - temp[XX] * tempDotTempV * invNormTemp * invNormTemp);
520 v[YY] = vi[YY] + c * (tempV[YY] - temp[YY] * tempDotTempV * invNormTemp * invNormTemp);
521 v[ZZ] = vi[ZZ] + c * (tempV[ZZ] - temp[ZZ] * tempDotTempV * invNormTemp * invNormTemp);
525 template<VSiteCalculatePosition calculatePosition, VSiteCalculateVelocity calculateVelocity>
526 static void constr_vsite3FAD(const rvec xi,
537 { // Note: a = d * cos(theta)
538 // b = d * sin(theta)
541 pbc_rvec_sub(pbc, xj, xi, xij);
542 pbc_rvec_sub(pbc, xk, xj, xjk);
545 const real invdij = inverseNorm(xij);
546 const real xijDotXjk = iprod(xij, xjk);
547 const real c1 = invdij * invdij * xijDotXjk;
548 xp[XX] = xjk[XX] - c1 * xij[XX];
549 xp[YY] = xjk[YY] - c1 * xij[YY];
550 xp[ZZ] = xjk[ZZ] - c1 * xij[ZZ];
551 const real a1 = a * invdij;
552 const real invNormXp = inverseNorm(xp);
553 const real b1 = b * invNormXp;
556 if (calculatePosition == VSiteCalculatePosition::Yes)
558 x[XX] = xi[XX] + a1 * xij[XX] + b1 * xp[XX];
559 x[YY] = xi[YY] + a1 * xij[YY] + b1 * xp[YY];
560 x[ZZ] = xi[ZZ] + a1 * xij[ZZ] + b1 * xp[ZZ];
562 /* TOTAL: 63 flops */
565 if (calculateVelocity == VSiteCalculateVelocity::Yes)
569 rvec_sub(vj, vi, vij);
570 rvec_sub(vk, vj, vjk);
572 const real vijDotXjkPlusXijDotVjk = iprod(vij, xjk) + iprod(xij, vjk);
573 const real xijDotVij = iprod(xij, vij);
574 const real invNormXij2 = invdij * invdij;
578 - xij[XX] * invNormXij2
579 * (vijDotXjkPlusXijDotVjk - invNormXij2 * xijDotXjk * xijDotVij * 2)
580 - vij[XX] * xijDotXjk * invNormXij2;
582 - xij[YY] * invNormXij2
583 * (vijDotXjkPlusXijDotVjk - invNormXij2 * xijDotXjk * xijDotVij * 2)
584 - vij[YY] * xijDotXjk * invNormXij2;
586 - xij[ZZ] * invNormXij2
587 * (vijDotXjkPlusXijDotVjk - invNormXij2 * xijDotXjk * xijDotVij * 2)
588 - vij[ZZ] * xijDotXjk * invNormXij2;
590 const real xpDotVp = iprod(xp, vp);
592 v[XX] = vi[XX] + a1 * (vij[XX] - xij[XX] * xijDotVij * invdij * invdij)
593 + b1 * (vp[XX] - xp[XX] * xpDotVp * invNormXp * invNormXp);
594 v[YY] = vi[YY] + a1 * (vij[YY] - xij[YY] * xijDotVij * invdij * invdij)
595 + b1 * (vp[YY] - xp[YY] * xpDotVp * invNormXp * invNormXp);
596 v[ZZ] = vi[ZZ] + a1 * (vij[ZZ] - xij[ZZ] * xijDotVij * invdij * invdij)
597 + b1 * (vp[ZZ] - xp[ZZ] * xpDotVp * invNormXp * invNormXp);
601 template<VSiteCalculatePosition calculatePosition, VSiteCalculateVelocity calculateVelocity>
602 static void constr_vsite3OUT(const rvec xi,
617 pbc_rvec_sub(pbc, xj, xi, xij);
618 pbc_rvec_sub(pbc, xk, xi, xik);
619 cprod(xij, xik, temp);
622 if (calculatePosition == VSiteCalculatePosition::Yes)
624 x[XX] = xi[XX] + a * xij[XX] + b * xik[XX] + c * temp[XX];
625 x[YY] = xi[YY] + a * xij[YY] + b * xik[YY] + c * temp[YY];
626 x[ZZ] = xi[ZZ] + a * xij[ZZ] + b * xik[ZZ] + c * temp[ZZ];
628 /* TOTAL: 33 flops */
631 if (calculateVelocity == VSiteCalculateVelocity::Yes)
635 rvec_sub(vj, vi, vij);
636 rvec_sub(vk, vi, vik);
640 cprod(vij, xik, temp1);
641 cprod(xij, vik, temp2);
643 v[XX] = vi[XX] + a * vij[XX] + b * vik[XX] + c * (temp1[XX] + temp2[XX]);
644 v[YY] = vi[YY] + a * vij[YY] + b * vik[YY] + c * (temp1[YY] + temp2[YY]);
645 v[ZZ] = vi[ZZ] + a * vij[ZZ] + b * vik[ZZ] + c * (temp1[ZZ] + temp2[ZZ]);
649 template<VSiteCalculatePosition calculatePosition, VSiteCalculateVelocity calculateVelocity>
650 static void constr_vsite4FD(const rvec xi,
665 rvec xij, xjk, xjl, temp;
668 pbc_rvec_sub(pbc, xj, xi, xij);
669 pbc_rvec_sub(pbc, xk, xj, xjk);
670 pbc_rvec_sub(pbc, xl, xj, xjl);
673 /* temp goes from i to a point on the plane jkl */
674 temp[XX] = xij[XX] + a * xjk[XX] + b * xjl[XX];
675 temp[YY] = xij[YY] + a * xjk[YY] + b * xjl[YY];
676 temp[ZZ] = xij[ZZ] + a * xjk[ZZ] + b * xjl[ZZ];
679 const real invRm = inverseNorm(temp);
683 if (calculatePosition == VSiteCalculatePosition::Yes)
685 x[XX] = xi[XX] + d * temp[XX];
686 x[YY] = xi[YY] + d * temp[YY];
687 x[ZZ] = xi[ZZ] + d * temp[ZZ];
689 /* TOTAL: 43 flops */
691 if (calculateVelocity == VSiteCalculateVelocity::Yes)
697 rvec_sub(vj, vi, vij);
698 rvec_sub(vk, vj, vjk);
699 rvec_sub(vl, vj, vjl);
702 vm[XX] = vij[XX] + a * vjk[XX] + b * vjl[XX];
703 vm[YY] = vij[YY] + a * vjk[YY] + b * vjl[YY];
704 vm[ZZ] = vij[ZZ] + a * vjk[ZZ] + b * vjl[ZZ];
706 const real vmDotRm = iprod(vm, temp);
707 v[XX] = vi[XX] + d * (vm[XX] - temp[XX] * vmDotRm * invRm * invRm);
708 v[YY] = vi[YY] + d * (vm[YY] - temp[YY] * vmDotRm * invRm * invRm);
709 v[ZZ] = vi[ZZ] + d * (vm[ZZ] - temp[ZZ] * vmDotRm * invRm * invRm);
713 template<VSiteCalculatePosition calculatePosition, VSiteCalculateVelocity calculateVelocity>
714 static void constr_vsite4FDN(const rvec xi,
729 rvec xij, xik, xil, ra, rb, rja, rjb, rm;
732 pbc_rvec_sub(pbc, xj, xi, xij);
733 pbc_rvec_sub(pbc, xk, xi, xik);
734 pbc_rvec_sub(pbc, xl, xi, xil);
737 ra[XX] = a * xik[XX];
738 ra[YY] = a * xik[YY];
739 ra[ZZ] = a * xik[ZZ];
741 rb[XX] = b * xil[XX];
742 rb[YY] = b * xil[YY];
743 rb[ZZ] = b * xil[ZZ];
747 rvec_sub(ra, xij, rja);
748 rvec_sub(rb, xij, rjb);
754 const real invNormRm = inverseNorm(rm);
758 if (calculatePosition == VSiteCalculatePosition::Yes)
760 x[XX] = xi[XX] + d * rm[XX];
761 x[YY] = xi[YY] + d * rm[YY];
762 x[ZZ] = xi[ZZ] + d * rm[ZZ];
764 /* TOTAL: 47 flops */
767 if (calculateVelocity == VSiteCalculateVelocity::Yes)
772 rvec_sub(vj, vi, vij);
773 rvec_sub(vk, vi, vik);
774 rvec_sub(vl, vi, vil);
779 vja[XX] = a * vik[XX] - vij[XX];
780 vja[YY] = a * vik[YY] - vij[YY];
781 vja[ZZ] = a * vik[ZZ] - vij[ZZ];
782 vjb[XX] = b * vil[XX] - vij[XX];
783 vjb[YY] = b * vil[YY] - vij[YY];
784 vjb[ZZ] = b * vil[ZZ] - vij[ZZ];
788 cprod(vja, rjb, temp1);
789 cprod(rja, vjb, temp2);
792 vm[XX] = temp1[XX] + temp2[XX];
793 vm[YY] = temp1[YY] + temp2[YY];
794 vm[ZZ] = temp1[ZZ] + temp2[ZZ];
796 const real rmDotVm = iprod(rm, vm);
797 v[XX] = vi[XX] + d * (vm[XX] - rm[XX] * rmDotVm * invNormRm * invNormRm);
798 v[YY] = vi[YY] + d * (vm[YY] - rm[YY] * rmDotVm * invNormRm * invNormRm);
799 v[ZZ] = vi[ZZ] + d * (vm[ZZ] - rm[ZZ] * rmDotVm * invNormRm * invNormRm);
803 template<VSiteCalculatePosition calculatePosition, VSiteCalculateVelocity calculateVelocity>
804 static int constr_vsiten(const t_iatom* ia,
805 ArrayRef<const t_iparams> ip,
816 const int n3 = 3 * ip[ia[0]].vsiten.n;
817 const int av = ia[1];
819 copy_rvec(x[ai], x1);
820 copy_rvec(v[ai], v1);
822 for (int i = 3; i < n3; i += 3)
825 a = ip[ia[i]].vsiten.a;
826 if (calculatePosition == VSiteCalculatePosition::Yes)
830 pbc_dx_aiuc(pbc, x[ai], x1, dx);
834 rvec_sub(x[ai], x1, dx);
836 dsum[XX] += a * dx[XX];
837 dsum[YY] += a * dx[YY];
838 dsum[ZZ] += a * dx[ZZ];
841 if (calculateVelocity == VSiteCalculateVelocity::Yes)
843 rvec_sub(v[ai], v1, dx);
844 dvsum[XX] += a * dx[XX];
845 dvsum[YY] += a * dx[YY];
846 dvsum[ZZ] += a * dx[ZZ];
851 if (calculatePosition == VSiteCalculatePosition::Yes)
853 x[av][XX] = x1[XX] + dsum[XX];
854 x[av][YY] = x1[YY] + dsum[YY];
855 x[av][ZZ] = x1[ZZ] + dsum[ZZ];
858 if (calculateVelocity == VSiteCalculateVelocity::Yes)
860 v[av][XX] = v1[XX] + dvsum[XX];
861 v[av][YY] = v1[YY] + dvsum[YY];
862 v[av][ZZ] = v1[ZZ] + dvsum[ZZ];
872 //! PBC modes for vsite construction and spreading
875 all, //!< Apply normal, simple PBC for all vsites
876 none //!< No PBC treatment needed
879 /*! \brief Returns the PBC mode based on the system PBC and vsite properties
881 * \param[in] pbcPtr A pointer to a PBC struct or nullptr when no PBC treatment is required
883 static PbcMode getPbcMode(const t_pbc* pbcPtr)
885 if (pbcPtr == nullptr)
887 return PbcMode::none;
895 /*! \brief Executes the vsite construction task for a single thread
897 * \tparam operation Whether we are calculating positions, velocities, or both
898 * \param[in,out] x Coordinates to construct vsites for
899 * \param[in,out] v Velocities are generated for virtual sites if `operation` requires it
900 * \param[in] ip Interaction parameters for all interaction, only vsite parameters are used
901 * \param[in] ilist The interaction lists, only vsites are usesd
902 * \param[in] pbc_null PBC struct, used for PBC distance calculations when !=nullptr
904 template<VSiteCalculatePosition calculatePosition, VSiteCalculateVelocity calculateVelocity>
905 static void construct_vsites_thread(ArrayRef<RVec> x,
907 ArrayRef<const t_iparams> ip,
908 ArrayRef<const InteractionList> ilist,
909 const t_pbc* pbc_null)
911 if (calculateVelocity == VSiteCalculateVelocity::Yes)
913 GMX_RELEASE_ASSERT(!v.empty(),
914 "Can't calculate velocities without access to velocity vector.");
917 // Work around clang bug (unfixed as of Feb 2021)
918 // https://bugs.llvm.org/show_bug.cgi?id=35450
920 CLANG_DIAGNOSTIC_IGNORE(-Wunused-lambda-capture)
922 // GCC 8 falsely flags unused variables if constexpr prunes a code path, fixed in GCC 9
923 // https://gcc.gnu.org/bugzilla/show_bug.cgi?id=85827
925 GCC_DIAGNOSTIC_IGNORE(-Wunused-but-set-parameter)
927 // getVOrNull returns a velocity rvec if we need it, nullptr otherwise.
928 auto getVOrNull = [v](int idx) -> real* {
929 if (calculateVelocity == VSiteCalculateVelocity::Yes)
931 return v[idx].as_vec();
939 CLANG_DIAGNOSTIC_RESET
941 const PbcMode pbcMode = getPbcMode(pbc_null);
942 /* We need another pbc pointer, as with charge groups we switch per vsite */
943 const t_pbc* pbc_null2 = pbc_null;
945 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
947 if (ilist[ftype].empty())
953 int nra = interaction_function[ftype].nratoms;
955 int nr = ilist[ftype].size();
957 const t_iatom* ia = ilist[ftype].iatoms.data();
959 for (int i = 0; i < nr;)
962 /* The vsite and constructing atoms */
965 /* Constants for constructing vsites */
966 real a1 = ip[tp].vsite.a;
967 /* Copy the old position */
969 copy_rvec(x[avsite], xv);
971 /* Construct the vsite depending on type */
977 constr_vsite1<calculatePosition, calculateVelocity>(
978 x[ai], x[avsite], getVOrNull(ai), getVOrNull(avsite));
982 constr_vsite2<calculatePosition, calculateVelocity>(x[ai],
993 constr_vsite2FD<calculatePosition, calculateVelocity>(x[ai],
1000 getVOrNull(avsite));
1005 b1 = ip[tp].vsite.b;
1006 constr_vsite3<calculatePosition, calculateVelocity>(x[ai],
1016 getVOrNull(avsite));
1021 b1 = ip[tp].vsite.b;
1022 constr_vsite3FD<calculatePosition, calculateVelocity>(x[ai],
1032 getVOrNull(avsite));
1037 b1 = ip[tp].vsite.b;
1038 constr_vsite3FAD<calculatePosition, calculateVelocity>(x[ai],
1048 getVOrNull(avsite));
1053 b1 = ip[tp].vsite.b;
1054 c1 = ip[tp].vsite.c;
1055 constr_vsite3OUT<calculatePosition, calculateVelocity>(x[ai],
1066 getVOrNull(avsite));
1072 b1 = ip[tp].vsite.b;
1073 c1 = ip[tp].vsite.c;
1074 constr_vsite4FD<calculatePosition, calculateVelocity>(x[ai],
1087 getVOrNull(avsite));
1093 b1 = ip[tp].vsite.b;
1094 c1 = ip[tp].vsite.c;
1095 constr_vsite4FDN<calculatePosition, calculateVelocity>(x[ai],
1108 getVOrNull(avsite));
1111 inc = constr_vsiten<calculatePosition, calculateVelocity>(ia, ip, x, pbc_null2, v);
1114 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d", ftype, __FILE__, __LINE__);
1117 if (pbcMode == PbcMode::all)
1119 /* Keep the vsite in the same periodic image as before */
1121 int ishift = pbc_dx_aiuc(pbc_null, x[avsite], xv, dx);
1122 if (ishift != CENTRAL)
1124 rvec_add(xv, dx, x[avsite]);
1128 /* Increment loop variables */
1136 /*! \brief Dispatch the vsite construction tasks for all threads
1138 * \param[in] threadingInfo Used to divide work over threads when != nullptr
1139 * \param[in,out] x Coordinates to construct vsites for
1140 * \param[in,out] v When not empty, velocities are generated for virtual sites
1141 * \param[in] ip Interaction parameters for all interaction, only vsite parameters are used
1142 * \param[in] ilist The interaction lists, only vsites are usesd
1143 * \param[in] domainInfo Information about PBC and DD
1144 * \param[in] box Used for PBC when PBC is set in domainInfo
1146 template<VSiteCalculatePosition calculatePosition, VSiteCalculateVelocity calculateVelocity>
1147 static void construct_vsites(const ThreadingInfo* threadingInfo,
1150 ArrayRef<const t_iparams> ip,
1151 ArrayRef<const InteractionList> ilist,
1152 const DomainInfo& domainInfo,
1155 const bool useDomdec = domainInfo.useDomdec();
1157 t_pbc pbc, *pbc_null;
1159 /* We only need to do pbc when we have inter update-group vsites.
1160 * Note that with domain decomposition we do not need to apply PBC here
1161 * when we have at least 3 domains along each dimension. Currently we
1162 * do not optimize this case.
1164 if (domainInfo.pbcType_ != PbcType::No && domainInfo.useMolPbc_)
1166 /* This is wasting some CPU time as we now do this multiple times
1170 clear_ivec(null_ivec);
1171 pbc_null = set_pbc_dd(
1172 &pbc, domainInfo.pbcType_, useDomdec ? domainInfo.domdec_->numCells : null_ivec, FALSE, box);
1181 if (calculateVelocity == VSiteCalculateVelocity::Yes)
1183 dd_move_x_and_v_vsites(
1184 *domainInfo.domdec_, box, as_rvec_array(x.data()), as_rvec_array(v.data()));
1188 dd_move_x_vsites(*domainInfo.domdec_, box, as_rvec_array(x.data()));
1192 if (threadingInfo == nullptr || threadingInfo->numThreads() == 1)
1194 construct_vsites_thread<calculatePosition, calculateVelocity>(x, v, ip, ilist, pbc_null);
1198 #pragma omp parallel num_threads(threadingInfo->numThreads())
1202 const int th = gmx_omp_get_thread_num();
1203 const VsiteThread& tData = threadingInfo->threadData(th);
1204 GMX_ASSERT(tData.rangeStart >= 0,
1205 "The thread data should be initialized before calling construct_vsites");
1207 construct_vsites_thread<calculatePosition, calculateVelocity>(
1208 x, v, ip, tData.ilist, pbc_null);
1209 if (tData.useInterdependentTask)
1211 /* Here we don't need a barrier (unlike the spreading),
1212 * since both tasks only construct vsites from particles,
1213 * or local vsites, not from non-local vsites.
1215 construct_vsites_thread<calculatePosition, calculateVelocity>(
1216 x, v, ip, tData.idTask.ilist, pbc_null);
1219 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1221 /* Now we can construct the vsites that might depend on other vsites */
1222 construct_vsites_thread<calculatePosition, calculateVelocity>(
1223 x, v, ip, threadingInfo->threadDataNonLocalDependent().ilist, pbc_null);
1227 void VirtualSitesHandler::Impl::construct(ArrayRef<RVec> x,
1230 VSiteOperation operation) const
1234 case VSiteOperation::Positions:
1235 construct_vsites<VSiteCalculatePosition::Yes, VSiteCalculateVelocity::No>(
1236 &threadingInfo_, x, v, iparams_, ilists_, domainInfo_, box);
1238 case VSiteOperation::Velocities:
1239 construct_vsites<VSiteCalculatePosition::No, VSiteCalculateVelocity::Yes>(
1240 &threadingInfo_, x, v, iparams_, ilists_, domainInfo_, box);
1242 case VSiteOperation::PositionsAndVelocities:
1243 construct_vsites<VSiteCalculatePosition::Yes, VSiteCalculateVelocity::Yes>(
1244 &threadingInfo_, x, v, iparams_, ilists_, domainInfo_, box);
1246 default: gmx_fatal(FARGS, "Unknown virtual site operation");
1250 void VirtualSitesHandler::construct(ArrayRef<RVec> x, ArrayRef<RVec> v, const matrix box, VSiteOperation operation) const
1252 impl_->construct(x, v, box, operation);
1255 void constructVirtualSites(ArrayRef<RVec> x, ArrayRef<const t_iparams> ip, ArrayRef<const InteractionList> ilist)
1259 const DomainInfo domainInfo;
1260 construct_vsites<VSiteCalculatePosition::Yes, VSiteCalculateVelocity::No>(
1261 nullptr, x, {}, ip, ilist, domainInfo, nullptr);
1265 /* Force spreading routines */
1267 static void spread_vsite1(const t_iatom ia[], ArrayRef<RVec> f)
1269 const int av = ia[1];
1270 const int ai = ia[2];
1275 template<VirialHandling virialHandling>
1276 static void spread_vsite2(const t_iatom ia[],
1278 ArrayRef<const RVec> x,
1280 ArrayRef<RVec> fshift,
1290 svmul(1 - a, f[av], fi);
1291 svmul(a, f[av], fj);
1294 rvec_inc(f[ai], fi);
1295 rvec_inc(f[aj], fj);
1298 if (virialHandling == VirialHandling::Pbc)
1304 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
1305 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
1313 if (siv != CENTRAL || sij != CENTRAL)
1315 rvec_inc(fshift[siv], f[av]);
1316 rvec_dec(fshift[CENTRAL], fi);
1317 rvec_dec(fshift[sij], fj);
1321 /* TOTAL: 13 flops */
1324 void constructVirtualSitesGlobal(const gmx_mtop_t& mtop, gmx::ArrayRef<gmx::RVec> x)
1326 GMX_ASSERT(x.ssize() >= mtop.natoms, "x should contain the whole system");
1327 GMX_ASSERT(!mtop.moleculeBlockIndices.empty(),
1328 "molblock indices are needed in constructVsitesGlobal");
1330 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
1332 const gmx_molblock_t& molb = mtop.molblock[mb];
1333 const gmx_moltype_t& molt = mtop.moltype[molb.type];
1334 if (vsiteIlistNrCount(molt.ilist) > 0)
1336 int atomOffset = mtop.moleculeBlockIndices[mb].globalAtomStart;
1337 for (int mol = 0; mol < molb.nmol; mol++)
1339 constructVirtualSites(
1340 x.subArray(atomOffset, molt.atoms.nr), mtop.ffparams.iparams, molt.ilist);
1341 atomOffset += molt.atoms.nr;
1347 template<VirialHandling virialHandling>
1348 static void spread_vsite2FD(const t_iatom ia[],
1350 ArrayRef<const RVec> x,
1352 ArrayRef<RVec> fshift,
1356 const int av = ia[1];
1357 const int ai = ia[2];
1358 const int aj = ia[3];
1360 copy_rvec(f[av], fv);
1363 int sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1366 const real invDistance = inverseNorm(xij);
1367 const real b = a * invDistance;
1368 /* 4 + ?10? flops */
1370 const real fproj = iprod(xij, fv) * invDistance * invDistance;
1373 fj[XX] = b * (fv[XX] - fproj * xij[XX]);
1374 fj[YY] = b * (fv[YY] - fproj * xij[YY]);
1375 fj[ZZ] = b * (fv[ZZ] - fproj * xij[ZZ]);
1378 /* b is already calculated in constr_vsite2FD
1379 storing b somewhere will save flops. */
1381 f[ai][XX] += fv[XX] - fj[XX];
1382 f[ai][YY] += fv[YY] - fj[YY];
1383 f[ai][ZZ] += fv[ZZ] - fj[ZZ];
1384 f[aj][XX] += fj[XX];
1385 f[aj][YY] += fj[YY];
1386 f[aj][ZZ] += fj[ZZ];
1389 if (virialHandling == VirialHandling::Pbc)
1395 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1402 if (svi != CENTRAL || sji != CENTRAL)
1404 rvec_dec(fshift[svi], fv);
1405 fshift[CENTRAL][XX] += fv[XX] - fj[XX];
1406 fshift[CENTRAL][YY] += fv[YY] - fj[YY];
1407 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ];
1408 fshift[sji][XX] += fj[XX];
1409 fshift[sji][YY] += fj[YY];
1410 fshift[sji][ZZ] += fj[ZZ];
1414 if (virialHandling == VirialHandling::NonLinear)
1416 /* Under this condition, the virial for the current forces is not
1417 * calculated from the redistributed forces. This means that
1418 * the effect of non-linear virtual site constructions on the virial
1419 * needs to be added separately. This contribution can be calculated
1420 * in many ways, but the simplest and cheapest way is to use
1421 * the first constructing atom ai as a reference position in space:
1422 * subtract (xv-xi)*fv and add (xj-xi)*fj.
1426 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1428 for (int i = 0; i < DIM; i++)
1430 for (int j = 0; j < DIM; j++)
1432 /* As xix is a linear combination of j and k, use that here */
1433 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * fj[j];
1438 /* TOTAL: 38 flops */
1441 template<VirialHandling virialHandling>
1442 static void spread_vsite3(const t_iatom ia[],
1445 ArrayRef<const RVec> x,
1447 ArrayRef<RVec> fshift,
1450 rvec fi, fj, fk, dx;
1458 svmul(1 - a - b, f[av], fi);
1459 svmul(a, f[av], fj);
1460 svmul(b, f[av], fk);
1463 rvec_inc(f[ai], fi);
1464 rvec_inc(f[aj], fj);
1465 rvec_inc(f[ak], fk);
1468 if (virialHandling == VirialHandling::Pbc)
1475 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
1476 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
1477 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
1486 if (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL)
1488 rvec_inc(fshift[siv], f[av]);
1489 rvec_dec(fshift[CENTRAL], fi);
1490 rvec_dec(fshift[sij], fj);
1491 rvec_dec(fshift[sik], fk);
1495 /* TOTAL: 20 flops */
1498 template<VirialHandling virialHandling>
1499 static void spread_vsite3FD(const t_iatom ia[],
1502 ArrayRef<const RVec> x,
1504 ArrayRef<RVec> fshift,
1509 rvec xvi, xij, xjk, xix, fv, temp;
1510 t_iatom av, ai, aj, ak;
1517 copy_rvec(f[av], fv);
1519 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1520 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1523 /* xix goes from i to point x on the line jk */
1524 xix[XX] = xij[XX] + a * xjk[XX];
1525 xix[YY] = xij[YY] + a * xjk[YY];
1526 xix[ZZ] = xij[ZZ] + a * xjk[ZZ];
1529 const real invDistance = inverseNorm(xix);
1530 const real c = b * invDistance;
1531 /* 4 + ?10? flops */
1533 fproj = iprod(xix, fv) * invDistance * invDistance; /* = (xix . f)/(xix . xix) */
1535 temp[XX] = c * (fv[XX] - fproj * xix[XX]);
1536 temp[YY] = c * (fv[YY] - fproj * xix[YY]);
1537 temp[ZZ] = c * (fv[ZZ] - fproj * xix[ZZ]);
1540 /* c is already calculated in constr_vsite3FD
1541 storing c somewhere will save 26 flops! */
1544 f[ai][XX] += fv[XX] - temp[XX];
1545 f[ai][YY] += fv[YY] - temp[YY];
1546 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
1547 f[aj][XX] += a1 * temp[XX];
1548 f[aj][YY] += a1 * temp[YY];
1549 f[aj][ZZ] += a1 * temp[ZZ];
1550 f[ak][XX] += a * temp[XX];
1551 f[ak][YY] += a * temp[YY];
1552 f[ak][ZZ] += a * temp[ZZ];
1555 if (virialHandling == VirialHandling::Pbc)
1560 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1567 if (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL)
1569 rvec_dec(fshift[svi], fv);
1570 fshift[CENTRAL][XX] += fv[XX] - (1 + a) * temp[XX];
1571 fshift[CENTRAL][YY] += fv[YY] - (1 + a) * temp[YY];
1572 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a) * temp[ZZ];
1573 fshift[sji][XX] += temp[XX];
1574 fshift[sji][YY] += temp[YY];
1575 fshift[sji][ZZ] += temp[ZZ];
1576 fshift[skj][XX] += a * temp[XX];
1577 fshift[skj][YY] += a * temp[YY];
1578 fshift[skj][ZZ] += a * temp[ZZ];
1582 if (virialHandling == VirialHandling::NonLinear)
1584 /* Under this condition, the virial for the current forces is not
1585 * calculated from the redistributed forces. This means that
1586 * the effect of non-linear virtual site constructions on the virial
1587 * needs to be added separately. This contribution can be calculated
1588 * in many ways, but the simplest and cheapest way is to use
1589 * the first constructing atom ai as a reference position in space:
1590 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
1594 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1596 for (int i = 0; i < DIM; i++)
1598 for (int j = 0; j < DIM; j++)
1600 /* As xix is a linear combination of j and k, use that here */
1601 dxdf[i][j] += -xiv[i] * fv[j] + xix[i] * temp[j];
1606 /* TOTAL: 61 flops */
1609 template<VirialHandling virialHandling>
1610 static void spread_vsite3FAD(const t_iatom ia[],
1613 ArrayRef<const RVec> x,
1615 ArrayRef<RVec> fshift,
1619 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
1620 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
1621 t_iatom av, ai, aj, ak;
1628 copy_rvec(f[ia[1]], fv);
1630 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1631 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1634 invdij = inverseNorm(xij);
1635 invdij2 = invdij * invdij;
1636 c1 = iprod(xij, xjk) * invdij2;
1637 xperp[XX] = xjk[XX] - c1 * xij[XX];
1638 xperp[YY] = xjk[YY] - c1 * xij[YY];
1639 xperp[ZZ] = xjk[ZZ] - c1 * xij[ZZ];
1640 /* xperp in plane ijk, perp. to ij */
1641 invdp = inverseNorm(xperp);
1646 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
1647 storing them somewhere will save 45 flops! */
1649 fproj = iprod(xij, fv) * invdij2;
1650 svmul(fproj, xij, Fpij); /* proj. f on xij */
1651 svmul(iprod(xperp, fv) * invdp * invdp, xperp, Fppp); /* proj. f on xperp */
1652 svmul(b1 * fproj, xperp, f3);
1655 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
1656 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
1657 for (int d = 0; d < DIM; d++)
1665 f[ai][XX] += fv[XX] - f1[XX] + c1 * f2[XX] + f3[XX];
1666 f[ai][YY] += fv[YY] - f1[YY] + c1 * f2[YY] + f3[YY];
1667 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1 * f2[ZZ] + f3[ZZ];
1668 f[aj][XX] += f1[XX] - c2 * f2[XX] - f3[XX];
1669 f[aj][YY] += f1[YY] - c2 * f2[YY] - f3[YY];
1670 f[aj][ZZ] += f1[ZZ] - c2 * f2[ZZ] - f3[ZZ];
1671 f[ak][XX] += f2[XX];
1672 f[ak][YY] += f2[YY];
1673 f[ak][ZZ] += f2[ZZ];
1676 if (virialHandling == VirialHandling::Pbc)
1682 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1689 if (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL)
1691 rvec_dec(fshift[svi], fv);
1692 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1 - c1) * f2[XX] + f3[XX];
1693 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1 - c1) * f2[YY] + f3[YY];
1694 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1 - c1) * f2[ZZ] + f3[ZZ];
1695 fshift[sji][XX] += f1[XX] - c1 * f2[XX] - f3[XX];
1696 fshift[sji][YY] += f1[YY] - c1 * f2[YY] - f3[YY];
1697 fshift[sji][ZZ] += f1[ZZ] - c1 * f2[ZZ] - f3[ZZ];
1698 fshift[skj][XX] += f2[XX];
1699 fshift[skj][YY] += f2[YY];
1700 fshift[skj][ZZ] += f2[ZZ];
1704 if (virialHandling == VirialHandling::NonLinear)
1707 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1709 for (int i = 0; i < DIM; i++)
1711 for (int j = 0; j < DIM; j++)
1713 /* Note that xik=xij+xjk, so we have to add xij*f2 */
1714 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * (f1[j] + (1 - c2) * f2[j] - f3[j])
1720 /* TOTAL: 113 flops */
1723 template<VirialHandling virialHandling>
1724 static void spread_vsite3OUT(const t_iatom ia[],
1728 ArrayRef<const RVec> x,
1730 ArrayRef<RVec> fshift,
1734 rvec xvi, xij, xik, fv, fj, fk;
1744 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1745 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1748 copy_rvec(f[av], fv);
1755 fj[XX] = a * fv[XX] - xik[ZZ] * cfy + xik[YY] * cfz;
1756 fj[YY] = xik[ZZ] * cfx + a * fv[YY] - xik[XX] * cfz;
1757 fj[ZZ] = -xik[YY] * cfx + xik[XX] * cfy + a * fv[ZZ];
1759 fk[XX] = b * fv[XX] + xij[ZZ] * cfy - xij[YY] * cfz;
1760 fk[YY] = -xij[ZZ] * cfx + b * fv[YY] + xij[XX] * cfz;
1761 fk[ZZ] = xij[YY] * cfx - xij[XX] * cfy + b * fv[ZZ];
1764 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
1765 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
1766 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1767 rvec_inc(f[aj], fj);
1768 rvec_inc(f[ak], fk);
1771 if (virialHandling == VirialHandling::Pbc)
1776 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1783 if (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL)
1785 rvec_dec(fshift[svi], fv);
1786 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
1787 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
1788 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1789 rvec_inc(fshift[sji], fj);
1790 rvec_inc(fshift[ski], fk);
1794 if (virialHandling == VirialHandling::NonLinear)
1798 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1800 for (int i = 0; i < DIM; i++)
1802 for (int j = 0; j < DIM; j++)
1804 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * fj[j] + xik[i] * fk[j];
1809 /* TOTAL: 54 flops */
1812 template<VirialHandling virialHandling>
1813 static void spread_vsite4FD(const t_iatom ia[],
1817 ArrayRef<const RVec> x,
1819 ArrayRef<RVec> fshift,
1824 rvec xvi, xij, xjk, xjl, xix, fv, temp;
1825 int av, ai, aj, ak, al;
1826 int sji, skj, slj, m;
1834 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1835 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1836 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1839 /* xix goes from i to point x on the plane jkl */
1840 for (m = 0; m < DIM; m++)
1842 xix[m] = xij[m] + a * xjk[m] + b * xjl[m];
1846 const real invDistance = inverseNorm(xix);
1847 const real d = c * invDistance;
1848 /* 4 + ?10? flops */
1850 copy_rvec(f[av], fv);
1852 fproj = iprod(xix, fv) * invDistance * invDistance; /* = (xix . f)/(xix . xix) */
1854 for (m = 0; m < DIM; m++)
1856 temp[m] = d * (fv[m] - fproj * xix[m]);
1860 /* c is already calculated in constr_vsite3FD
1861 storing c somewhere will save 35 flops! */
1864 for (m = 0; m < DIM; m++)
1866 f[ai][m] += fv[m] - temp[m];
1867 f[aj][m] += a1 * temp[m];
1868 f[ak][m] += a * temp[m];
1869 f[al][m] += b * temp[m];
1873 if (virialHandling == VirialHandling::Pbc)
1878 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1885 if (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL)
1887 rvec_dec(fshift[svi], fv);
1888 for (m = 0; m < DIM; m++)
1890 fshift[CENTRAL][m] += fv[m] - (1 + a + b) * temp[m];
1891 fshift[sji][m] += temp[m];
1892 fshift[skj][m] += a * temp[m];
1893 fshift[slj][m] += b * temp[m];
1898 if (virialHandling == VirialHandling::NonLinear)
1903 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1905 for (i = 0; i < DIM; i++)
1907 for (j = 0; j < DIM; j++)
1909 dxdf[i][j] += -xiv[i] * fv[j] + xix[i] * temp[j];
1914 /* TOTAL: 77 flops */
1917 template<VirialHandling virialHandling>
1918 static void spread_vsite4FDN(const t_iatom ia[],
1922 ArrayRef<const RVec> x,
1924 ArrayRef<RVec> fshift,
1928 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1929 rvec fv, fj, fk, fl;
1932 int av, ai, aj, ak, al;
1935 /* DEBUG: check atom indices */
1942 copy_rvec(f[av], fv);
1944 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1945 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1946 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1949 ra[XX] = a * xik[XX];
1950 ra[YY] = a * xik[YY];
1951 ra[ZZ] = a * xik[ZZ];
1953 rb[XX] = b * xil[XX];
1954 rb[YY] = b * xil[YY];
1955 rb[ZZ] = b * xil[ZZ];
1959 rvec_sub(ra, xij, rja);
1960 rvec_sub(rb, xij, rjb);
1961 rvec_sub(rb, ra, rab);
1964 cprod(rja, rjb, rm);
1967 invrm = inverseNorm(rm);
1968 denom = invrm * invrm;
1971 cfx = c * invrm * fv[XX];
1972 cfy = c * invrm * fv[YY];
1973 cfz = c * invrm * fv[ZZ];
1984 fj[XX] = (-rm[XX] * rt[XX]) * cfx + (rab[ZZ] - rm[YY] * rt[XX]) * cfy
1985 + (-rab[YY] - rm[ZZ] * rt[XX]) * cfz;
1986 fj[YY] = (-rab[ZZ] - rm[XX] * rt[YY]) * cfx + (-rm[YY] * rt[YY]) * cfy
1987 + (rab[XX] - rm[ZZ] * rt[YY]) * cfz;
1988 fj[ZZ] = (rab[YY] - rm[XX] * rt[ZZ]) * cfx + (-rab[XX] - rm[YY] * rt[ZZ]) * cfy
1989 + (-rm[ZZ] * rt[ZZ]) * cfz;
1995 rt[XX] *= denom * a;
1996 rt[YY] *= denom * a;
1997 rt[ZZ] *= denom * a;
2000 fk[XX] = (-rm[XX] * rt[XX]) * cfx + (-a * rjb[ZZ] - rm[YY] * rt[XX]) * cfy
2001 + (a * rjb[YY] - rm[ZZ] * rt[XX]) * cfz;
2002 fk[YY] = (a * rjb[ZZ] - rm[XX] * rt[YY]) * cfx + (-rm[YY] * rt[YY]) * cfy
2003 + (-a * rjb[XX] - rm[ZZ] * rt[YY]) * cfz;
2004 fk[ZZ] = (-a * rjb[YY] - rm[XX] * rt[ZZ]) * cfx + (a * rjb[XX] - rm[YY] * rt[ZZ]) * cfy
2005 + (-rm[ZZ] * rt[ZZ]) * cfz;
2011 rt[XX] *= denom * b;
2012 rt[YY] *= denom * b;
2013 rt[ZZ] *= denom * b;
2016 fl[XX] = (-rm[XX] * rt[XX]) * cfx + (b * rja[ZZ] - rm[YY] * rt[XX]) * cfy
2017 + (-b * rja[YY] - rm[ZZ] * rt[XX]) * cfz;
2018 fl[YY] = (-b * rja[ZZ] - rm[XX] * rt[YY]) * cfx + (-rm[YY] * rt[YY]) * cfy
2019 + (b * rja[XX] - rm[ZZ] * rt[YY]) * cfz;
2020 fl[ZZ] = (b * rja[YY] - rm[XX] * rt[ZZ]) * cfx + (-b * rja[XX] - rm[YY] * rt[ZZ]) * cfy
2021 + (-rm[ZZ] * rt[ZZ]) * cfz;
2024 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
2025 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
2026 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
2027 rvec_inc(f[aj], fj);
2028 rvec_inc(f[ak], fk);
2029 rvec_inc(f[al], fl);
2032 if (virialHandling == VirialHandling::Pbc)
2037 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
2044 if (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL)
2046 rvec_dec(fshift[svi], fv);
2047 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
2048 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
2049 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
2050 rvec_inc(fshift[sij], fj);
2051 rvec_inc(fshift[sik], fk);
2052 rvec_inc(fshift[sil], fl);
2056 if (virialHandling == VirialHandling::NonLinear)
2061 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
2063 for (i = 0; i < DIM; i++)
2065 for (j = 0; j < DIM; j++)
2067 dxdf[i][j] += -xiv[i] * fv[j] + xij[i] * fj[j] + xik[i] * fk[j] + xil[i] * fl[j];
2072 /* Total: 207 flops (Yuck!) */
2075 template<VirialHandling virialHandling>
2076 static int spread_vsiten(const t_iatom ia[],
2077 ArrayRef<const t_iparams> ip,
2078 ArrayRef<const RVec> x,
2080 ArrayRef<RVec> fshift,
2088 n3 = 3 * ip[ia[0]].vsiten.n;
2090 copy_rvec(x[av], xv);
2092 for (i = 0; i < n3; i += 3)
2097 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
2103 a = ip[ia[i]].vsiten.a;
2104 svmul(a, f[av], fi);
2105 rvec_inc(f[ai], fi);
2107 if (virialHandling == VirialHandling::Pbc && siv != CENTRAL)
2109 rvec_inc(fshift[siv], fi);
2110 rvec_dec(fshift[CENTRAL], fi);
2120 //! Returns the number of virtual sites in the interaction list, for VSITEN the number of atoms
2121 static int vsite_count(ArrayRef<const InteractionList> ilist, int ftype)
2123 if (ftype == F_VSITEN)
2125 return ilist[ftype].size() / 3;
2129 return ilist[ftype].size() / (1 + interaction_function[ftype].nratoms);
2133 //! Executes the force spreading task for a single thread
2134 template<VirialHandling virialHandling>
2135 static void spreadForceForThread(ArrayRef<const RVec> x,
2137 ArrayRef<RVec> fshift,
2139 ArrayRef<const t_iparams> ip,
2140 ArrayRef<const InteractionList> ilist,
2141 const t_pbc* pbc_null)
2143 const PbcMode pbcMode = getPbcMode(pbc_null);
2144 /* We need another pbc pointer, as with charge groups we switch per vsite */
2145 const t_pbc* pbc_null2 = pbc_null;
2146 gmx::ArrayRef<const int> vsite_pbc;
2148 /* this loop goes backwards to be able to build *
2149 * higher type vsites from lower types */
2150 for (int ftype = c_ftypeVsiteEnd - 1; ftype >= c_ftypeVsiteStart; ftype--)
2152 if (ilist[ftype].empty())
2158 int nra = interaction_function[ftype].nratoms;
2160 int nr = ilist[ftype].size();
2162 const t_iatom* ia = ilist[ftype].iatoms.data();
2164 if (pbcMode == PbcMode::all)
2166 pbc_null2 = pbc_null;
2169 for (int i = 0; i < nr;)
2173 /* Constants for constructing */
2175 a1 = ip[tp].vsite.a;
2176 /* Construct the vsite depending on type */
2179 case F_VSITE1: spread_vsite1(ia, f); break;
2181 spread_vsite2<virialHandling>(ia, a1, x, f, fshift, pbc_null2);
2184 spread_vsite2FD<virialHandling>(ia, a1, x, f, fshift, dxdf, pbc_null2);
2187 b1 = ip[tp].vsite.b;
2188 spread_vsite3<virialHandling>(ia, a1, b1, x, f, fshift, pbc_null2);
2191 b1 = ip[tp].vsite.b;
2192 spread_vsite3FD<virialHandling>(ia, a1, b1, x, f, fshift, dxdf, pbc_null2);
2195 b1 = ip[tp].vsite.b;
2196 spread_vsite3FAD<virialHandling>(ia, a1, b1, x, f, fshift, dxdf, pbc_null2);
2199 b1 = ip[tp].vsite.b;
2200 c1 = ip[tp].vsite.c;
2201 spread_vsite3OUT<virialHandling>(ia, a1, b1, c1, x, f, fshift, dxdf, pbc_null2);
2204 b1 = ip[tp].vsite.b;
2205 c1 = ip[tp].vsite.c;
2206 spread_vsite4FD<virialHandling>(ia, a1, b1, c1, x, f, fshift, dxdf, pbc_null2);
2209 b1 = ip[tp].vsite.b;
2210 c1 = ip[tp].vsite.c;
2211 spread_vsite4FDN<virialHandling>(ia, a1, b1, c1, x, f, fshift, dxdf, pbc_null2);
2214 inc = spread_vsiten<virialHandling>(ia, ip, x, f, fshift, pbc_null2);
2217 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d", ftype, __FILE__, __LINE__);
2219 clear_rvec(f[ia[1]]);
2221 /* Increment loop variables */
2229 //! Wrapper function for calling the templated thread-local spread function
2230 static void spreadForceWrapper(ArrayRef<const RVec> x,
2232 const VirialHandling virialHandling,
2233 ArrayRef<RVec> fshift,
2235 const bool clearDxdf,
2236 ArrayRef<const t_iparams> ip,
2237 ArrayRef<const InteractionList> ilist,
2238 const t_pbc* pbc_null)
2240 if (virialHandling == VirialHandling::NonLinear && clearDxdf)
2245 switch (virialHandling)
2247 case VirialHandling::None:
2248 spreadForceForThread<VirialHandling::None>(x, f, fshift, dxdf, ip, ilist, pbc_null);
2250 case VirialHandling::Pbc:
2251 spreadForceForThread<VirialHandling::Pbc>(x, f, fshift, dxdf, ip, ilist, pbc_null);
2253 case VirialHandling::NonLinear:
2254 spreadForceForThread<VirialHandling::NonLinear>(x, f, fshift, dxdf, ip, ilist, pbc_null);
2259 //! Clears the task force buffer elements that are written by task idTask
2260 static void clearTaskForceBufferUsedElements(InterdependentTask* idTask)
2262 int ntask = idTask->spreadTask.size();
2263 for (int ti = 0; ti < ntask; ti++)
2265 const AtomIndex* atomList = &idTask->atomIndex[idTask->spreadTask[ti]];
2266 int natom = atomList->atom.size();
2267 RVec* force = idTask->force.data();
2268 for (int i = 0; i < natom; i++)
2270 clear_rvec(force[atomList->atom[i]]);
2275 void VirtualSitesHandler::Impl::spreadForces(ArrayRef<const RVec> x,
2277 const VirialHandling virialHandling,
2278 ArrayRef<RVec> fshift,
2282 gmx_wallcycle* wcycle)
2284 wallcycle_start(wcycle, ewcVSITESPREAD);
2286 const bool useDomdec = domainInfo_.useDomdec();
2288 t_pbc pbc, *pbc_null;
2290 if (domainInfo_.useMolPbc_)
2292 /* This is wasting some CPU time as we now do this multiple times
2295 pbc_null = set_pbc_dd(
2296 &pbc, domainInfo_.pbcType_, useDomdec ? domainInfo_.domdec_->numCells : nullptr, FALSE, box);
2305 dd_clear_f_vsites(*domainInfo_.domdec_, f);
2308 const int numThreads = threadingInfo_.numThreads();
2310 if (numThreads == 1)
2313 spreadForceWrapper(x, f, virialHandling, fshift, dxdf, true, iparams_, ilists_, pbc_null);
2315 if (virialHandling == VirialHandling::NonLinear)
2317 for (int i = 0; i < DIM; i++)
2319 for (int j = 0; j < DIM; j++)
2321 virial[i][j] += -0.5 * dxdf[i][j];
2328 /* First spread the vsites that might depend on non-local vsites */
2329 auto& nlDependentVSites = threadingInfo_.threadDataNonLocalDependent();
2330 spreadForceWrapper(x,
2334 nlDependentVSites.dxdf,
2337 nlDependentVSites.ilist,
2340 #pragma omp parallel num_threads(numThreads)
2344 int thread = gmx_omp_get_thread_num();
2345 VsiteThread& tData = threadingInfo_.threadData(thread);
2347 ArrayRef<RVec> fshift_t;
2348 if (virialHandling == VirialHandling::Pbc)
2356 fshift_t = tData.fshift;
2358 for (int i = 0; i < SHIFTS; i++)
2360 clear_rvec(fshift_t[i]);
2365 if (tData.useInterdependentTask)
2367 /* Spread the vsites that spread outside our local range.
2368 * This is done using a thread-local force buffer force.
2369 * First we need to copy the input vsite forces to force.
2371 InterdependentTask* idTask = &tData.idTask;
2373 /* Clear the buffer elements set by our task during
2374 * the last call to spread_vsite_f.
2376 clearTaskForceBufferUsedElements(idTask);
2378 int nvsite = idTask->vsite.size();
2379 for (int i = 0; i < nvsite; i++)
2381 copy_rvec(f[idTask->vsite[i]], idTask->force[idTask->vsite[i]]);
2383 spreadForceWrapper(x,
2393 /* We need a barrier before reducing forces below
2394 * that have been produced by a different thread above.
2398 /* Loop over all thread task and reduce forces they
2399 * produced on atoms that fall in our range.
2400 * Note that atomic reduction would be a simpler solution,
2401 * but that might not have good support on all platforms.
2403 int ntask = idTask->reduceTask.size();
2404 for (int ti = 0; ti < ntask; ti++)
2406 const InterdependentTask& idt_foreign =
2407 threadingInfo_.threadData(idTask->reduceTask[ti]).idTask;
2408 const AtomIndex& atomList = idt_foreign.atomIndex[thread];
2409 const RVec* f_foreign = idt_foreign.force.data();
2411 for (int ind : atomList.atom)
2413 rvec_inc(f[ind], f_foreign[ind]);
2414 /* Clearing of f_foreign is done at the next step */
2417 /* Clear the vsite forces, both in f and force */
2418 for (int i = 0; i < nvsite; i++)
2420 int ind = tData.idTask.vsite[i];
2422 clear_rvec(tData.idTask.force[ind]);
2426 /* Spread the vsites that spread locally only */
2428 x, f, virialHandling, fshift_t, tData.dxdf, false, iparams_, tData.ilist, pbc_null);
2430 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2433 if (virialHandling == VirialHandling::Pbc)
2435 for (int th = 1; th < numThreads; th++)
2437 for (int i = 0; i < SHIFTS; i++)
2439 rvec_inc(fshift[i], threadingInfo_.threadData(th).fshift[i]);
2444 if (virialHandling == VirialHandling::NonLinear)
2446 for (int th = 0; th < numThreads + 1; th++)
2448 /* MSVC doesn't like matrix references, so we use a pointer */
2449 const matrix& dxdf = threadingInfo_.threadData(th).dxdf;
2451 for (int i = 0; i < DIM; i++)
2453 for (int j = 0; j < DIM; j++)
2455 virial[i][j] += -0.5 * dxdf[i][j];
2464 dd_move_f_vsites(*domainInfo_.domdec_, f, fshift);
2467 inc_nrnb(nrnb, eNR_VSITE1, vsite_count(ilists_, F_VSITE1));
2468 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(ilists_, F_VSITE2));
2469 inc_nrnb(nrnb, eNR_VSITE2FD, vsite_count(ilists_, F_VSITE2FD));
2470 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(ilists_, F_VSITE3));
2471 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(ilists_, F_VSITE3FD));
2472 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(ilists_, F_VSITE3FAD));
2473 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(ilists_, F_VSITE3OUT));
2474 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(ilists_, F_VSITE4FD));
2475 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(ilists_, F_VSITE4FDN));
2476 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(ilists_, F_VSITEN));
2478 wallcycle_stop(wcycle, ewcVSITESPREAD);
2481 /*! \brief Returns the an array with group indices for each atom
2483 * \param[in] grouping The paritioning of the atom range into atom groups
2485 static std::vector<int> makeAtomToGroupMapping(const gmx::RangePartitioning& grouping)
2487 std::vector<int> atomToGroup(grouping.fullRange().end(), 0);
2489 for (int group = 0; group < grouping.numBlocks(); group++)
2491 auto block = grouping.block(group);
2492 std::fill(atomToGroup.begin() + block.begin(), atomToGroup.begin() + block.end(), group);
2498 int countNonlinearVsites(const gmx_mtop_t& mtop)
2500 int numNonlinearVsites = 0;
2501 for (const gmx_molblock_t& molb : mtop.molblock)
2503 const gmx_moltype_t& molt = mtop.moltype[molb.type];
2505 for (const auto& ilist : extractILists(molt.ilist, IF_VSITE))
2507 if (ilist.functionType != F_VSITE2 && ilist.functionType != F_VSITE3
2508 && ilist.functionType != F_VSITEN)
2510 numNonlinearVsites += molb.nmol * ilist.iatoms.size() / (1 + NRAL(ilist.functionType));
2515 return numNonlinearVsites;
2518 void VirtualSitesHandler::spreadForces(ArrayRef<const RVec> x,
2520 const VirialHandling virialHandling,
2521 ArrayRef<RVec> fshift,
2525 gmx_wallcycle* wcycle)
2527 impl_->spreadForces(x, f, virialHandling, fshift, virial, nrnb, box, wcycle);
2530 int countInterUpdategroupVsites(const gmx_mtop_t& mtop,
2531 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingPerMoleculetype)
2533 int n_intercg_vsite = 0;
2534 for (const gmx_molblock_t& molb : mtop.molblock)
2536 const gmx_moltype_t& molt = mtop.moltype[molb.type];
2538 std::vector<int> atomToGroup;
2539 if (!updateGroupingPerMoleculetype.empty())
2541 atomToGroup = makeAtomToGroupMapping(updateGroupingPerMoleculetype[molb.type]);
2543 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2545 const int nral = NRAL(ftype);
2546 const InteractionList& il = molt.ilist[ftype];
2547 for (int i = 0; i < il.size(); i += 1 + nral)
2549 bool isInterGroup = atomToGroup.empty();
2552 const int group = atomToGroup[il.iatoms[1 + i]];
2553 for (int a = 1; a < nral; a++)
2555 if (atomToGroup[il.iatoms[1 + a]] != group)
2557 isInterGroup = true;
2564 n_intercg_vsite += molb.nmol;
2570 return n_intercg_vsite;
2573 std::unique_ptr<VirtualSitesHandler> makeVirtualSitesHandler(const gmx_mtop_t& mtop,
2574 const t_commrec* cr,
2577 GMX_RELEASE_ASSERT(cr != nullptr, "We need a valid commrec");
2579 std::unique_ptr<VirtualSitesHandler> vsite;
2581 /* check if there are vsites */
2583 for (int ftype = 0; ftype < F_NRE; ftype++)
2585 if (interaction_function[ftype].flags & IF_VSITE)
2587 GMX_ASSERT(ftype >= c_ftypeVsiteStart && ftype < c_ftypeVsiteEnd,
2588 "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
2590 nvsite += gmx_mtop_ftype_count(mtop, ftype);
2594 GMX_ASSERT(ftype < c_ftypeVsiteStart || ftype >= c_ftypeVsiteEnd,
2595 "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
2604 return std::make_unique<VirtualSitesHandler>(mtop, cr->dd, pbcType);
2607 ThreadingInfo::ThreadingInfo() : numThreads_(gmx_omp_nthreads_get(emntVSITE))
2609 if (numThreads_ > 1)
2611 /* We need one extra thread data structure for the overlap vsites */
2612 tData_.resize(numThreads_ + 1);
2613 #pragma omp parallel for num_threads(numThreads_) schedule(static)
2614 for (int thread = 0; thread < numThreads_; thread++)
2618 tData_[thread] = std::make_unique<VsiteThread>();
2620 InterdependentTask& idTask = tData_[thread]->idTask;
2622 idTask.atomIndex.resize(numThreads_);
2624 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2626 if (numThreads_ > 1)
2628 tData_[numThreads_] = std::make_unique<VsiteThread>();
2633 //! Returns the number of inter update-group vsites
2634 static int getNumInterUpdategroupVsites(const gmx_mtop_t& mtop, const gmx_domdec_t* domdec)
2636 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingPerMoleculetype;
2639 updateGroupingPerMoleculetype = getUpdateGroupingPerMoleculetype(*domdec);
2642 return countInterUpdategroupVsites(mtop, updateGroupingPerMoleculetype);
2645 VirtualSitesHandler::Impl::Impl(const gmx_mtop_t& mtop, gmx_domdec_t* domdec, const PbcType pbcType) :
2646 numInterUpdategroupVirtualSites_(getNumInterUpdategroupVsites(mtop, domdec)),
2647 domainInfo_({ pbcType, pbcType != PbcType::No && numInterUpdategroupVirtualSites_ > 0, domdec }),
2648 iparams_(mtop.ffparams.iparams)
2652 VirtualSitesHandler::VirtualSitesHandler(const gmx_mtop_t& mtop, gmx_domdec_t* domdec, const PbcType pbcType) :
2653 impl_(new Impl(mtop, domdec, pbcType))
2657 //! Flag that atom \p atom which is home in another task, if it has not already been added before
2658 static inline void flagAtom(InterdependentTask* idTask, const int atom, const int numThreads, const int numAtomsPerThread)
2660 if (!idTask->use[atom])
2662 idTask->use[atom] = true;
2663 int thread = atom / numAtomsPerThread;
2664 /* Assign all non-local atom force writes to thread 0 */
2665 if (thread >= numThreads)
2669 idTask->atomIndex[thread].atom.push_back(atom);
2673 /*! \brief Here we try to assign all vsites that are in our local range.
2675 * Our task local atom range is tData->rangeStart - tData->rangeEnd.
2676 * Vsites that depend only on local atoms, as indicated by taskIndex[]==thread,
2677 * are assigned to task tData->ilist. Vsites that depend on non-local atoms
2678 * but not on other vsites are assigned to task tData->id_task.ilist.
2679 * taskIndex[] is set for all vsites in our range, either to our local tasks
2680 * or to the single last task as taskIndex[]=2*nthreads.
2682 static void assignVsitesToThread(VsiteThread* tData,
2686 gmx::ArrayRef<int> taskIndex,
2687 ArrayRef<const InteractionList> ilist,
2688 ArrayRef<const t_iparams> ip,
2689 const ParticleType* ptype)
2691 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2693 tData->ilist[ftype].clear();
2694 tData->idTask.ilist[ftype].clear();
2696 const int nral1 = 1 + NRAL(ftype);
2697 const int* iat = ilist[ftype].iatoms.data();
2698 for (int i = 0; i < ilist[ftype].size();)
2700 /* Get the number of iatom entries in this virtual site.
2701 * The 3 below for F_VSITEN is from 1+NRAL(ftype)=3
2703 const int numIAtoms = (ftype == F_VSITEN ? ip[iat[i]].vsiten.n * 3 : nral1);
2705 if (iat[1 + i] < tData->rangeStart || iat[1 + i] >= tData->rangeEnd)
2707 /* This vsite belongs to a different thread */
2712 /* We would like to assign this vsite to task thread,
2713 * but it might depend on atoms outside the atom range of thread
2714 * or on another vsite not assigned to task thread.
2717 if (ftype != F_VSITEN)
2719 for (int j = i + 2; j < i + nral1; j++)
2721 /* Do a range check to avoid a harmless race on taskIndex */
2722 if (iat[j] < tData->rangeStart || iat[j] >= tData->rangeEnd || taskIndex[iat[j]] != thread)
2724 if (!tData->useInterdependentTask || ptype[iat[j]] == ParticleType::VSite)
2726 /* At least one constructing atom is a vsite
2727 * that is not assigned to the same thread.
2728 * Put this vsite into a separate task.
2734 /* There are constructing atoms outside our range,
2735 * put this vsite into a second task to be executed
2736 * on the same thread. During construction no barrier
2737 * is needed between the two tasks on the same thread.
2738 * During spreading we need to run this task with
2739 * an additional thread-local intermediate force buffer
2740 * (or atomic reduction) and a barrier between the two
2743 task = nthread + thread;
2749 for (int j = i + 2; j < i + numIAtoms; j += 3)
2751 /* Do a range check to avoid a harmless race on taskIndex */
2752 if (iat[j] < tData->rangeStart || iat[j] >= tData->rangeEnd || taskIndex[iat[j]] != thread)
2754 GMX_ASSERT(ptype[iat[j]] != ParticleType::VSite,
2755 "A vsite to be assigned in assignVsitesToThread has a vsite as "
2756 "a constructing atom that does not belong to our task, such "
2757 "vsites should be assigned to the single 'master' task");
2759 task = nthread + thread;
2764 /* Update this vsite's thread index entry */
2765 taskIndex[iat[1 + i]] = task;
2767 if (task == thread || task == nthread + thread)
2769 /* Copy this vsite to the thread data struct of thread */
2770 InteractionList* il_task;
2773 il_task = &tData->ilist[ftype];
2777 il_task = &tData->idTask.ilist[ftype];
2779 /* Copy the vsite data to the thread-task local array */
2780 il_task->push_back(iat[i], numIAtoms - 1, iat + i + 1);
2781 if (task == nthread + thread)
2783 /* This vsite writes outside our own task force block.
2784 * Put it into the interdependent task list and flag
2785 * the atoms involved for reduction.
2787 tData->idTask.vsite.push_back(iat[i + 1]);
2788 if (ftype != F_VSITEN)
2790 for (int j = i + 2; j < i + nral1; j++)
2792 flagAtom(&tData->idTask, iat[j], nthread, natperthread);
2797 for (int j = i + 2; j < i + numIAtoms; j += 3)
2799 flagAtom(&tData->idTask, iat[j], nthread, natperthread);
2810 /*! \brief Assign all vsites with taskIndex[]==task to task tData */
2811 static void assignVsitesToSingleTask(VsiteThread* tData,
2813 gmx::ArrayRef<const int> taskIndex,
2814 ArrayRef<const InteractionList> ilist,
2815 ArrayRef<const t_iparams> ip)
2817 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2819 tData->ilist[ftype].clear();
2820 tData->idTask.ilist[ftype].clear();
2822 int nral1 = 1 + NRAL(ftype);
2824 const int* iat = ilist[ftype].iatoms.data();
2825 InteractionList* il_task = &tData->ilist[ftype];
2827 for (int i = 0; i < ilist[ftype].size();)
2829 if (ftype == F_VSITEN)
2831 /* The 3 below is from 1+NRAL(ftype)=3 */
2832 inc = ip[iat[i]].vsiten.n * 3;
2834 /* Check if the vsite is assigned to our task */
2835 if (taskIndex[iat[1 + i]] == task)
2837 /* Copy the vsite data to the thread-task local array */
2838 il_task->push_back(iat[i], inc - 1, iat + i + 1);
2846 void ThreadingInfo::setVirtualSites(ArrayRef<const InteractionList> ilists,
2847 ArrayRef<const t_iparams> iparams,
2848 const t_mdatoms& mdatoms,
2849 const bool useDomdec)
2851 if (numThreads_ <= 1)
2857 /* The current way of distributing the vsites over threads in primitive.
2858 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
2859 * without taking into account how the vsites are distributed.
2860 * Without domain decomposition we at least tighten the upper bound
2861 * of the range (useful for common systems such as a vsite-protein
2863 * With domain decomposition, as long as the vsites are distributed
2864 * uniformly in each domain along the major dimension, usually x,
2865 * it will also perform well.
2867 int vsite_atom_range;
2871 vsite_atom_range = -1;
2872 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2875 if (ftype != F_VSITEN)
2877 int nral1 = 1 + NRAL(ftype);
2878 ArrayRef<const int> iat = ilists[ftype].iatoms;
2879 for (int i = 0; i < ilists[ftype].size(); i += nral1)
2881 for (int j = i + 1; j < i + nral1; j++)
2883 vsite_atom_range = std::max(vsite_atom_range, iat[j]);
2891 ArrayRef<const int> iat = ilists[ftype].iatoms;
2894 while (i < ilists[ftype].size())
2896 /* The 3 below is from 1+NRAL(ftype)=3 */
2897 vs_ind_end = i + iparams[iat[i]].vsiten.n * 3;
2899 vsite_atom_range = std::max(vsite_atom_range, iat[i + 1]);
2900 while (i < vs_ind_end)
2902 vsite_atom_range = std::max(vsite_atom_range, iat[i + 2]);
2910 natperthread = (vsite_atom_range + numThreads_ - 1) / numThreads_;
2914 /* Any local or not local atom could be involved in virtual sites.
2915 * But since we usually have very few non-local virtual sites
2916 * (only non-local vsites that depend on local vsites),
2917 * we distribute the local atom range equally over the threads.
2918 * When assigning vsites to threads, we should take care that the last
2919 * threads also covers the non-local range.
2921 vsite_atom_range = mdatoms.nr;
2922 natperthread = (mdatoms.homenr + numThreads_ - 1) / numThreads_;
2928 "virtual site thread dist: natoms %d, range %d, natperthread %d\n",
2934 /* To simplify the vsite assignment, we make an index which tells us
2935 * to which task particles, both non-vsites and vsites, are assigned.
2937 taskIndex_.resize(mdatoms.nr);
2939 /* Initialize the task index array. Here we assign the non-vsite
2940 * particles to task=thread, so we easily figure out if vsites
2941 * depend on local and/or non-local particles in assignVsitesToThread.
2945 for (int i = 0; i < mdatoms.nr; i++)
2947 if (mdatoms.ptype[i] == ParticleType::VSite)
2949 /* vsites are not assigned to a task yet */
2954 /* assign non-vsite particles to task thread */
2955 taskIndex_[i] = thread;
2957 if (i == (thread + 1) * natperthread && thread < numThreads_)
2964 #pragma omp parallel num_threads(numThreads_)
2968 int thread = gmx_omp_get_thread_num();
2969 VsiteThread& tData = *tData_[thread];
2971 /* Clear the buffer use flags that were set before */
2972 if (tData.useInterdependentTask)
2974 InterdependentTask& idTask = tData.idTask;
2976 /* To avoid an extra OpenMP barrier in spread_vsite_f,
2977 * we clear the force buffer at the next step,
2978 * so we need to do it here as well.
2980 clearTaskForceBufferUsedElements(&idTask);
2982 idTask.vsite.resize(0);
2983 for (int t = 0; t < numThreads_; t++)
2985 AtomIndex& atomIndex = idTask.atomIndex[t];
2986 int natom = atomIndex.atom.size();
2987 for (int i = 0; i < natom; i++)
2989 idTask.use[atomIndex.atom[i]] = false;
2991 atomIndex.atom.resize(0);
2996 /* To avoid large f_buf allocations of #threads*vsite_atom_range
2997 * we don't use task2 with more than 200000 atoms. This doesn't
2998 * affect performance, since with such a large range relatively few
2999 * vsites will end up in the separate task.
3000 * Note that useTask2 should be the same for all threads.
3002 tData.useInterdependentTask = (vsite_atom_range <= 200000);
3003 if (tData.useInterdependentTask)
3005 size_t natoms_use_in_vsites = vsite_atom_range;
3006 InterdependentTask& idTask = tData.idTask;
3007 /* To avoid resizing and re-clearing every nstlist steps,
3008 * we never down size the force buffer.
3010 if (natoms_use_in_vsites > idTask.force.size() || natoms_use_in_vsites > idTask.use.size())
3012 idTask.force.resize(natoms_use_in_vsites, { 0, 0, 0 });
3013 idTask.use.resize(natoms_use_in_vsites, false);
3017 /* Assign all vsites that can execute independently on threads */
3018 tData.rangeStart = thread * natperthread;
3019 if (thread < numThreads_ - 1)
3021 tData.rangeEnd = (thread + 1) * natperthread;
3025 /* The last thread should cover up to the end of the range */
3026 tData.rangeEnd = mdatoms.nr;
3028 assignVsitesToThread(
3029 &tData, thread, numThreads_, natperthread, taskIndex_, ilists, iparams, mdatoms.ptype);
3031 if (tData.useInterdependentTask)
3033 /* In the worst case, all tasks write to force ranges of
3034 * all other tasks, leading to #tasks^2 scaling (this is only
3035 * the overhead, the actual flops remain constant).
3036 * But in most cases there is far less coupling. To improve
3037 * scaling at high thread counts we therefore construct
3038 * an index to only loop over the actually affected tasks.
3040 InterdependentTask& idTask = tData.idTask;
3042 /* Ensure assignVsitesToThread finished on other threads */
3045 idTask.spreadTask.resize(0);
3046 idTask.reduceTask.resize(0);
3047 for (int t = 0; t < numThreads_; t++)
3049 /* Do we write to the force buffer of task t? */
3050 if (!idTask.atomIndex[t].atom.empty())
3052 idTask.spreadTask.push_back(t);
3054 /* Does task t write to our force buffer? */
3055 if (!tData_[t]->idTask.atomIndex[thread].atom.empty())
3057 idTask.reduceTask.push_back(t);
3062 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
3064 /* Assign all remaining vsites, that will have taskIndex[]=2*vsite->nthreads,
3065 * to a single task that will not run in parallel with other tasks.
3067 assignVsitesToSingleTask(tData_[numThreads_].get(), 2 * numThreads_, taskIndex_, ilists, iparams);
3069 if (debug && numThreads_ > 1)
3072 "virtual site useInterdependentTask %d, nuse:\n",
3073 static_cast<int>(tData_[0]->useInterdependentTask));
3074 for (int th = 0; th < numThreads_ + 1; th++)
3076 fprintf(debug, " %4d", tData_[th]->idTask.nuse);
3078 fprintf(debug, "\n");
3080 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
3082 if (!ilists[ftype].empty())
3084 fprintf(debug, "%-20s thread dist:", interaction_function[ftype].longname);
3085 for (int th = 0; th < numThreads_ + 1; th++)
3089 tData_[th]->ilist[ftype].size(),
3090 tData_[th]->idTask.ilist[ftype].size());
3092 fprintf(debug, "\n");
3098 int nrOrig = vsiteIlistNrCount(ilists);
3100 for (int th = 0; th < numThreads_ + 1; th++)
3102 nrThreaded += vsiteIlistNrCount(tData_[th]->ilist) + vsiteIlistNrCount(tData_[th]->idTask.ilist);
3104 GMX_ASSERT(nrThreaded == nrOrig,
3105 "The number of virtual sites assigned to all thread task has to match the total "
3106 "number of virtual sites");
3110 void VirtualSitesHandler::Impl::setVirtualSites(ArrayRef<const InteractionList> ilists,
3111 const t_mdatoms& mdatoms)
3115 threadingInfo_.setVirtualSites(ilists, iparams_, mdatoms, domainInfo_.useDomdec());
3118 void VirtualSitesHandler::setVirtualSites(ArrayRef<const InteractionList> ilists, const t_mdatoms& mdatoms)
3120 impl_->setVirtualSites(ilists, mdatoms);