Remove cppcheck
[gromacs.git] / src / gromacs / mdlib / vsite.cpp
blobe5085f769c7477f7fce9a6bdaae7add765494ee9
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "gmxpre.h"
39 #include "vsite.h"
41 #include <stdio.h>
43 #include <algorithm>
44 #include <vector>
46 #include "gromacs/domdec/domdec.h"
47 #include "gromacs/domdec/domdec_struct.h"
48 #include "gromacs/gmxlib/network.h"
49 #include "gromacs/gmxlib/nrnb.h"
50 #include "gromacs/math/functions.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdlib/gmx_omp_nthreads.h"
53 #include "gromacs/mdtypes/commrec.h"
54 #include "gromacs/mdtypes/mdatom.h"
55 #include "gromacs/pbcutil/ishift.h"
56 #include "gromacs/pbcutil/mshift.h"
57 #include "gromacs/pbcutil/pbc.h"
58 #include "gromacs/timing/wallcycle.h"
59 #include "gromacs/topology/ifunc.h"
60 #include "gromacs/topology/mtop_util.h"
61 #include "gromacs/topology/topology.h"
62 #include "gromacs/utility/exceptions.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/gmxassert.h"
65 #include "gromacs/utility/gmxomp.h"
66 #include "gromacs/utility/smalloc.h"
68 /* The strategy used here for assigning virtual sites to (thread-)tasks
69 * is as follows:
71 * We divide the atom range that vsites operate on (natoms_local with DD,
72 * 0 - last atom involved in vsites without DD) equally over all threads.
74 * Vsites in the local range constructed from atoms in the local range
75 * and/or other vsites that are fully local are assigned to a simple,
76 * independent task.
78 * Vsites that are not assigned after using the above criterion get assigned
79 * to a so called "interdependent" thread task when none of the constructing
80 * atoms is a vsite. These tasks are called interdependent, because one task
81 * accesses atoms assigned to a different task/thread.
82 * Note that this option is turned off with large (local) atom counts
83 * to avoid high memory usage.
85 * Any remaining vsites are assigned to a separate master thread task.
88 using gmx::RVec;
90 static void init_ilist(t_ilist *ilist)
92 for (int i = 0; i < F_NRE; i++)
94 ilist[i].nr = 0;
95 ilist[i].nalloc = 0;
96 ilist[i].iatoms = nullptr;
100 /*! \brief List of atom indices belonging to a task */
101 struct AtomIndex {
102 //! List of atom indices
103 std::vector<int> atom;
106 /*! \brief Data structure for thread tasks that use constructing atoms outside their own atom range */
107 struct InterdependentTask
109 //! The interaction lists, only vsite entries are used
110 t_ilist ilist[F_NRE];
111 //! Thread/task-local force buffer
112 std::vector<RVec> force;
113 //! The atom indices of the vsites of our task
114 std::vector<int> vsite;
115 //! Flags if elements in force are spread to or not
116 std::vector<bool> use;
117 //! The number of entries set to true in use
118 int nuse;
119 //! Array of atoms indices, size nthreads, covering all nuse set elements in use
120 std::vector<AtomIndex> atomIndex;
121 //! List of tasks (force blocks) this task spread forces to
122 std::vector<int> spreadTask;
123 //! List of tasks that write to this tasks force block range
124 std::vector<int> reduceTask;
126 InterdependentTask()
128 init_ilist(ilist);
129 nuse = 0;
133 /*! \brief Vsite thread task data structure */
134 struct VsiteThread {
135 //! Start of atom range of this task
136 int rangeStart;
137 //! End of atom range of this task
138 int rangeEnd;
139 //! The interaction lists, only vsite entries are used
140 t_ilist ilist[F_NRE];
141 //! Local fshift accumulation buffer
142 rvec fshift[SHIFTS];
143 //! Local virial dx*df accumulation buffer
144 matrix dxdf;
145 //! 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
146 bool useInterdependentTask;
147 //! Data for vsites that involve constructing atoms in the atom range of other threads/tasks
148 InterdependentTask idTask;
150 /*! \brief Constructor */
151 VsiteThread()
153 rangeStart = -1;
154 rangeEnd = -1;
155 init_ilist(ilist);
156 clear_rvecs(SHIFTS, fshift);
157 clear_mat(dxdf);
158 useInterdependentTask = false;
163 /* The start and end values of for the vsite indices in the ftype enum.
164 * The validity of these values is checked in init_vsite.
165 * This is used to avoid loops over all ftypes just to get the vsite entries.
166 * (We should replace the fixed ilist array by only the used entries.)
168 static const int c_ftypeVsiteStart = F_VSITE2;
169 static const int c_ftypeVsiteEnd = F_VSITEN + 1;
172 /*! \brief Returns the sum of the vsite ilist sizes over all vsite types
174 * \param[in] ilist The interaction list
176 static int vsiteIlistNrCount(const t_ilist *ilist)
178 int nr = 0;
179 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
181 nr += ilist[ftype].nr;
184 return nr;
187 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
189 if (pbc)
191 return pbc_dx_aiuc(pbc, xi, xj, dx);
193 else
195 rvec_sub(xi, xj, dx);
196 return CENTRAL;
200 /* Vsite construction routines */
202 static void constr_vsite2(const rvec xi, const rvec xj, rvec x, real a, const t_pbc *pbc)
204 real b = 1 - a;
205 /* 1 flop */
207 if (pbc)
209 rvec dx;
210 pbc_dx_aiuc(pbc, xj, xi, dx);
211 x[XX] = xi[XX] + a*dx[XX];
212 x[YY] = xi[YY] + a*dx[YY];
213 x[ZZ] = xi[ZZ] + a*dx[ZZ];
215 else
217 x[XX] = b*xi[XX] + a*xj[XX];
218 x[YY] = b*xi[YY] + a*xj[YY];
219 x[ZZ] = b*xi[ZZ] + a*xj[ZZ];
220 /* 9 Flops */
223 /* TOTAL: 10 flops */
226 static void constr_vsite3(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b,
227 const t_pbc *pbc)
229 real c = 1 - a - b;
230 /* 2 flops */
232 if (pbc)
234 rvec dxj, dxk;
236 pbc_dx_aiuc(pbc, xj, xi, dxj);
237 pbc_dx_aiuc(pbc, xk, xi, dxk);
238 x[XX] = xi[XX] + a*dxj[XX] + b*dxk[XX];
239 x[YY] = xi[YY] + a*dxj[YY] + b*dxk[YY];
240 x[ZZ] = xi[ZZ] + a*dxj[ZZ] + b*dxk[ZZ];
242 else
244 x[XX] = c*xi[XX] + a*xj[XX] + b*xk[XX];
245 x[YY] = c*xi[YY] + a*xj[YY] + b*xk[YY];
246 x[ZZ] = c*xi[ZZ] + a*xj[ZZ] + b*xk[ZZ];
247 /* 15 Flops */
250 /* TOTAL: 17 flops */
253 static void constr_vsite3FD(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b,
254 const t_pbc *pbc)
256 rvec xij, xjk, temp;
257 real c;
259 pbc_rvec_sub(pbc, xj, xi, xij);
260 pbc_rvec_sub(pbc, xk, xj, xjk);
261 /* 6 flops */
263 /* temp goes from i to a point on the line jk */
264 temp[XX] = xij[XX] + a*xjk[XX];
265 temp[YY] = xij[YY] + a*xjk[YY];
266 temp[ZZ] = xij[ZZ] + a*xjk[ZZ];
267 /* 6 flops */
269 c = b*gmx::invsqrt(iprod(temp, temp));
270 /* 6 + 10 flops */
272 x[XX] = xi[XX] + c*temp[XX];
273 x[YY] = xi[YY] + c*temp[YY];
274 x[ZZ] = xi[ZZ] + c*temp[ZZ];
275 /* 6 Flops */
277 /* TOTAL: 34 flops */
280 static void constr_vsite3FAD(const rvec xi, const rvec xj, const rvec xk, rvec x, real a, real b, const t_pbc *pbc)
282 rvec xij, xjk, xp;
283 real a1, b1, c1, invdij;
285 pbc_rvec_sub(pbc, xj, xi, xij);
286 pbc_rvec_sub(pbc, xk, xj, xjk);
287 /* 6 flops */
289 invdij = gmx::invsqrt(iprod(xij, xij));
290 c1 = invdij * invdij * iprod(xij, xjk);
291 xp[XX] = xjk[XX] - c1*xij[XX];
292 xp[YY] = xjk[YY] - c1*xij[YY];
293 xp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
294 a1 = a*invdij;
295 b1 = b*gmx::invsqrt(iprod(xp, xp));
296 /* 45 */
298 x[XX] = xi[XX] + a1*xij[XX] + b1*xp[XX];
299 x[YY] = xi[YY] + a1*xij[YY] + b1*xp[YY];
300 x[ZZ] = xi[ZZ] + a1*xij[ZZ] + b1*xp[ZZ];
301 /* 12 Flops */
303 /* TOTAL: 63 flops */
306 static void constr_vsite3OUT(const rvec xi, const rvec xj, const rvec xk, rvec x,
307 real a, real b, real c, const t_pbc *pbc)
309 rvec xij, xik, temp;
311 pbc_rvec_sub(pbc, xj, xi, xij);
312 pbc_rvec_sub(pbc, xk, xi, xik);
313 cprod(xij, xik, temp);
314 /* 15 Flops */
316 x[XX] = xi[XX] + a*xij[XX] + b*xik[XX] + c*temp[XX];
317 x[YY] = xi[YY] + a*xij[YY] + b*xik[YY] + c*temp[YY];
318 x[ZZ] = xi[ZZ] + a*xij[ZZ] + b*xik[ZZ] + c*temp[ZZ];
319 /* 18 Flops */
321 /* TOTAL: 33 flops */
324 static void constr_vsite4FD(const rvec xi, const rvec xj, const rvec xk, const rvec xl, rvec x,
325 real a, real b, real c, const t_pbc *pbc)
327 rvec xij, xjk, xjl, temp;
328 real d;
330 pbc_rvec_sub(pbc, xj, xi, xij);
331 pbc_rvec_sub(pbc, xk, xj, xjk);
332 pbc_rvec_sub(pbc, xl, xj, xjl);
333 /* 9 flops */
335 /* temp goes from i to a point on the plane jkl */
336 temp[XX] = xij[XX] + a*xjk[XX] + b*xjl[XX];
337 temp[YY] = xij[YY] + a*xjk[YY] + b*xjl[YY];
338 temp[ZZ] = xij[ZZ] + a*xjk[ZZ] + b*xjl[ZZ];
339 /* 12 flops */
341 d = c*gmx::invsqrt(iprod(temp, temp));
342 /* 6 + 10 flops */
344 x[XX] = xi[XX] + d*temp[XX];
345 x[YY] = xi[YY] + d*temp[YY];
346 x[ZZ] = xi[ZZ] + d*temp[ZZ];
347 /* 6 Flops */
349 /* TOTAL: 43 flops */
352 static void constr_vsite4FDN(const rvec xi, const rvec xj, const rvec xk, const rvec xl, rvec x,
353 real a, real b, real c, const t_pbc *pbc)
355 rvec xij, xik, xil, ra, rb, rja, rjb, rm;
356 real d;
358 pbc_rvec_sub(pbc, xj, xi, xij);
359 pbc_rvec_sub(pbc, xk, xi, xik);
360 pbc_rvec_sub(pbc, xl, xi, xil);
361 /* 9 flops */
363 ra[XX] = a*xik[XX];
364 ra[YY] = a*xik[YY];
365 ra[ZZ] = a*xik[ZZ];
367 rb[XX] = b*xil[XX];
368 rb[YY] = b*xil[YY];
369 rb[ZZ] = b*xil[ZZ];
371 /* 6 flops */
373 rvec_sub(ra, xij, rja);
374 rvec_sub(rb, xij, rjb);
375 /* 6 flops */
377 cprod(rja, rjb, rm);
378 /* 9 flops */
380 d = c*gmx::invsqrt(norm2(rm));
381 /* 5+5+1 flops */
383 x[XX] = xi[XX] + d*rm[XX];
384 x[YY] = xi[YY] + d*rm[YY];
385 x[ZZ] = xi[ZZ] + d*rm[ZZ];
386 /* 6 Flops */
388 /* TOTAL: 47 flops */
392 static int constr_vsiten(const t_iatom *ia, const t_iparams ip[],
393 rvec *x, const t_pbc *pbc)
395 rvec x1, dx;
396 dvec dsum;
397 int n3, av, ai;
398 real a;
400 n3 = 3*ip[ia[0]].vsiten.n;
401 av = ia[1];
402 ai = ia[2];
403 copy_rvec(x[ai], x1);
404 clear_dvec(dsum);
405 for (int i = 3; i < n3; i += 3)
407 ai = ia[i+2];
408 a = ip[ia[i]].vsiten.a;
409 if (pbc)
411 pbc_dx_aiuc(pbc, x[ai], x1, dx);
413 else
415 rvec_sub(x[ai], x1, dx);
417 dsum[XX] += a*dx[XX];
418 dsum[YY] += a*dx[YY];
419 dsum[ZZ] += a*dx[ZZ];
420 /* 9 Flops */
423 x[av][XX] = x1[XX] + dsum[XX];
424 x[av][YY] = x1[YY] + dsum[YY];
425 x[av][ZZ] = x1[ZZ] + dsum[ZZ];
427 return n3;
430 /*! \brief PBC modes for vsite construction and spreading */
431 enum class PbcMode
433 all, // Apply normal, simple PBC for all vsites
434 chargeGroup, // Keep vsite in the same periodic image as the rest of it's charge group
435 none // No PBC treatment needed
438 /*! \brief Returns the PBC mode based on the system PBC and vsite properties
440 * \param[in] pbcPtr A pointer to a PBC struct or nullptr when no PBC treatment is required
441 * \param[in] vsite A pointer to the vsite struct, can be nullptr
443 static PbcMode getPbcMode(const t_pbc *pbcPtr,
444 const gmx_vsite_t *vsite)
446 if (pbcPtr == nullptr)
448 return PbcMode::none;
450 else if (vsite != nullptr && vsite->bHaveChargeGroups)
452 return PbcMode::chargeGroup;
454 else
456 return PbcMode::all;
460 static void construct_vsites_thread(const gmx_vsite_t *vsite,
461 rvec x[],
462 real dt, rvec *v,
463 const t_iparams ip[], const t_ilist ilist[],
464 const t_pbc *pbc_null)
466 real inv_dt;
467 if (v != nullptr)
469 inv_dt = 1.0/dt;
471 else
473 inv_dt = 1.0;
476 const PbcMode pbcMode = getPbcMode(pbc_null, vsite);
477 /* We need another pbc pointer, as with charge groups we switch per vsite */
478 const t_pbc *pbc_null2 = pbc_null;
479 const int *vsite_pbc = nullptr;
481 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
483 if (ilist[ftype].nr == 0)
485 continue;
488 { // TODO remove me
489 int nra = interaction_function[ftype].nratoms;
490 int inc = 1 + nra;
491 int nr = ilist[ftype].nr;
493 const t_iatom *ia = ilist[ftype].iatoms;
495 if (pbcMode == PbcMode::chargeGroup)
497 vsite_pbc = vsite->vsite_pbc_loc[ftype - c_ftypeVsiteStart];
500 for (int i = 0; i < nr; )
502 int tp = ia[0];
503 /* The vsite and constructing atoms */
504 int avsite = ia[1];
505 int ai = ia[2];
506 /* Constants for constructing vsites */
507 real a1 = ip[tp].vsite.a;
508 /* Check what kind of pbc we need to use */
509 int pbc_atom;
510 rvec xpbc;
511 if (pbcMode == PbcMode::all)
513 /* No charge groups, vsite follows its own pbc */
514 pbc_atom = avsite;
515 copy_rvec(x[avsite], xpbc);
517 else if (pbcMode == PbcMode::chargeGroup)
519 pbc_atom = vsite_pbc[i/(1 + nra)];
520 if (pbc_atom > -2)
522 if (pbc_atom >= 0)
524 /* We need to copy the coordinates here,
525 * single for single atom cg's pbc_atom
526 * is the vsite itself.
528 copy_rvec(x[pbc_atom], xpbc);
530 pbc_null2 = pbc_null;
532 else
534 pbc_null2 = nullptr;
537 else
539 pbc_atom = -2;
541 /* Copy the old position */
542 rvec xv;
543 copy_rvec(x[avsite], xv);
545 /* Construct the vsite depending on type */
546 int aj, ak, al;
547 real b1, c1;
548 switch (ftype)
550 case F_VSITE2:
551 aj = ia[3];
552 constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
553 break;
554 case F_VSITE3:
555 aj = ia[3];
556 ak = ia[4];
557 b1 = ip[tp].vsite.b;
558 constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
559 break;
560 case F_VSITE3FD:
561 aj = ia[3];
562 ak = ia[4];
563 b1 = ip[tp].vsite.b;
564 constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
565 break;
566 case F_VSITE3FAD:
567 aj = ia[3];
568 ak = ia[4];
569 b1 = ip[tp].vsite.b;
570 constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
571 break;
572 case F_VSITE3OUT:
573 aj = ia[3];
574 ak = ia[4];
575 b1 = ip[tp].vsite.b;
576 c1 = ip[tp].vsite.c;
577 constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
578 break;
579 case F_VSITE4FD:
580 aj = ia[3];
581 ak = ia[4];
582 al = ia[5];
583 b1 = ip[tp].vsite.b;
584 c1 = ip[tp].vsite.c;
585 constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
586 pbc_null2);
587 break;
588 case F_VSITE4FDN:
589 aj = ia[3];
590 ak = ia[4];
591 al = ia[5];
592 b1 = ip[tp].vsite.b;
593 c1 = ip[tp].vsite.c;
594 constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
595 pbc_null2);
596 break;
597 case F_VSITEN:
598 inc = constr_vsiten(ia, ip, x, pbc_null2);
599 break;
600 default:
601 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
602 ftype, __FILE__, __LINE__);
605 if (pbc_atom >= 0)
607 /* Match the pbc of this vsite to the rest of its charge group */
608 rvec dx;
609 int ishift = pbc_dx_aiuc(pbc_null, x[avsite], xpbc, dx);
610 if (ishift != CENTRAL)
612 rvec_add(xpbc, dx, x[avsite]);
615 if (v != nullptr)
617 /* Calculate velocity of vsite... */
618 rvec vv;
619 rvec_sub(x[avsite], xv, vv);
620 svmul(inv_dt, vv, v[avsite]);
623 /* Increment loop variables */
624 i += inc;
625 ia += inc;
631 void construct_vsites(const gmx_vsite_t *vsite,
632 rvec x[],
633 real dt, rvec *v,
634 const t_iparams ip[], const t_ilist ilist[],
635 int ePBC, gmx_bool bMolPBC,
636 const t_commrec *cr,
637 const matrix box)
639 const bool useDomdec = (vsite != nullptr && vsite->useDomdec);
640 GMX_ASSERT(!useDomdec || (cr != nullptr && DOMAINDECOMP(cr)), "When vsites are set up with domain decomposition, we need a valid commrec");
641 // TODO: Remove this assertion when we remove charge groups
642 GMX_ASSERT(vsite != nullptr || ePBC == epbcNONE, "Without a vsite struct we can not do PBC (in case we have charge groups)");
644 t_pbc pbc, *pbc_null;
646 /* We only need to do pbc when we have inter-cg vsites.
647 * Note that with domain decomposition we do not need to apply PBC here
648 * when we have at least 3 domains along each dimension. Currently we
649 * do not optimize this case.
651 if (ePBC != epbcNONE && (useDomdec || bMolPBC) &&
652 !(vsite != nullptr && vsite->n_intercg_vsite == 0))
654 /* This is wasting some CPU time as we now do this multiple times
655 * per MD step.
657 ivec null_ivec;
658 clear_ivec(null_ivec);
659 pbc_null = set_pbc_dd(&pbc, ePBC,
660 useDomdec ? cr->dd->nc : null_ivec,
661 FALSE, box);
663 else
665 pbc_null = nullptr;
668 if (useDomdec)
670 dd_move_x_vsites(cr->dd, box, x);
673 if (vsite == nullptr || vsite->nthreads == 1)
675 construct_vsites_thread(vsite,
676 x, dt, v,
677 ip, ilist,
678 pbc_null);
680 else
682 #pragma omp parallel num_threads(vsite->nthreads)
686 const int th = gmx_omp_get_thread_num();
687 const VsiteThread &tData = *vsite->tData[th];
688 GMX_ASSERT(tData.rangeStart >= 0, "The thread data should be initialized before calling construct_vsites");
690 construct_vsites_thread(vsite,
691 x, dt, v,
692 ip, tData.ilist,
693 pbc_null);
694 if (tData.useInterdependentTask)
696 /* Here we don't need a barrier (unlike the spreading),
697 * since both tasks only construct vsites from particles,
698 * or local vsites, not from non-local vsites.
700 construct_vsites_thread(vsite,
701 x, dt, v,
702 ip, tData.idTask.ilist,
703 pbc_null);
706 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
708 /* Now we can construct the vsites that might depend on other vsites */
709 construct_vsites_thread(vsite,
710 x, dt, v,
711 ip, vsite->tData[vsite->nthreads]->ilist,
712 pbc_null);
716 static void spread_vsite2(const t_iatom ia[], real a,
717 const rvec x[],
718 rvec f[], rvec fshift[],
719 const t_pbc *pbc, const t_graph *g)
721 rvec fi, fj, dx;
722 t_iatom av, ai, aj;
723 ivec di;
724 int siv, sij;
726 av = ia[1];
727 ai = ia[2];
728 aj = ia[3];
730 svmul(1 - a, f[av], fi);
731 svmul( a, f[av], fj);
732 /* 7 flop */
734 rvec_inc(f[ai], fi);
735 rvec_inc(f[aj], fj);
736 /* 6 Flops */
738 if (g)
740 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
741 siv = IVEC2IS(di);
742 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
743 sij = IVEC2IS(di);
745 else if (pbc)
747 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
748 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
750 else
752 siv = CENTRAL;
753 sij = CENTRAL;
756 if (fshift && (siv != CENTRAL || sij != CENTRAL))
758 rvec_inc(fshift[siv], f[av]);
759 rvec_dec(fshift[CENTRAL], fi);
760 rvec_dec(fshift[sij], fj);
763 /* TOTAL: 13 flops */
766 void constructVsitesGlobal(const gmx_mtop_t &mtop,
767 gmx::ArrayRef<gmx::RVec> x)
769 GMX_ASSERT(x.size() >= static_cast<gmx::index>(mtop.natoms), "x should contain the whole system");
770 GMX_ASSERT(!mtop.moleculeBlockIndices.empty(), "molblock indices are needed in constructVsitesGlobal");
772 for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
774 const gmx_molblock_t &molb = mtop.molblock[mb];
775 const gmx_moltype_t &molt = mtop.moltype[molb.type];
776 if (vsiteIlistNrCount(molt.ilist) > 0)
778 int atomOffset = mtop.moleculeBlockIndices[mb].globalAtomStart;
779 for (int mol = 0; mol < molb.nmol; mol++)
781 construct_vsites(nullptr, as_rvec_array(x.data()) + atomOffset,
782 0.0, nullptr,
783 mtop.ffparams.iparams, molt.ilist,
784 epbcNONE, TRUE, nullptr, nullptr);
785 atomOffset += molt.atoms.nr;
791 static void spread_vsite3(const t_iatom ia[], real a, real b,
792 const rvec x[],
793 rvec f[], rvec fshift[],
794 const t_pbc *pbc, const t_graph *g)
796 rvec fi, fj, fk, dx;
797 int av, ai, aj, ak;
798 ivec di;
799 int siv, sij, sik;
801 av = ia[1];
802 ai = ia[2];
803 aj = ia[3];
804 ak = ia[4];
806 svmul(1 - a - b, f[av], fi);
807 svmul( a, f[av], fj);
808 svmul( b, f[av], fk);
809 /* 11 flops */
811 rvec_inc(f[ai], fi);
812 rvec_inc(f[aj], fj);
813 rvec_inc(f[ak], fk);
814 /* 9 Flops */
816 if (g)
818 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ia[1]), di);
819 siv = IVEC2IS(di);
820 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
821 sij = IVEC2IS(di);
822 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), di);
823 sik = IVEC2IS(di);
825 else if (pbc)
827 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
828 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
829 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
831 else
833 siv = CENTRAL;
834 sij = CENTRAL;
835 sik = CENTRAL;
838 if (fshift && (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL))
840 rvec_inc(fshift[siv], f[av]);
841 rvec_dec(fshift[CENTRAL], fi);
842 rvec_dec(fshift[sij], fj);
843 rvec_dec(fshift[sik], fk);
846 /* TOTAL: 20 flops */
849 static void spread_vsite3FD(const t_iatom ia[], real a, real b,
850 const rvec x[],
851 rvec f[], rvec fshift[],
852 gmx_bool VirCorr, matrix dxdf,
853 const t_pbc *pbc, const t_graph *g)
855 real c, invl, fproj, a1;
856 rvec xvi, xij, xjk, xix, fv, temp;
857 t_iatom av, ai, aj, ak;
858 int svi, sji, skj;
859 ivec di;
861 av = ia[1];
862 ai = ia[2];
863 aj = ia[3];
864 ak = ia[4];
865 copy_rvec(f[av], fv);
867 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
868 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
869 /* 6 flops */
871 /* xix goes from i to point x on the line jk */
872 xix[XX] = xij[XX]+a*xjk[XX];
873 xix[YY] = xij[YY]+a*xjk[YY];
874 xix[ZZ] = xij[ZZ]+a*xjk[ZZ];
875 /* 6 flops */
877 invl = gmx::invsqrt(iprod(xix, xix));
878 c = b*invl;
879 /* 4 + ?10? flops */
881 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
883 temp[XX] = c*(fv[XX]-fproj*xix[XX]);
884 temp[YY] = c*(fv[YY]-fproj*xix[YY]);
885 temp[ZZ] = c*(fv[ZZ]-fproj*xix[ZZ]);
886 /* 16 */
888 /* c is already calculated in constr_vsite3FD
889 storing c somewhere will save 26 flops! */
891 a1 = 1 - a;
892 f[ai][XX] += fv[XX] - temp[XX];
893 f[ai][YY] += fv[YY] - temp[YY];
894 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
895 f[aj][XX] += a1*temp[XX];
896 f[aj][YY] += a1*temp[YY];
897 f[aj][ZZ] += a1*temp[ZZ];
898 f[ak][XX] += a*temp[XX];
899 f[ak][YY] += a*temp[YY];
900 f[ak][ZZ] += a*temp[ZZ];
901 /* 19 Flops */
903 if (g)
905 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
906 svi = IVEC2IS(di);
907 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
908 sji = IVEC2IS(di);
909 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
910 skj = IVEC2IS(di);
912 else if (pbc)
914 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
916 else
918 svi = CENTRAL;
921 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
923 rvec_dec(fshift[svi], fv);
924 fshift[CENTRAL][XX] += fv[XX] - (1 + a)*temp[XX];
925 fshift[CENTRAL][YY] += fv[YY] - (1 + a)*temp[YY];
926 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a)*temp[ZZ];
927 fshift[ sji][XX] += temp[XX];
928 fshift[ sji][YY] += temp[YY];
929 fshift[ sji][ZZ] += temp[ZZ];
930 fshift[ skj][XX] += a*temp[XX];
931 fshift[ skj][YY] += a*temp[YY];
932 fshift[ skj][ZZ] += a*temp[ZZ];
935 if (VirCorr)
937 /* When VirCorr=TRUE, the virial for the current forces is not
938 * calculated from the redistributed forces. This means that
939 * the effect of non-linear virtual site constructions on the virial
940 * needs to be added separately. This contribution can be calculated
941 * in many ways, but the simplest and cheapest way is to use
942 * the first constructing atom ai as a reference position in space:
943 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
945 rvec xiv;
947 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
949 for (int i = 0; i < DIM; i++)
951 for (int j = 0; j < DIM; j++)
953 /* As xix is a linear combination of j and k, use that here */
954 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
959 /* TOTAL: 61 flops */
962 static void spread_vsite3FAD(const t_iatom ia[], real a, real b,
963 const rvec x[],
964 rvec f[], rvec fshift[],
965 gmx_bool VirCorr, matrix dxdf,
966 const t_pbc *pbc, const t_graph *g)
968 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
969 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
970 t_iatom av, ai, aj, ak;
971 int svi, sji, skj, d;
972 ivec di;
974 av = ia[1];
975 ai = ia[2];
976 aj = ia[3];
977 ak = ia[4];
978 copy_rvec(f[ia[1]], fv);
980 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
981 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
982 /* 6 flops */
984 invdij = gmx::invsqrt(iprod(xij, xij));
985 invdij2 = invdij * invdij;
986 c1 = iprod(xij, xjk) * invdij2;
987 xperp[XX] = xjk[XX] - c1*xij[XX];
988 xperp[YY] = xjk[YY] - c1*xij[YY];
989 xperp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
990 /* xperp in plane ijk, perp. to ij */
991 invdp = gmx::invsqrt(iprod(xperp, xperp));
992 a1 = a*invdij;
993 b1 = b*invdp;
994 /* 45 flops */
996 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
997 storing them somewhere will save 45 flops! */
999 fproj = iprod(xij, fv)*invdij2;
1000 svmul(fproj, xij, Fpij); /* proj. f on xij */
1001 svmul(iprod(xperp, fv)*invdp*invdp, xperp, Fppp); /* proj. f on xperp */
1002 svmul(b1*fproj, xperp, f3);
1003 /* 23 flops */
1005 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
1006 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
1007 for (d = 0; (d < DIM); d++)
1009 f1[d] *= a1;
1010 f2[d] *= b1;
1012 /* 12 flops */
1014 c2 = 1 + c1;
1015 f[ai][XX] += fv[XX] - f1[XX] + c1*f2[XX] + f3[XX];
1016 f[ai][YY] += fv[YY] - f1[YY] + c1*f2[YY] + f3[YY];
1017 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1*f2[ZZ] + f3[ZZ];
1018 f[aj][XX] += f1[XX] - c2*f2[XX] - f3[XX];
1019 f[aj][YY] += f1[YY] - c2*f2[YY] - f3[YY];
1020 f[aj][ZZ] += f1[ZZ] - c2*f2[ZZ] - f3[ZZ];
1021 f[ak][XX] += f2[XX];
1022 f[ak][YY] += f2[YY];
1023 f[ak][ZZ] += f2[ZZ];
1024 /* 30 Flops */
1026 if (g)
1028 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1029 svi = IVEC2IS(di);
1030 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1031 sji = IVEC2IS(di);
1032 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
1033 skj = IVEC2IS(di);
1035 else if (pbc)
1037 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1039 else
1041 svi = CENTRAL;
1044 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
1046 rvec_dec(fshift[svi], fv);
1047 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1-c1)*f2[XX] + f3[XX];
1048 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1-c1)*f2[YY] + f3[YY];
1049 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1-c1)*f2[ZZ] + f3[ZZ];
1050 fshift[ sji][XX] += f1[XX] - c1 *f2[XX] - f3[XX];
1051 fshift[ sji][YY] += f1[YY] - c1 *f2[YY] - f3[YY];
1052 fshift[ sji][ZZ] += f1[ZZ] - c1 *f2[ZZ] - f3[ZZ];
1053 fshift[ skj][XX] += f2[XX];
1054 fshift[ skj][YY] += f2[YY];
1055 fshift[ skj][ZZ] += f2[ZZ];
1058 if (VirCorr)
1060 rvec xiv;
1061 int i, j;
1063 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1065 for (i = 0; i < DIM; i++)
1067 for (j = 0; j < DIM; j++)
1069 /* Note that xik=xij+xjk, so we have to add xij*f2 */
1070 dxdf[i][j] +=
1071 -xiv[i]*fv[j]
1072 + xij[i]*(f1[j] + (1 - c2)*f2[j] - f3[j])
1073 + xjk[i]*f2[j];
1078 /* TOTAL: 113 flops */
1081 static void spread_vsite3OUT(const t_iatom ia[], real a, real b, real c,
1082 const rvec x[],
1083 rvec f[], rvec fshift[],
1084 gmx_bool VirCorr, matrix dxdf,
1085 const t_pbc *pbc, const t_graph *g)
1087 rvec xvi, xij, xik, fv, fj, fk;
1088 real cfx, cfy, cfz;
1089 int av, ai, aj, ak;
1090 ivec di;
1091 int svi, sji, ski;
1093 av = ia[1];
1094 ai = ia[2];
1095 aj = ia[3];
1096 ak = ia[4];
1098 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1099 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1100 /* 6 Flops */
1102 copy_rvec(f[av], fv);
1104 cfx = c*fv[XX];
1105 cfy = c*fv[YY];
1106 cfz = c*fv[ZZ];
1107 /* 3 Flops */
1109 fj[XX] = a*fv[XX] - xik[ZZ]*cfy + xik[YY]*cfz;
1110 fj[YY] = xik[ZZ]*cfx + a*fv[YY] - xik[XX]*cfz;
1111 fj[ZZ] = -xik[YY]*cfx + xik[XX]*cfy + a*fv[ZZ];
1113 fk[XX] = b*fv[XX] + xij[ZZ]*cfy - xij[YY]*cfz;
1114 fk[YY] = -xij[ZZ]*cfx + b*fv[YY] + xij[XX]*cfz;
1115 fk[ZZ] = xij[YY]*cfx - xij[XX]*cfy + b*fv[ZZ];
1116 /* 30 Flops */
1118 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
1119 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
1120 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1121 rvec_inc(f[aj], fj);
1122 rvec_inc(f[ak], fk);
1123 /* 15 Flops */
1125 if (g)
1127 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1128 svi = IVEC2IS(di);
1129 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1130 sji = IVEC2IS(di);
1131 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1132 ski = IVEC2IS(di);
1134 else if (pbc)
1136 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1138 else
1140 svi = CENTRAL;
1143 if (fshift && (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL))
1145 rvec_dec(fshift[svi], fv);
1146 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
1147 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
1148 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1149 rvec_inc(fshift[sji], fj);
1150 rvec_inc(fshift[ski], fk);
1153 if (VirCorr)
1155 rvec xiv;
1157 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1159 for (int i = 0; i < DIM; i++)
1161 for (int j = 0; j < DIM; j++)
1163 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j];
1168 /* TOTAL: 54 flops */
1171 static void spread_vsite4FD(const t_iatom ia[], real a, real b, real c,
1172 const rvec x[],
1173 rvec f[], rvec fshift[],
1174 gmx_bool VirCorr, matrix dxdf,
1175 const t_pbc *pbc, const t_graph *g)
1177 real d, invl, fproj, a1;
1178 rvec xvi, xij, xjk, xjl, xix, fv, temp;
1179 int av, ai, aj, ak, al;
1180 ivec di;
1181 int svi, sji, skj, slj, m;
1183 av = ia[1];
1184 ai = ia[2];
1185 aj = ia[3];
1186 ak = ia[4];
1187 al = ia[5];
1189 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1190 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1191 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1192 /* 9 flops */
1194 /* xix goes from i to point x on the plane jkl */
1195 for (m = 0; m < DIM; m++)
1197 xix[m] = xij[m] + a*xjk[m] + b*xjl[m];
1199 /* 12 flops */
1201 invl = gmx::invsqrt(iprod(xix, xix));
1202 d = c*invl;
1203 /* 4 + ?10? flops */
1205 copy_rvec(f[av], fv);
1207 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
1209 for (m = 0; m < DIM; m++)
1211 temp[m] = d*(fv[m] - fproj*xix[m]);
1213 /* 16 */
1215 /* c is already calculated in constr_vsite3FD
1216 storing c somewhere will save 35 flops! */
1218 a1 = 1 - a - b;
1219 for (m = 0; m < DIM; m++)
1221 f[ai][m] += fv[m] - temp[m];
1222 f[aj][m] += a1*temp[m];
1223 f[ak][m] += a*temp[m];
1224 f[al][m] += b*temp[m];
1226 /* 26 Flops */
1228 if (g)
1230 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1231 svi = IVEC2IS(di);
1232 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1233 sji = IVEC2IS(di);
1234 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
1235 skj = IVEC2IS(di);
1236 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, aj), di);
1237 slj = IVEC2IS(di);
1239 else if (pbc)
1241 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1243 else
1245 svi = CENTRAL;
1248 if (fshift &&
1249 (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL))
1251 rvec_dec(fshift[svi], fv);
1252 for (m = 0; m < DIM; m++)
1254 fshift[CENTRAL][m] += fv[m] - (1 + a + b)*temp[m];
1255 fshift[ sji][m] += temp[m];
1256 fshift[ skj][m] += a*temp[m];
1257 fshift[ slj][m] += b*temp[m];
1261 if (VirCorr)
1263 rvec xiv;
1264 int i, j;
1266 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1268 for (i = 0; i < DIM; i++)
1270 for (j = 0; j < DIM; j++)
1272 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
1277 /* TOTAL: 77 flops */
1281 static void spread_vsite4FDN(const t_iatom ia[], real a, real b, real c,
1282 const rvec x[],
1283 rvec f[], rvec fshift[],
1284 gmx_bool VirCorr, matrix dxdf,
1285 const t_pbc *pbc, const t_graph *g)
1287 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1288 rvec fv, fj, fk, fl;
1289 real invrm, denom;
1290 real cfx, cfy, cfz;
1291 ivec di;
1292 int av, ai, aj, ak, al;
1293 int svi, sij, sik, sil;
1295 /* DEBUG: check atom indices */
1296 av = ia[1];
1297 ai = ia[2];
1298 aj = ia[3];
1299 ak = ia[4];
1300 al = ia[5];
1302 copy_rvec(f[av], fv);
1304 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1305 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1306 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1307 /* 9 flops */
1309 ra[XX] = a*xik[XX];
1310 ra[YY] = a*xik[YY];
1311 ra[ZZ] = a*xik[ZZ];
1313 rb[XX] = b*xil[XX];
1314 rb[YY] = b*xil[YY];
1315 rb[ZZ] = b*xil[ZZ];
1317 /* 6 flops */
1319 rvec_sub(ra, xij, rja);
1320 rvec_sub(rb, xij, rjb);
1321 rvec_sub(rb, ra, rab);
1322 /* 9 flops */
1324 cprod(rja, rjb, rm);
1325 /* 9 flops */
1327 invrm = gmx::invsqrt(norm2(rm));
1328 denom = invrm*invrm;
1329 /* 5+5+2 flops */
1331 cfx = c*invrm*fv[XX];
1332 cfy = c*invrm*fv[YY];
1333 cfz = c*invrm*fv[ZZ];
1334 /* 6 Flops */
1336 cprod(rm, rab, rt);
1337 /* 9 flops */
1339 rt[XX] *= denom;
1340 rt[YY] *= denom;
1341 rt[ZZ] *= denom;
1342 /* 3flops */
1344 fj[XX] = ( -rm[XX]*rt[XX]) * cfx + ( rab[ZZ]-rm[YY]*rt[XX]) * cfy + (-rab[YY]-rm[ZZ]*rt[XX]) * cfz;
1345 fj[YY] = (-rab[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( rab[XX]-rm[ZZ]*rt[YY]) * cfz;
1346 fj[ZZ] = ( rab[YY]-rm[XX]*rt[ZZ]) * cfx + (-rab[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1347 /* 30 flops */
1349 cprod(rjb, rm, rt);
1350 /* 9 flops */
1352 rt[XX] *= denom*a;
1353 rt[YY] *= denom*a;
1354 rt[ZZ] *= denom*a;
1355 /* 3flops */
1357 fk[XX] = ( -rm[XX]*rt[XX]) * cfx + (-a*rjb[ZZ]-rm[YY]*rt[XX]) * cfy + ( a*rjb[YY]-rm[ZZ]*rt[XX]) * cfz;
1358 fk[YY] = ( a*rjb[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + (-a*rjb[XX]-rm[ZZ]*rt[YY]) * cfz;
1359 fk[ZZ] = (-a*rjb[YY]-rm[XX]*rt[ZZ]) * cfx + ( a*rjb[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1360 /* 36 flops */
1362 cprod(rm, rja, rt);
1363 /* 9 flops */
1365 rt[XX] *= denom*b;
1366 rt[YY] *= denom*b;
1367 rt[ZZ] *= denom*b;
1368 /* 3flops */
1370 fl[XX] = ( -rm[XX]*rt[XX]) * cfx + ( b*rja[ZZ]-rm[YY]*rt[XX]) * cfy + (-b*rja[YY]-rm[ZZ]*rt[XX]) * cfz;
1371 fl[YY] = (-b*rja[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( b*rja[XX]-rm[ZZ]*rt[YY]) * cfz;
1372 fl[ZZ] = ( b*rja[YY]-rm[XX]*rt[ZZ]) * cfx + (-b*rja[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1373 /* 36 flops */
1375 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1376 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1377 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1378 rvec_inc(f[aj], fj);
1379 rvec_inc(f[ak], fk);
1380 rvec_inc(f[al], fl);
1381 /* 21 flops */
1383 if (g)
1385 ivec_sub(SHIFT_IVEC(g, av), SHIFT_IVEC(g, ai), di);
1386 svi = IVEC2IS(di);
1387 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1388 sij = IVEC2IS(di);
1389 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1390 sik = IVEC2IS(di);
1391 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, ai), di);
1392 sil = IVEC2IS(di);
1394 else if (pbc)
1396 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1398 else
1400 svi = CENTRAL;
1403 if (fshift && (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL))
1405 rvec_dec(fshift[svi], fv);
1406 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1407 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1408 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1409 rvec_inc(fshift[sij], fj);
1410 rvec_inc(fshift[sik], fk);
1411 rvec_inc(fshift[sil], fl);
1414 if (VirCorr)
1416 rvec xiv;
1417 int i, j;
1419 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1421 for (i = 0; i < DIM; i++)
1423 for (j = 0; j < DIM; j++)
1425 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j] + xil[i]*fl[j];
1430 /* Total: 207 flops (Yuck!) */
1434 static int spread_vsiten(const t_iatom ia[], const t_iparams ip[],
1435 const rvec x[],
1436 rvec f[], rvec fshift[],
1437 const t_pbc *pbc, const t_graph *g)
1439 rvec xv, dx, fi;
1440 int n3, av, i, ai;
1441 real a;
1442 ivec di;
1443 int siv;
1445 n3 = 3*ip[ia[0]].vsiten.n;
1446 av = ia[1];
1447 copy_rvec(x[av], xv);
1449 for (i = 0; i < n3; i += 3)
1451 ai = ia[i+2];
1452 if (g)
1454 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
1455 siv = IVEC2IS(di);
1457 else if (pbc)
1459 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1461 else
1463 siv = CENTRAL;
1465 a = ip[ia[i]].vsiten.a;
1466 svmul(a, f[av], fi);
1467 rvec_inc(f[ai], fi);
1468 if (fshift && siv != CENTRAL)
1470 rvec_inc(fshift[siv], fi);
1471 rvec_dec(fshift[CENTRAL], fi);
1473 /* 6 Flops */
1476 return n3;
1480 static int vsite_count(const t_ilist *ilist, int ftype)
1482 if (ftype == F_VSITEN)
1484 return ilist[ftype].nr/3;
1486 else
1488 return ilist[ftype].nr/(1 + interaction_function[ftype].nratoms);
1492 static void spread_vsite_f_thread(const gmx_vsite_t *vsite,
1493 const rvec x[],
1494 rvec f[], rvec *fshift,
1495 gmx_bool VirCorr, matrix dxdf,
1496 t_iparams ip[], const t_ilist ilist[],
1497 const t_graph *g, const t_pbc *pbc_null)
1499 const PbcMode pbcMode = getPbcMode(pbc_null, vsite);
1500 /* We need another pbc pointer, as with charge groups we switch per vsite */
1501 const t_pbc *pbc_null2 = pbc_null;
1502 const int *vsite_pbc = nullptr;
1504 /* this loop goes backwards to be able to build *
1505 * higher type vsites from lower types */
1506 for (int ftype = c_ftypeVsiteEnd - 1; ftype >= c_ftypeVsiteStart; ftype--)
1508 if (ilist[ftype].nr == 0)
1510 continue;
1513 { // TODO remove me
1514 int nra = interaction_function[ftype].nratoms;
1515 int inc = 1 + nra;
1516 int nr = ilist[ftype].nr;
1518 const t_iatom *ia = ilist[ftype].iatoms;
1520 if (pbcMode == PbcMode::all)
1522 pbc_null2 = pbc_null;
1524 else if (pbcMode == PbcMode::chargeGroup)
1526 vsite_pbc = vsite->vsite_pbc_loc[ftype - c_ftypeVsiteStart];
1529 for (int i = 0; i < nr; )
1531 if (vsite_pbc != nullptr)
1533 if (vsite_pbc[i/(1 + nra)] > -2)
1535 pbc_null2 = pbc_null;
1537 else
1539 pbc_null2 = nullptr;
1543 int tp = ia[0];
1545 /* Constants for constructing */
1546 real a1, b1, c1;
1547 a1 = ip[tp].vsite.a;
1548 /* Construct the vsite depending on type */
1549 switch (ftype)
1551 case F_VSITE2:
1552 spread_vsite2(ia, a1, x, f, fshift, pbc_null2, g);
1553 break;
1554 case F_VSITE3:
1555 b1 = ip[tp].vsite.b;
1556 spread_vsite3(ia, a1, b1, x, f, fshift, pbc_null2, g);
1557 break;
1558 case F_VSITE3FD:
1559 b1 = ip[tp].vsite.b;
1560 spread_vsite3FD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1561 break;
1562 case F_VSITE3FAD:
1563 b1 = ip[tp].vsite.b;
1564 spread_vsite3FAD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1565 break;
1566 case F_VSITE3OUT:
1567 b1 = ip[tp].vsite.b;
1568 c1 = ip[tp].vsite.c;
1569 spread_vsite3OUT(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1570 break;
1571 case F_VSITE4FD:
1572 b1 = ip[tp].vsite.b;
1573 c1 = ip[tp].vsite.c;
1574 spread_vsite4FD(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1575 break;
1576 case F_VSITE4FDN:
1577 b1 = ip[tp].vsite.b;
1578 c1 = ip[tp].vsite.c;
1579 spread_vsite4FDN(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1580 break;
1581 case F_VSITEN:
1582 inc = spread_vsiten(ia, ip, x, f, fshift, pbc_null2, g);
1583 break;
1584 default:
1585 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
1586 ftype, __FILE__, __LINE__);
1588 clear_rvec(f[ia[1]]);
1590 /* Increment loop variables */
1591 i += inc;
1592 ia += inc;
1598 /*! \brief Clears the task force buffer elements that are written by task idTask */
1599 static void clearTaskForceBufferUsedElements(InterdependentTask *idTask)
1601 int ntask = idTask->spreadTask.size();
1602 for (int ti = 0; ti < ntask; ti++)
1604 const AtomIndex *atomList = &idTask->atomIndex[idTask->spreadTask[ti]];
1605 int natom = atomList->atom.size();
1606 RVec *force = idTask->force.data();
1607 for (int i = 0; i < natom; i++)
1609 clear_rvec(force[atomList->atom[i]]);
1614 void spread_vsite_f(const gmx_vsite_t *vsite,
1615 const rvec * gmx_restrict x,
1616 rvec * gmx_restrict f, rvec * gmx_restrict fshift,
1617 gmx_bool VirCorr, matrix vir,
1618 t_nrnb *nrnb, const t_idef *idef,
1619 int ePBC, gmx_bool bMolPBC, const t_graph *g, const matrix box,
1620 const t_commrec *cr, gmx_wallcycle *wcycle)
1622 wallcycle_start(wcycle, ewcVSITESPREAD);
1623 const bool useDomdec = vsite->useDomdec;
1624 GMX_ASSERT(!useDomdec || (cr != nullptr && DOMAINDECOMP(cr)), "When vsites are set up with domain decomposition, we need a valid commrec");
1626 t_pbc pbc, *pbc_null;
1628 /* We only need to do pbc when we have inter-cg vsites */
1629 if ((useDomdec || bMolPBC) && vsite->n_intercg_vsite)
1631 /* This is wasting some CPU time as we now do this multiple times
1632 * per MD step.
1634 pbc_null = set_pbc_dd(&pbc, ePBC, useDomdec ? cr->dd->nc : nullptr, FALSE, box);
1636 else
1638 pbc_null = nullptr;
1641 if (useDomdec)
1643 dd_clear_f_vsites(cr->dd, f);
1646 if (vsite->nthreads == 1)
1648 matrix dxdf;
1649 if (VirCorr)
1651 clear_mat(dxdf);
1653 spread_vsite_f_thread(vsite,
1654 x, f, fshift,
1655 VirCorr, dxdf,
1656 idef->iparams, idef->il,
1657 g, pbc_null);
1659 if (VirCorr)
1661 for (int i = 0; i < DIM; i++)
1663 for (int j = 0; j < DIM; j++)
1665 vir[i][j] += -0.5*dxdf[i][j];
1670 else
1672 /* First spread the vsites that might depend on non-local vsites */
1673 if (VirCorr)
1675 clear_mat(vsite->tData[vsite->nthreads]->dxdf);
1677 spread_vsite_f_thread(vsite,
1678 x, f, fshift,
1679 VirCorr, vsite->tData[vsite->nthreads]->dxdf,
1680 idef->iparams,
1681 vsite->tData[vsite->nthreads]->ilist,
1682 g, pbc_null);
1684 #pragma omp parallel num_threads(vsite->nthreads)
1688 int thread = gmx_omp_get_thread_num();
1689 VsiteThread *tData = vsite->tData[thread];
1691 rvec *fshift_t;
1692 if (thread == 0 || fshift == nullptr)
1694 fshift_t = fshift;
1696 else
1698 fshift_t = tData->fshift;
1700 for (int i = 0; i < SHIFTS; i++)
1702 clear_rvec(fshift_t[i]);
1705 if (VirCorr)
1707 clear_mat(tData->dxdf);
1710 if (tData->useInterdependentTask)
1712 /* Spread the vsites that spread outside our local range.
1713 * This is done using a thread-local force buffer force.
1714 * First we need to copy the input vsite forces to force.
1716 InterdependentTask *idTask = &tData->idTask;
1718 /* Clear the buffer elements set by our task during
1719 * the last call to spread_vsite_f.
1721 clearTaskForceBufferUsedElements(idTask);
1723 int nvsite = idTask->vsite.size();
1724 for (int i = 0; i < nvsite; i++)
1726 copy_rvec(f[idTask->vsite[i]],
1727 idTask->force[idTask->vsite[i]]);
1729 spread_vsite_f_thread(vsite,
1730 x, as_rvec_array(idTask->force.data()), fshift_t,
1731 VirCorr, tData->dxdf,
1732 idef->iparams,
1733 tData->idTask.ilist,
1734 g, pbc_null);
1736 /* We need a barrier before reducing forces below
1737 * that have been produced by a different thread above.
1739 #pragma omp barrier
1741 /* Loop over all thread task and reduce forces they
1742 * produced on atoms that fall in our range.
1743 * Note that atomic reduction would be a simpler solution,
1744 * but that might not have good support on all platforms.
1746 int ntask = idTask->reduceTask.size();
1747 for (int ti = 0; ti < ntask; ti++)
1749 const InterdependentTask *idt_foreign = &vsite->tData[idTask->reduceTask[ti]]->idTask;
1750 const AtomIndex *atomList = &idt_foreign->atomIndex[thread];
1751 const RVec *f_foreign = idt_foreign->force.data();
1753 int natom = atomList->atom.size();
1754 for (int i = 0; i < natom; i++)
1756 int ind = atomList->atom[i];
1757 rvec_inc(f[ind], f_foreign[ind]);
1758 /* Clearing of f_foreign is done at the next step */
1761 /* Clear the vsite forces, both in f and force */
1762 for (int i = 0; i < nvsite; i++)
1764 int ind = tData->idTask.vsite[i];
1765 clear_rvec(f[ind]);
1766 clear_rvec(tData->idTask.force[ind]);
1770 /* Spread the vsites that spread locally only */
1771 spread_vsite_f_thread(vsite,
1772 x, f, fshift_t,
1773 VirCorr, tData->dxdf,
1774 idef->iparams,
1775 tData->ilist,
1776 g, pbc_null);
1778 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1781 if (fshift != nullptr)
1783 for (int th = 1; th < vsite->nthreads; th++)
1785 for (int i = 0; i < SHIFTS; i++)
1787 rvec_inc(fshift[i], vsite->tData[th]->fshift[i]);
1792 if (VirCorr)
1794 for (int th = 0; th < vsite->nthreads + 1; th++)
1796 /* MSVC doesn't like matrix references, so we use a pointer */
1797 const matrix *dxdf = &vsite->tData[th]->dxdf;
1799 for (int i = 0; i < DIM; i++)
1801 for (int j = 0; j < DIM; j++)
1803 vir[i][j] += -0.5*(*dxdf)[i][j];
1810 if (useDomdec)
1812 dd_move_f_vsites(cr->dd, f, fshift);
1815 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(idef->il, F_VSITE2));
1816 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(idef->il, F_VSITE3));
1817 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(idef->il, F_VSITE3FD));
1818 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(idef->il, F_VSITE3FAD));
1819 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(idef->il, F_VSITE3OUT));
1820 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(idef->il, F_VSITE4FD));
1821 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(idef->il, F_VSITE4FDN));
1822 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(idef->il, F_VSITEN));
1824 wallcycle_stop(wcycle, ewcVSITESPREAD);
1827 /*! \brief Returns the an array with charge-group indices for each atom
1829 * \param[in] chargeGroups The charge group block struct
1831 static std::vector<int> atom2cg(const t_block &chargeGroups)
1833 std::vector<int> a2cg(chargeGroups.index[chargeGroups.nr], 0);
1835 for (int chargeGroup = 0; chargeGroup < chargeGroups.nr; chargeGroup++)
1837 std::fill(a2cg.begin() + chargeGroups.index[chargeGroup],
1838 a2cg.begin() + chargeGroups.index[chargeGroup + 1],
1839 chargeGroup);
1842 return a2cg;
1845 int count_intercg_vsites(const gmx_mtop_t *mtop)
1847 int n_intercg_vsite = 0;
1848 for (const gmx_molblock_t &molb : mtop->molblock)
1850 const gmx_moltype_t &molt = mtop->moltype[molb.type];
1852 std::vector<int> a2cg = atom2cg(molt.cgs);
1853 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
1855 int nral = NRAL(ftype);
1856 const t_ilist &il = molt.ilist[ftype];
1857 const t_iatom *ia = il.iatoms;
1858 for (int i = 0; i < il.nr; i += 1 + nral)
1860 int cg = a2cg[ia[1+i]];
1861 for (int a = 1; a < nral; a++)
1863 if (a2cg[ia[1+a]] != cg)
1865 n_intercg_vsite += molb.nmol;
1866 break;
1873 return n_intercg_vsite;
1876 static int **get_vsite_pbc(const t_iparams *iparams, const t_ilist *ilist,
1877 const t_atom *atom, const t_mdatoms *md,
1878 const t_block &cgs)
1880 /* Make an atom to charge group index */
1881 std::vector<int> a2cg = atom2cg(cgs);
1883 /* Make an array that tells if the pbc of an atom is set */
1884 std::vector<bool> pbc_set(cgs.index[cgs.nr], false);
1885 /* PBC is set for all non vsites */
1886 for (int a = 0; a < cgs.index[cgs.nr]; a++)
1888 if ((atom && atom[a].ptype != eptVSite) ||
1889 (md && md->ptype[a] != eptVSite))
1891 pbc_set[a] = true;
1895 int **vsite_pbc;
1896 snew(vsite_pbc, F_VSITEN-F_VSITE2+1);
1898 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
1900 { // TODO remove me
1901 int nral = NRAL(ftype);
1902 const t_ilist *il = &ilist[ftype];
1903 const t_iatom *ia = il->iatoms;
1904 int *vsite_pbc_f;
1906 snew(vsite_pbc[ftype-F_VSITE2], il->nr/(1 + nral));
1907 vsite_pbc_f = vsite_pbc[ftype-F_VSITE2];
1909 int i = 0;
1910 while (i < il->nr)
1912 int vsi = i/(1 + nral);
1913 t_iatom vsite = ia[i+1];
1914 int cg_v = a2cg[vsite];
1915 /* A value of -2 signals that this vsite and its contructing
1916 * atoms are all within the same cg, so no pbc is required.
1918 vsite_pbc_f[vsi] = -2;
1919 /* Check if constructing atoms are outside the vsite's cg */
1920 int nc3 = 0;
1921 if (ftype == F_VSITEN)
1923 nc3 = 3*iparams[ia[i]].vsiten.n;
1924 for (int j = 0; j < nc3; j += 3)
1926 if (a2cg[ia[i+j+2]] != cg_v)
1928 vsite_pbc_f[vsi] = -1;
1932 else
1934 for (int a = 1; a < nral; a++)
1936 if (a2cg[ia[i+1+a]] != cg_v)
1938 vsite_pbc_f[vsi] = -1;
1942 if (vsite_pbc_f[vsi] == -1)
1944 /* Check if this is the first processed atom of a vsite only cg */
1945 gmx_bool bViteOnlyCG_and_FirstAtom = TRUE;
1946 for (int a = cgs.index[cg_v]; a < cgs.index[cg_v + 1]; a++)
1948 /* Non-vsites already have pbc set, so simply check for pbc_set */
1949 if (pbc_set[a])
1951 bViteOnlyCG_and_FirstAtom = FALSE;
1952 break;
1955 if (bViteOnlyCG_and_FirstAtom)
1957 /* First processed atom of a vsite only charge group.
1958 * The pbc of the input coordinates to construct_vsites
1959 * should be preserved.
1961 vsite_pbc_f[vsi] = vsite;
1963 else if (cg_v != a2cg[ia[1+i+1]])
1965 /* This vsite has a different charge group index
1966 * than it's first constructing atom
1967 * and the charge group has more than one atom,
1968 * search for the first normal particle
1969 * or vsite that already had its pbc defined.
1970 * If nothing is found, use full pbc for this vsite.
1972 for (int a = cgs.index[cg_v]; a < cgs.index[cg_v + 1]; a++)
1974 if (a != vsite && pbc_set[a])
1976 vsite_pbc_f[vsi] = a;
1977 if (gmx_debug_at)
1979 fprintf(debug, "vsite %d match pbc with atom %d\n",
1980 vsite+1, a+1);
1982 break;
1985 if (gmx_debug_at)
1987 fprintf(debug, "vsite atom %d cg %d - %d pbc atom %d\n",
1988 vsite+1, cgs.index[cg_v] + 1, cgs.index[cg_v + 1],
1989 vsite_pbc_f[vsi] + 1);
1993 if (ftype == F_VSITEN)
1995 /* The other entries in vsite_pbc_f are not used for center vsites */
1996 i += nc3;
1998 else
2000 i += 1 + nral;
2003 /* This vsite now has its pbc defined */
2004 pbc_set[vsite] = true;
2009 return vsite_pbc;
2013 gmx_vsite_t *initVsite(const gmx_mtop_t &mtop,
2014 const t_commrec *cr)
2016 GMX_RELEASE_ASSERT(cr != nullptr, "We need a valid commrec");
2018 /* check if there are vsites */
2019 int nvsite = 0;
2020 for (int ftype = 0; ftype < F_NRE; ftype++)
2022 if (interaction_function[ftype].flags & IF_VSITE)
2024 GMX_ASSERT(ftype >= c_ftypeVsiteStart && ftype < c_ftypeVsiteEnd, "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
2026 nvsite += gmx_mtop_ftype_count(&mtop, ftype);
2028 else
2030 GMX_ASSERT(ftype < c_ftypeVsiteStart || ftype >= c_ftypeVsiteEnd, "c_ftypeVsiteStart and/or c_ftypeVsiteEnd do not have correct values");
2034 if (nvsite == 0)
2036 return nullptr;
2039 gmx_vsite_t *vsite = new(gmx_vsite_t);
2041 vsite->n_intercg_vsite = count_intercg_vsites(&mtop);
2043 vsite->bHaveChargeGroups = (ncg_mtop(&mtop) < mtop.natoms);
2045 vsite->useDomdec = (DOMAINDECOMP(cr) && cr->dd->nnodes > 1);
2047 /* If we don't have charge groups, the vsite follows its own pbc.
2049 * With charge groups, each vsite needs to follow the pbc of the charge
2050 * group. Thus we need to keep track of PBC. Here we assume that without
2051 * domain decomposition all molecules are whole (which will not be
2052 * the case with periodic molecules).
2054 if (vsite->bHaveChargeGroups &&
2055 vsite->n_intercg_vsite > 0 &&
2056 DOMAINDECOMP(cr))
2058 vsite->nvsite_pbc_molt = mtop.moltype.size();
2059 snew(vsite->vsite_pbc_molt, vsite->nvsite_pbc_molt);
2060 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
2062 const gmx_moltype_t &molt = mtop.moltype[mt];
2063 vsite->vsite_pbc_molt[mt] = get_vsite_pbc(mtop.ffparams.iparams,
2064 molt.ilist,
2065 molt.atoms.atom, nullptr,
2066 molt.cgs);
2069 snew(vsite->vsite_pbc_loc_nalloc, c_ftypeVsiteEnd - c_ftypeVsiteStart);
2070 snew(vsite->vsite_pbc_loc, c_ftypeVsiteEnd - c_ftypeVsiteStart);
2072 else
2074 vsite->vsite_pbc_molt = nullptr;
2075 vsite->vsite_pbc_loc = nullptr;
2078 vsite->nthreads = gmx_omp_nthreads_get(emntVSITE);
2080 if (vsite->nthreads > 1)
2082 /* We need one extra thread data structure for the overlap vsites */
2083 snew(vsite->tData, vsite->nthreads + 1);
2084 #pragma omp parallel for num_threads(vsite->nthreads) schedule(static)
2085 for (int thread = 0; thread < vsite->nthreads; thread++)
2089 vsite->tData[thread] = new VsiteThread;
2091 InterdependentTask *idTask = &vsite->tData[thread]->idTask;
2092 idTask->nuse = 0;
2093 idTask->atomIndex.resize(vsite->nthreads);
2095 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2097 if (vsite->nthreads > 1)
2099 vsite->tData[vsite->nthreads] = new VsiteThread;
2103 vsite->taskIndex = nullptr;
2104 vsite->taskIndexNalloc = 0;
2106 return vsite;
2109 static inline void flagAtom(InterdependentTask *idTask, int atom,
2110 int thread, int nthread, int natperthread)
2112 if (!idTask->use[atom])
2114 idTask->use[atom] = true;
2115 thread = atom/natperthread;
2116 /* Assign all non-local atom force writes to thread 0 */
2117 if (thread >= nthread)
2119 thread = 0;
2121 idTask->atomIndex[thread].atom.push_back(atom);
2125 /*\brief Here we try to assign all vsites that are in our local range.
2127 * Our task local atom range is tData->rangeStart - tData->rangeEnd.
2128 * Vsites that depend only on local atoms, as indicated by taskIndex[]==thread,
2129 * are assigned to task tData->ilist. Vsites that depend on non-local atoms
2130 * but not on other vsites are assigned to task tData->id_task.ilist.
2131 * taskIndex[] is set for all vsites in our range, either to our local tasks
2132 * or to the single last task as taskIndex[]=2*nthreads.
2134 static void assignVsitesToThread(VsiteThread *tData,
2135 int thread,
2136 int nthread,
2137 int natperthread,
2138 int *taskIndex,
2139 const t_ilist *ilist,
2140 const t_iparams *ip,
2141 const unsigned short *ptype)
2143 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2145 tData->ilist[ftype].nr = 0;
2146 tData->idTask.ilist[ftype].nr = 0;
2148 int nral1 = 1 + NRAL(ftype);
2149 int inc = nral1;
2150 t_iatom *iat = ilist[ftype].iatoms;
2151 for (int i = 0; i < ilist[ftype].nr; )
2153 if (ftype == F_VSITEN)
2155 /* The 3 below is from 1+NRAL(ftype)=3 */
2156 inc = ip[iat[i]].vsiten.n*3;
2159 if (iat[1 + i] < tData->rangeStart ||
2160 iat[1 + i] >= tData->rangeEnd)
2162 /* This vsite belongs to a different thread */
2163 i += inc;
2164 continue;
2167 /* We would like to assign this vsite to task thread,
2168 * but it might depend on atoms outside the atom range of thread
2169 * or on another vsite not assigned to task thread.
2171 int task = thread;
2172 if (ftype != F_VSITEN)
2174 for (int j = i + 2; j < i + nral1; j++)
2176 /* Do a range check to avoid a harmless race on taskIndex */
2177 if (iat[j] < tData->rangeStart ||
2178 iat[j] >= tData->rangeEnd ||
2179 taskIndex[iat[j]] != thread)
2181 if (!tData->useInterdependentTask ||
2182 ptype[iat[j]] == eptVSite)
2184 /* At least one constructing atom is a vsite
2185 * that is not assigned to the same thread.
2186 * Put this vsite into a separate task.
2188 task = 2*nthread;
2189 break;
2192 /* There are constructing atoms outside our range,
2193 * put this vsite into a second task to be executed
2194 * on the same thread. During construction no barrier
2195 * is needed between the two tasks on the same thread.
2196 * During spreading we need to run this task with
2197 * an additional thread-local intermediate force buffer
2198 * (or atomic reduction) and a barrier between the two
2199 * tasks.
2201 task = nthread + thread;
2205 else
2207 for (int j = i + 2; j < i + inc; j += 3)
2209 /* Do a range check to avoid a harmless race on taskIndex */
2210 if (iat[j] < tData->rangeStart ||
2211 iat[j] >= tData->rangeEnd ||
2212 taskIndex[iat[j]] != thread)
2214 GMX_ASSERT(ptype[iat[j]] != eptVSite, "A vsite to be assigned in assignVsitesToThread has a vsite as a constructing atom that does not belong to our task, such vsites should be assigned to the single 'master' task");
2216 task = nthread + thread;
2221 /* Update this vsite's thread index entry */
2222 taskIndex[iat[1+i]] = task;
2224 if (task == thread || task == nthread + thread)
2226 /* Copy this vsite to the thread data struct of thread */
2227 t_ilist *il_task;
2228 if (task == thread)
2230 il_task = &tData->ilist[ftype];
2232 else
2234 il_task = &tData->idTask.ilist[ftype];
2236 /* Ensure we have sufficient memory allocated */
2237 if (il_task->nr + inc > il_task->nalloc)
2239 il_task->nalloc = over_alloc_large(il_task->nr + inc);
2240 srenew(il_task->iatoms, il_task->nalloc);
2242 /* Copy the vsite data to the thread-task local array */
2243 for (int j = i; j < i + inc; j++)
2245 il_task->iatoms[il_task->nr++] = iat[j];
2247 if (task == nthread + thread)
2249 /* This vsite write outside our own task force block.
2250 * Put it into the interdependent task list and flag
2251 * the atoms involved for reduction.
2253 tData->idTask.vsite.push_back(iat[i + 1]);
2254 if (ftype != F_VSITEN)
2256 for (int j = i + 2; j < i + nral1; j++)
2258 flagAtom(&tData->idTask, iat[j],
2259 thread, nthread, natperthread);
2262 else
2264 for (int j = i + 2; j < i + inc; j += 3)
2266 flagAtom(&tData->idTask, iat[j],
2267 thread, nthread, natperthread);
2273 i += inc;
2278 /*! \brief Assign all vsites with taskIndex[]==task to task tData */
2279 static void assignVsitesToSingleTask(VsiteThread *tData,
2280 int task,
2281 const int *taskIndex,
2282 const t_ilist *ilist,
2283 const t_iparams *ip)
2285 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2287 tData->ilist[ftype].nr = 0;
2288 tData->idTask.ilist[ftype].nr = 0;
2290 int nral1 = 1 + NRAL(ftype);
2291 int inc = nral1;
2292 t_iatom *iat = ilist[ftype].iatoms;
2293 t_ilist *il_task = &tData->ilist[ftype];
2295 for (int i = 0; i < ilist[ftype].nr; )
2297 if (ftype == F_VSITEN)
2299 /* The 3 below is from 1+NRAL(ftype)=3 */
2300 inc = ip[iat[i]].vsiten.n*3;
2302 /* Check if the vsite is assigned to our task */
2303 if (taskIndex[iat[1 + i]] == task)
2305 /* Ensure we have sufficient memory allocated */
2306 if (il_task->nr + inc > il_task->nalloc)
2308 il_task->nalloc = over_alloc_large(il_task->nr + inc);
2309 srenew(il_task->iatoms, il_task->nalloc);
2311 /* Copy the vsite data to the thread-task local array */
2312 for (int j = i; j < i + inc; j++)
2314 il_task->iatoms[il_task->nr++] = iat[j];
2318 i += inc;
2323 void split_vsites_over_threads(const t_ilist *ilist,
2324 const t_iparams *ip,
2325 const t_mdatoms *mdatoms,
2326 gmx_vsite_t *vsite)
2328 int vsite_atom_range, natperthread;
2330 if (vsite->nthreads == 1)
2332 /* Nothing to do */
2333 return;
2336 /* The current way of distributing the vsites over threads in primitive.
2337 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
2338 * without taking into account how the vsites are distributed.
2339 * Without domain decomposition we at least tighten the upper bound
2340 * of the range (useful for common systems such as a vsite-protein
2341 * in 3-site water).
2342 * With domain decomposition, as long as the vsites are distributed
2343 * uniformly in each domain along the major dimension, usually x,
2344 * it will also perform well.
2346 if (!vsite->useDomdec)
2348 vsite_atom_range = -1;
2349 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2351 { // TODO remove me
2352 if (ftype != F_VSITEN)
2354 int nral1 = 1 + NRAL(ftype);
2355 const t_iatom *iat = ilist[ftype].iatoms;
2356 for (int i = 0; i < ilist[ftype].nr; i += nral1)
2358 for (int j = i + 1; j < i + nral1; j++)
2360 vsite_atom_range = std::max(vsite_atom_range, iat[j]);
2364 else
2366 int vs_ind_end;
2368 const t_iatom *iat = ilist[ftype].iatoms;
2370 int i = 0;
2371 while (i < ilist[ftype].nr)
2373 /* The 3 below is from 1+NRAL(ftype)=3 */
2374 vs_ind_end = i + ip[iat[i]].vsiten.n*3;
2376 vsite_atom_range = std::max(vsite_atom_range, iat[i+1]);
2377 while (i < vs_ind_end)
2379 vsite_atom_range = std::max(vsite_atom_range, iat[i+2]);
2380 i += 3;
2386 vsite_atom_range++;
2387 natperthread = (vsite_atom_range + vsite->nthreads - 1)/vsite->nthreads;
2389 else
2391 /* Any local or not local atom could be involved in virtual sites.
2392 * But since we usually have very few non-local virtual sites
2393 * (only non-local vsites that depend on local vsites),
2394 * we distribute the local atom range equally over the threads.
2395 * When assigning vsites to threads, we should take care that the last
2396 * threads also covers the non-local range.
2398 vsite_atom_range = mdatoms->nr;
2399 natperthread = (mdatoms->homenr + vsite->nthreads - 1)/vsite->nthreads;
2402 if (debug)
2404 fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n", mdatoms->nr, vsite_atom_range, natperthread);
2407 /* To simplify the vsite assignment, we make an index which tells us
2408 * to which task particles, both non-vsites and vsites, are assigned.
2410 if (mdatoms->nr > vsite->taskIndexNalloc)
2412 vsite->taskIndexNalloc = over_alloc_large(mdatoms->nr);
2413 srenew(vsite->taskIndex, vsite->taskIndexNalloc);
2416 /* Initialize the task index array. Here we assign the non-vsite
2417 * particles to task=thread, so we easily figure out if vsites
2418 * depend on local and/or non-local particles in assignVsitesToThread.
2420 int *taskIndex = vsite->taskIndex;
2422 int thread = 0;
2423 for (int i = 0; i < mdatoms->nr; i++)
2425 if (mdatoms->ptype[i] == eptVSite)
2427 /* vsites are not assigned to a task yet */
2428 taskIndex[i] = -1;
2430 else
2432 /* assign non-vsite particles to task thread */
2433 taskIndex[i] = thread;
2435 if (i == (thread + 1)*natperthread && thread < vsite->nthreads)
2437 thread++;
2442 #pragma omp parallel num_threads(vsite->nthreads)
2446 int thread = gmx_omp_get_thread_num();
2447 VsiteThread *tData = vsite->tData[thread];
2449 /* Clear the buffer use flags that were set before */
2450 if (tData->useInterdependentTask)
2452 InterdependentTask *idTask = &tData->idTask;
2454 /* To avoid an extra OpenMP barrier in spread_vsite_f,
2455 * we clear the force buffer at the next step,
2456 * so we need to do it here as well.
2458 clearTaskForceBufferUsedElements(idTask);
2460 idTask->vsite.resize(0);
2461 for (int t = 0; t < vsite->nthreads; t++)
2463 AtomIndex *atomIndex = &idTask->atomIndex[t];
2464 int natom = atomIndex->atom.size();
2465 for (int i = 0; i < natom; i++)
2467 idTask->use[atomIndex->atom[i]] = false;
2469 atomIndex->atom.resize(0);
2471 idTask->nuse = 0;
2474 /* To avoid large f_buf allocations of #threads*vsite_atom_range
2475 * we don't use task2 with more than 200000 atoms. This doesn't
2476 * affect performance, since with such a large range relatively few
2477 * vsites will end up in the separate task.
2478 * Note that useTask2 should be the same for all threads.
2480 tData->useInterdependentTask = (vsite_atom_range <= 200000);
2481 if (tData->useInterdependentTask)
2483 size_t natoms_use_in_vsites = vsite_atom_range;
2484 InterdependentTask *idTask = &tData->idTask;
2485 /* To avoid resizing and re-clearing every nstlist steps,
2486 * we never down size the force buffer.
2488 if (natoms_use_in_vsites > idTask->force.size() ||
2489 natoms_use_in_vsites > idTask->use.size())
2491 idTask->force.resize(natoms_use_in_vsites, { 0, 0, 0 });
2492 idTask->use.resize(natoms_use_in_vsites, false);
2496 /* Assign all vsites that can execute independently on threads */
2497 tData->rangeStart = thread *natperthread;
2498 if (thread < vsite->nthreads - 1)
2500 tData->rangeEnd = (thread + 1)*natperthread;
2502 else
2504 /* The last thread should cover up to the end of the range */
2505 tData->rangeEnd = mdatoms->nr;
2507 assignVsitesToThread(tData,
2508 thread, vsite->nthreads,
2509 natperthread,
2510 taskIndex,
2511 ilist, ip, mdatoms->ptype);
2513 if (tData->useInterdependentTask)
2515 /* In the worst case, all tasks write to force ranges of
2516 * all other tasks, leading to #tasks^2 scaling (this is only
2517 * the overhead, the actual flops remain constant).
2518 * But in most cases there is far less coupling. To improve
2519 * scaling at high thread counts we therefore construct
2520 * an index to only loop over the actually affected tasks.
2522 InterdependentTask *idTask = &tData->idTask;
2524 /* Ensure assignVsitesToThread finished on other threads */
2525 #pragma omp barrier
2527 idTask->spreadTask.resize(0);
2528 idTask->reduceTask.resize(0);
2529 for (int t = 0; t < vsite->nthreads; t++)
2531 /* Do we write to the force buffer of task t? */
2532 if (!idTask->atomIndex[t].atom.empty())
2534 idTask->spreadTask.push_back(t);
2536 /* Does task t write to our force buffer? */
2537 if (!vsite->tData[t]->idTask.atomIndex[thread].atom.empty())
2539 idTask->reduceTask.push_back(t);
2544 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2546 /* Assign all remaining vsites, that will have taskIndex[]=2*vsite->nthreads,
2547 * to a single task that will not run in parallel with other tasks.
2549 assignVsitesToSingleTask(vsite->tData[vsite->nthreads],
2550 2*vsite->nthreads,
2551 taskIndex,
2552 ilist, ip);
2554 if (debug && vsite->nthreads > 1)
2556 fprintf(debug, "virtual site useInterdependentTask %d, nuse:\n",
2557 int{vsite->tData[0]->useInterdependentTask});
2558 for (int th = 0; th < vsite->nthreads + 1; th++)
2560 fprintf(debug, " %4d", vsite->tData[th]->idTask.nuse);
2562 fprintf(debug, "\n");
2564 for (int ftype = c_ftypeVsiteStart; ftype < c_ftypeVsiteEnd; ftype++)
2566 if (ilist[ftype].nr > 0)
2568 fprintf(debug, "%-20s thread dist:",
2569 interaction_function[ftype].longname);
2570 for (int th = 0; th < vsite->nthreads + 1; th++)
2572 fprintf(debug, " %4d %4d ",
2573 vsite->tData[th]->ilist[ftype].nr,
2574 vsite->tData[th]->idTask.ilist[ftype].nr);
2576 fprintf(debug, "\n");
2581 #ifndef NDEBUG
2582 int nrOrig = vsiteIlistNrCount(ilist);
2583 int nrThreaded = 0;
2584 for (int th = 0; th < vsite->nthreads + 1; th++)
2586 nrThreaded +=
2587 vsiteIlistNrCount(vsite->tData[th]->ilist) +
2588 vsiteIlistNrCount(vsite->tData[th]->idTask.ilist);
2590 GMX_ASSERT(nrThreaded == nrOrig, "The number of virtual sites assigned to all thread task has to match the total number of virtual sites");
2591 #endif
2594 void set_vsite_top(gmx_vsite_t *vsite,
2595 const gmx_localtop_t *top,
2596 const t_mdatoms *md)
2598 if (vsite->n_intercg_vsite > 0 && vsite->bHaveChargeGroups)
2600 vsite->vsite_pbc_loc = get_vsite_pbc(top->idef.iparams,
2601 top->idef.il, nullptr, md,
2602 top->cgs);
2605 if (vsite->nthreads > 1)
2607 if (vsite->bHaveChargeGroups)
2609 gmx_fatal(FARGS, "The combination of threading, virtual sites and charge groups is not implemented");
2612 split_vsites_over_threads(top->idef.il, top->idef.iparams,
2613 md, vsite);