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, 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.
41 #include "gromacs/legacyheaders/typedefs.h"
42 #include "gromacs/legacyheaders/types/commrec.h"
43 #include "gromacs/legacyheaders/vsite.h"
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/legacyheaders/nrnb.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/legacyheaders/network.h"
48 #include "gromacs/legacyheaders/domdec.h"
49 #include "gromacs/topology/mtop_util.h"
50 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
52 #include "gromacs/pbcutil/ishift.h"
53 #include "gromacs/pbcutil/mshift.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "gromacs/utility/gmxomp.h"
56 #include "gromacs/utility/smalloc.h"
58 /* Routines to send/recieve coordinates and force
59 * of constructing atoms.
62 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
66 return pbc_dx_aiuc(pbc, xi, xj, dx);
75 /* Vsite construction routines */
77 static void constr_vsite2(rvec xi, rvec xj, rvec x, real a, t_pbc *pbc)
87 pbc_dx_aiuc(pbc, xj, xi, dx);
88 x[XX] = xi[XX] + a*dx[XX];
89 x[YY] = xi[YY] + a*dx[YY];
90 x[ZZ] = xi[ZZ] + a*dx[ZZ];
94 x[XX] = b*xi[XX] + a*xj[XX];
95 x[YY] = b*xi[YY] + a*xj[YY];
96 x[ZZ] = b*xi[ZZ] + a*xj[ZZ];
100 /* TOTAL: 10 flops */
103 static void constr_vsite3(rvec xi, rvec xj, rvec xk, rvec x, real a, real b,
114 pbc_dx_aiuc(pbc, xj, xi, dxj);
115 pbc_dx_aiuc(pbc, xk, xi, dxk);
116 x[XX] = xi[XX] + a*dxj[XX] + b*dxk[XX];
117 x[YY] = xi[YY] + a*dxj[YY] + b*dxk[YY];
118 x[ZZ] = xi[ZZ] + a*dxj[ZZ] + b*dxk[ZZ];
122 x[XX] = c*xi[XX] + a*xj[XX] + b*xk[XX];
123 x[YY] = c*xi[YY] + a*xj[YY] + b*xk[YY];
124 x[ZZ] = c*xi[ZZ] + a*xj[ZZ] + b*xk[ZZ];
128 /* TOTAL: 17 flops */
131 static void constr_vsite3FD(rvec xi, rvec xj, rvec xk, rvec x, real a, real b,
137 pbc_rvec_sub(pbc, xj, xi, xij);
138 pbc_rvec_sub(pbc, xk, xj, xjk);
141 /* temp goes from i to a point on the line jk */
142 temp[XX] = xij[XX] + a*xjk[XX];
143 temp[YY] = xij[YY] + a*xjk[YY];
144 temp[ZZ] = xij[ZZ] + a*xjk[ZZ];
147 c = b*gmx_invsqrt(iprod(temp, temp));
150 x[XX] = xi[XX] + c*temp[XX];
151 x[YY] = xi[YY] + c*temp[YY];
152 x[ZZ] = xi[ZZ] + c*temp[ZZ];
155 /* TOTAL: 34 flops */
158 static void constr_vsite3FAD(rvec xi, rvec xj, rvec xk, rvec x, real a, real b, t_pbc *pbc)
161 real a1, b1, c1, invdij;
163 pbc_rvec_sub(pbc, xj, xi, xij);
164 pbc_rvec_sub(pbc, xk, xj, xjk);
167 invdij = gmx_invsqrt(iprod(xij, xij));
168 c1 = invdij * invdij * iprod(xij, xjk);
169 xp[XX] = xjk[XX] - c1*xij[XX];
170 xp[YY] = xjk[YY] - c1*xij[YY];
171 xp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
173 b1 = b*gmx_invsqrt(iprod(xp, xp));
176 x[XX] = xi[XX] + a1*xij[XX] + b1*xp[XX];
177 x[YY] = xi[YY] + a1*xij[YY] + b1*xp[YY];
178 x[ZZ] = xi[ZZ] + a1*xij[ZZ] + b1*xp[ZZ];
181 /* TOTAL: 63 flops */
184 static void constr_vsite3OUT(rvec xi, rvec xj, rvec xk, rvec x,
185 real a, real b, real c, t_pbc *pbc)
189 pbc_rvec_sub(pbc, xj, xi, xij);
190 pbc_rvec_sub(pbc, xk, xi, xik);
191 cprod(xij, xik, temp);
194 x[XX] = xi[XX] + a*xij[XX] + b*xik[XX] + c*temp[XX];
195 x[YY] = xi[YY] + a*xij[YY] + b*xik[YY] + c*temp[YY];
196 x[ZZ] = xi[ZZ] + a*xij[ZZ] + b*xik[ZZ] + c*temp[ZZ];
199 /* TOTAL: 33 flops */
202 static void constr_vsite4FD(rvec xi, rvec xj, rvec xk, rvec xl, rvec x,
203 real a, real b, real c, t_pbc *pbc)
205 rvec xij, xjk, xjl, temp;
208 pbc_rvec_sub(pbc, xj, xi, xij);
209 pbc_rvec_sub(pbc, xk, xj, xjk);
210 pbc_rvec_sub(pbc, xl, xj, xjl);
213 /* temp goes from i to a point on the plane jkl */
214 temp[XX] = xij[XX] + a*xjk[XX] + b*xjl[XX];
215 temp[YY] = xij[YY] + a*xjk[YY] + b*xjl[YY];
216 temp[ZZ] = xij[ZZ] + a*xjk[ZZ] + b*xjl[ZZ];
219 d = c*gmx_invsqrt(iprod(temp, temp));
222 x[XX] = xi[XX] + d*temp[XX];
223 x[YY] = xi[YY] + d*temp[YY];
224 x[ZZ] = xi[ZZ] + d*temp[ZZ];
227 /* TOTAL: 43 flops */
231 static void constr_vsite4FDN(rvec xi, rvec xj, rvec xk, rvec xl, rvec x,
232 real a, real b, real c, t_pbc *pbc)
234 rvec xij, xik, xil, ra, rb, rja, rjb, rm;
237 pbc_rvec_sub(pbc, xj, xi, xij);
238 pbc_rvec_sub(pbc, xk, xi, xik);
239 pbc_rvec_sub(pbc, xl, xi, xil);
252 rvec_sub(ra, xij, rja);
253 rvec_sub(rb, xij, rjb);
259 d = c*gmx_invsqrt(norm2(rm));
262 x[XX] = xi[XX] + d*rm[XX];
263 x[YY] = xi[YY] + d*rm[YY];
264 x[ZZ] = xi[ZZ] + d*rm[ZZ];
267 /* TOTAL: 47 flops */
271 static int constr_vsiten(t_iatom *ia, t_iparams ip[],
279 n3 = 3*ip[ia[0]].vsiten.n;
282 copy_rvec(x[ai], x1);
284 for (i = 3; i < n3; i += 3)
287 a = ip[ia[i]].vsiten.a;
290 pbc_dx_aiuc(pbc, x[ai], x1, dx);
294 rvec_sub(x[ai], x1, dx);
296 dsum[XX] += a*dx[XX];
297 dsum[YY] += a*dx[YY];
298 dsum[ZZ] += a*dx[ZZ];
302 x[av][XX] = x1[XX] + dsum[XX];
303 x[av][YY] = x1[YY] + dsum[YY];
304 x[av][ZZ] = x1[ZZ] + dsum[ZZ];
310 void construct_vsites_thread(gmx_vsite_t *vsite,
313 t_iparams ip[], t_ilist ilist[],
317 rvec xpbc, xv, vv, dx;
318 real a1, b1, c1, inv_dt;
319 int i, inc, ii, nra, nr, tp, ftype;
320 t_iatom avsite, ai, aj, ak, al, pbc_atom;
323 int *vsite_pbc, ishift;
324 rvec reftmp, vtmp, rtmp;
335 bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
339 for (ftype = 0; (ftype < F_NRE); ftype++)
341 if ((interaction_function[ftype].flags & IF_VSITE) &&
344 nra = interaction_function[ftype].nratoms;
346 nr = ilist[ftype].nr;
347 ia = ilist[ftype].iatoms;
351 pbc_null2 = pbc_null;
353 else if (pbc_null != NULL)
355 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
358 for (i = 0; i < nr; )
362 /* The vsite and constructing atoms */
366 /* Constants for constructing vsites */
368 /* Check what kind of pbc we need to use */
371 /* No charge groups, vsite follows its own pbc */
373 copy_rvec(x[avsite], xpbc);
375 else if (vsite_pbc != NULL)
377 pbc_atom = vsite_pbc[i/(1+nra)];
382 /* We need to copy the coordinates here,
383 * single for single atom cg's pbc_atom
384 * is the vsite itself.
386 copy_rvec(x[pbc_atom], xpbc);
388 pbc_null2 = pbc_null;
399 /* Copy the old position */
400 copy_rvec(x[avsite], xv);
402 /* Construct the vsite depending on type */
407 constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
413 constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
419 constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
425 constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
432 constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
440 constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
449 constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
453 inc = constr_vsiten(ia, ip, x, pbc_null2);
456 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
457 ftype, __FILE__, __LINE__);
462 /* Match the pbc of this vsite to the rest of its charge group */
463 ishift = pbc_dx_aiuc(pbc_null, x[avsite], xpbc, dx);
464 if (ishift != CENTRAL)
466 rvec_add(xpbc, dx, x[avsite]);
471 /* Calculate velocity of vsite... */
472 rvec_sub(x[avsite], xv, vv);
473 svmul(inv_dt, vv, v[avsite]);
476 /* Increment loop variables */
484 void construct_vsites(gmx_vsite_t *vsite,
487 t_iparams ip[], t_ilist ilist[],
488 int ePBC, gmx_bool bMolPBC,
489 t_commrec *cr, matrix box)
491 t_pbc pbc, *pbc_null;
495 bDomDec = cr && DOMAINDECOMP(cr);
497 /* We only need to do pbc when we have inter-cg vsites */
498 if (ePBC != epbcNONE && (bDomDec || bMolPBC) && vsite->n_intercg_vsite)
500 /* This is wasting some CPU time as we now do this multiple times
501 * per MD step. But how often do we have vsites with full pbc?
503 pbc_null = set_pbc_dd(&pbc, ePBC, cr != NULL ? cr->dd : NULL, FALSE, box);
512 dd_move_x_vsites(cr->dd, box, x);
515 if (vsite->nthreads == 1)
517 construct_vsites_thread(vsite,
524 #pragma omp parallel num_threads(vsite->nthreads)
526 construct_vsites_thread(vsite,
528 ip, vsite->tdata[gmx_omp_get_thread_num()].ilist,
531 /* Now we can construct the vsites that might depend on other vsites */
532 construct_vsites_thread(vsite,
534 ip, vsite->tdata[vsite->nthreads].ilist,
539 static void spread_vsite2(t_iatom ia[], real a,
540 rvec x[], rvec f[], rvec fshift[],
541 t_pbc *pbc, t_graph *g)
553 svmul(1-a, f[av], fi);
554 svmul( a, f[av], fj);
563 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
565 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
570 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
571 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
579 if (fshift && (siv != CENTRAL || sij != CENTRAL))
581 rvec_inc(fshift[siv], f[av]);
582 rvec_dec(fshift[CENTRAL], fi);
583 rvec_dec(fshift[sij], fj);
586 /* TOTAL: 13 flops */
589 void construct_vsites_mtop(gmx_vsite_t *vsite,
590 gmx_mtop_t *mtop, rvec x[])
593 gmx_molblock_t *molb;
597 for (mb = 0; mb < mtop->nmolblock; mb++)
599 molb = &mtop->molblock[mb];
600 molt = &mtop->moltype[molb->type];
601 for (mol = 0; mol < molb->nmol; mol++)
603 construct_vsites(vsite, x+as, 0.0, NULL,
604 mtop->ffparams.iparams, molt->ilist,
605 epbcNONE, TRUE, NULL, NULL);
606 as += molt->atoms.nr;
611 static void spread_vsite3(t_iatom ia[], real a, real b,
612 rvec x[], rvec f[], rvec fshift[],
613 t_pbc *pbc, t_graph *g)
616 atom_id av, ai, aj, ak;
625 svmul(1-a-b, f[av], fi);
626 svmul( a, f[av], fj);
627 svmul( b, f[av], fk);
637 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ia[1]), di);
639 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
641 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), di);
646 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
647 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
648 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
657 if (fshift && (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL))
659 rvec_inc(fshift[siv], f[av]);
660 rvec_dec(fshift[CENTRAL], fi);
661 rvec_dec(fshift[sij], fj);
662 rvec_dec(fshift[sik], fk);
665 /* TOTAL: 20 flops */
668 static void spread_vsite3FD(t_iatom ia[], real a, real b,
669 rvec x[], rvec f[], rvec fshift[],
670 gmx_bool VirCorr, matrix dxdf,
671 t_pbc *pbc, t_graph *g)
673 real fx, fy, fz, c, invl, fproj, a1;
674 rvec xvi, xij, xjk, xix, fv, temp;
675 t_iatom av, ai, aj, ak;
676 int svi, sji, skj, d;
683 copy_rvec(f[av], fv);
685 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
686 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
689 /* xix goes from i to point x on the line jk */
690 xix[XX] = xij[XX]+a*xjk[XX];
691 xix[YY] = xij[YY]+a*xjk[YY];
692 xix[ZZ] = xij[ZZ]+a*xjk[ZZ];
695 invl = gmx_invsqrt(iprod(xix, xix));
699 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
701 temp[XX] = c*(fv[XX]-fproj*xix[XX]);
702 temp[YY] = c*(fv[YY]-fproj*xix[YY]);
703 temp[ZZ] = c*(fv[ZZ]-fproj*xix[ZZ]);
706 /* c is already calculated in constr_vsite3FD
707 storing c somewhere will save 26 flops! */
710 f[ai][XX] += fv[XX] - temp[XX];
711 f[ai][YY] += fv[YY] - temp[YY];
712 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
713 f[aj][XX] += a1*temp[XX];
714 f[aj][YY] += a1*temp[YY];
715 f[aj][ZZ] += a1*temp[ZZ];
716 f[ak][XX] += a*temp[XX];
717 f[ak][YY] += a*temp[YY];
718 f[ak][ZZ] += a*temp[ZZ];
723 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
725 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
727 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
732 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
739 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
741 rvec_dec(fshift[svi], fv);
742 fshift[CENTRAL][XX] += fv[XX] - (1 + a)*temp[XX];
743 fshift[CENTRAL][YY] += fv[YY] - (1 + a)*temp[YY];
744 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a)*temp[ZZ];
745 fshift[ sji][XX] += temp[XX];
746 fshift[ sji][YY] += temp[YY];
747 fshift[ sji][ZZ] += temp[ZZ];
748 fshift[ skj][XX] += a*temp[XX];
749 fshift[ skj][YY] += a*temp[YY];
750 fshift[ skj][ZZ] += a*temp[ZZ];
755 /* When VirCorr=TRUE, the virial for the current forces is not
756 * calculated from the redistributed forces. This means that
757 * the effect of non-linear virtual site constructions on the virial
758 * needs to be added separately. This contribution can be calculated
759 * in many ways, but the simplest and cheapest way is to use
760 * the first constructing atom ai as a reference position in space:
761 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
766 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
768 for (i = 0; i < DIM; i++)
770 for (j = 0; j < DIM; j++)
772 /* As xix is a linear combination of j and k, use that here */
773 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
778 /* TOTAL: 61 flops */
781 static void spread_vsite3FAD(t_iatom ia[], real a, real b,
782 rvec x[], rvec f[], rvec fshift[],
783 gmx_bool VirCorr, matrix dxdf,
784 t_pbc *pbc, t_graph *g)
786 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
787 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
788 t_iatom av, ai, aj, ak;
789 int svi, sji, skj, d;
796 copy_rvec(f[ia[1]], fv);
798 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
799 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
802 invdij = gmx_invsqrt(iprod(xij, xij));
803 invdij2 = invdij * invdij;
804 c1 = iprod(xij, xjk) * invdij2;
805 xperp[XX] = xjk[XX] - c1*xij[XX];
806 xperp[YY] = xjk[YY] - c1*xij[YY];
807 xperp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
808 /* xperp in plane ijk, perp. to ij */
809 invdp = gmx_invsqrt(iprod(xperp, xperp));
814 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
815 storing them somewhere will save 45 flops! */
817 fproj = iprod(xij, fv)*invdij2;
818 svmul(fproj, xij, Fpij); /* proj. f on xij */
819 svmul(iprod(xperp, fv)*invdp*invdp, xperp, Fppp); /* proj. f on xperp */
820 svmul(b1*fproj, xperp, f3);
823 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
824 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
825 for (d = 0; (d < DIM); d++)
833 f[ai][XX] += fv[XX] - f1[XX] + c1*f2[XX] + f3[XX];
834 f[ai][YY] += fv[YY] - f1[YY] + c1*f2[YY] + f3[YY];
835 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1*f2[ZZ] + f3[ZZ];
836 f[aj][XX] += f1[XX] - c2*f2[XX] - f3[XX];
837 f[aj][YY] += f1[YY] - c2*f2[YY] - f3[YY];
838 f[aj][ZZ] += f1[ZZ] - c2*f2[ZZ] - f3[ZZ];
846 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
848 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
850 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
855 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
862 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
864 rvec_dec(fshift[svi], fv);
865 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1-c1)*f2[XX] + f3[XX];
866 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1-c1)*f2[YY] + f3[YY];
867 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1-c1)*f2[ZZ] + f3[ZZ];
868 fshift[ sji][XX] += f1[XX] - c1 *f2[XX] - f3[XX];
869 fshift[ sji][YY] += f1[YY] - c1 *f2[YY] - f3[YY];
870 fshift[ sji][ZZ] += f1[ZZ] - c1 *f2[ZZ] - f3[ZZ];
871 fshift[ skj][XX] += f2[XX];
872 fshift[ skj][YY] += f2[YY];
873 fshift[ skj][ZZ] += f2[ZZ];
881 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
883 for (i = 0; i < DIM; i++)
885 for (j = 0; j < DIM; j++)
887 /* Note that xik=xij+xjk, so we have to add xij*f2 */
890 + xij[i]*(f1[j] + (1 - c2)*f2[j] - f3[j])
896 /* TOTAL: 113 flops */
899 static void spread_vsite3OUT(t_iatom ia[], real a, real b, real c,
900 rvec x[], rvec f[], rvec fshift[],
901 gmx_bool VirCorr, matrix dxdf,
902 t_pbc *pbc, t_graph *g)
904 rvec xvi, xij, xik, fv, fj, fk;
906 atom_id av, ai, aj, ak;
915 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
916 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
919 copy_rvec(f[av], fv);
926 fj[XX] = a*fv[XX] - xik[ZZ]*cfy + xik[YY]*cfz;
927 fj[YY] = xik[ZZ]*cfx + a*fv[YY] - xik[XX]*cfz;
928 fj[ZZ] = -xik[YY]*cfx + xik[XX]*cfy + a*fv[ZZ];
930 fk[XX] = b*fv[XX] + xij[ZZ]*cfy - xij[YY]*cfz;
931 fk[YY] = -xij[ZZ]*cfx + b*fv[YY] + xij[XX]*cfz;
932 fk[ZZ] = xij[YY]*cfx - xij[XX]*cfy + b*fv[ZZ];
935 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
936 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
937 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
944 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
946 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
948 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
953 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
960 if (fshift && (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL))
962 rvec_dec(fshift[svi], fv);
963 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
964 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
965 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
966 rvec_inc(fshift[sji], fj);
967 rvec_inc(fshift[ski], fk);
975 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
977 for (i = 0; i < DIM; i++)
979 for (j = 0; j < DIM; j++)
981 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j];
986 /* TOTAL: 54 flops */
989 static void spread_vsite4FD(t_iatom ia[], real a, real b, real c,
990 rvec x[], rvec f[], rvec fshift[],
991 gmx_bool VirCorr, matrix dxdf,
992 t_pbc *pbc, t_graph *g)
994 real d, invl, fproj, a1;
995 rvec xvi, xij, xjk, xjl, xix, fv, temp;
996 atom_id av, ai, aj, ak, al;
998 int svi, sji, skj, slj, m;
1006 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1007 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1008 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1011 /* xix goes from i to point x on the plane jkl */
1012 for (m = 0; m < DIM; m++)
1014 xix[m] = xij[m] + a*xjk[m] + b*xjl[m];
1018 invl = gmx_invsqrt(iprod(xix, xix));
1020 /* 4 + ?10? flops */
1022 copy_rvec(f[av], fv);
1024 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
1026 for (m = 0; m < DIM; m++)
1028 temp[m] = d*(fv[m] - fproj*xix[m]);
1032 /* c is already calculated in constr_vsite3FD
1033 storing c somewhere will save 35 flops! */
1036 for (m = 0; m < DIM; m++)
1038 f[ai][m] += fv[m] - temp[m];
1039 f[aj][m] += a1*temp[m];
1040 f[ak][m] += a*temp[m];
1041 f[al][m] += b*temp[m];
1047 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1049 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1051 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
1053 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, aj), di);
1058 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1066 (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL))
1068 rvec_dec(fshift[svi], fv);
1069 for (m = 0; m < DIM; m++)
1071 fshift[CENTRAL][m] += fv[m] - (1 + a + b)*temp[m];
1072 fshift[ sji][m] += temp[m];
1073 fshift[ skj][m] += a*temp[m];
1074 fshift[ slj][m] += b*temp[m];
1083 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1085 for (i = 0; i < DIM; i++)
1087 for (j = 0; j < DIM; j++)
1089 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
1094 /* TOTAL: 77 flops */
1098 static void spread_vsite4FDN(t_iatom ia[], real a, real b, real c,
1099 rvec x[], rvec f[], rvec fshift[],
1100 gmx_bool VirCorr, matrix dxdf,
1101 t_pbc *pbc, t_graph *g)
1103 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1104 rvec fv, fj, fk, fl;
1108 int av, ai, aj, ak, al;
1109 int svi, sij, sik, sil;
1111 /* DEBUG: check atom indices */
1118 copy_rvec(f[av], fv);
1120 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1121 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1122 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1135 rvec_sub(ra, xij, rja);
1136 rvec_sub(rb, xij, rjb);
1137 rvec_sub(rb, ra, rab);
1140 cprod(rja, rjb, rm);
1143 invrm = gmx_invsqrt(norm2(rm));
1144 denom = invrm*invrm;
1147 cfx = c*invrm*fv[XX];
1148 cfy = c*invrm*fv[YY];
1149 cfz = c*invrm*fv[ZZ];
1160 fj[XX] = ( -rm[XX]*rt[XX]) * cfx + ( rab[ZZ]-rm[YY]*rt[XX]) * cfy + (-rab[YY]-rm[ZZ]*rt[XX]) * cfz;
1161 fj[YY] = (-rab[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( rab[XX]-rm[ZZ]*rt[YY]) * cfz;
1162 fj[ZZ] = ( rab[YY]-rm[XX]*rt[ZZ]) * cfx + (-rab[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1173 fk[XX] = ( -rm[XX]*rt[XX]) * cfx + (-a*rjb[ZZ]-rm[YY]*rt[XX]) * cfy + ( a*rjb[YY]-rm[ZZ]*rt[XX]) * cfz;
1174 fk[YY] = ( a*rjb[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + (-a*rjb[XX]-rm[ZZ]*rt[YY]) * cfz;
1175 fk[ZZ] = (-a*rjb[YY]-rm[XX]*rt[ZZ]) * cfx + ( a*rjb[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1186 fl[XX] = ( -rm[XX]*rt[XX]) * cfx + ( b*rja[ZZ]-rm[YY]*rt[XX]) * cfy + (-b*rja[YY]-rm[ZZ]*rt[XX]) * cfz;
1187 fl[YY] = (-b*rja[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( b*rja[XX]-rm[ZZ]*rt[YY]) * cfz;
1188 fl[ZZ] = ( b*rja[YY]-rm[XX]*rt[ZZ]) * cfx + (-b*rja[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1191 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1192 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1193 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1194 rvec_inc(f[aj], fj);
1195 rvec_inc(f[ak], fk);
1196 rvec_inc(f[al], fl);
1201 ivec_sub(SHIFT_IVEC(g, av), SHIFT_IVEC(g, ai), di);
1203 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1205 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1207 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, ai), di);
1212 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1219 if (fshift && (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL))
1221 rvec_dec(fshift[svi], fv);
1222 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1223 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1224 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1225 rvec_inc(fshift[sij], fj);
1226 rvec_inc(fshift[sik], fk);
1227 rvec_inc(fshift[sil], fl);
1235 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1237 for (i = 0; i < DIM; i++)
1239 for (j = 0; j < DIM; j++)
1241 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j] + xil[i]*fl[j];
1246 /* Total: 207 flops (Yuck!) */
1250 static int spread_vsiten(t_iatom ia[], t_iparams ip[],
1251 rvec x[], rvec f[], rvec fshift[],
1252 t_pbc *pbc, t_graph *g)
1260 n3 = 3*ip[ia[0]].vsiten.n;
1262 copy_rvec(x[av], xv);
1264 for (i = 0; i < n3; i += 3)
1269 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
1274 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1280 a = ip[ia[i]].vsiten.a;
1281 svmul(a, f[av], fi);
1282 rvec_inc(f[ai], fi);
1283 if (fshift && siv != CENTRAL)
1285 rvec_inc(fshift[siv], fi);
1286 rvec_dec(fshift[CENTRAL], fi);
1295 static int vsite_count(const t_ilist *ilist, int ftype)
1297 if (ftype == F_VSITEN)
1299 return ilist[ftype].nr/3;
1303 return ilist[ftype].nr/(1 + interaction_function[ftype].nratoms);
1307 static void spread_vsite_f_thread(gmx_vsite_t *vsite,
1308 rvec x[], rvec f[], rvec *fshift,
1309 gmx_bool VirCorr, matrix dxdf,
1310 t_iparams ip[], t_ilist ilist[],
1311 t_graph *g, t_pbc *pbc_null)
1315 int i, inc, m, nra, nr, tp, ftype;
1325 bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
1327 /* this loop goes backwards to be able to build *
1328 * higher type vsites from lower types */
1331 for (ftype = F_NRE-1; (ftype >= 0); ftype--)
1333 if ((interaction_function[ftype].flags & IF_VSITE) &&
1334 ilist[ftype].nr > 0)
1336 nra = interaction_function[ftype].nratoms;
1338 nr = ilist[ftype].nr;
1339 ia = ilist[ftype].iatoms;
1343 pbc_null2 = pbc_null;
1345 else if (pbc_null != NULL)
1347 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
1350 for (i = 0; i < nr; )
1352 if (vsite_pbc != NULL)
1354 if (vsite_pbc[i/(1+nra)] > -2)
1356 pbc_null2 = pbc_null;
1366 /* Constants for constructing */
1367 a1 = ip[tp].vsite.a;
1368 /* Construct the vsite depending on type */
1372 spread_vsite2(ia, a1, x, f, fshift, pbc_null2, g);
1375 b1 = ip[tp].vsite.b;
1376 spread_vsite3(ia, a1, b1, x, f, fshift, pbc_null2, g);
1379 b1 = ip[tp].vsite.b;
1380 spread_vsite3FD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1383 b1 = ip[tp].vsite.b;
1384 spread_vsite3FAD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1387 b1 = ip[tp].vsite.b;
1388 c1 = ip[tp].vsite.c;
1389 spread_vsite3OUT(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1392 b1 = ip[tp].vsite.b;
1393 c1 = ip[tp].vsite.c;
1394 spread_vsite4FD(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1397 b1 = ip[tp].vsite.b;
1398 c1 = ip[tp].vsite.c;
1399 spread_vsite4FDN(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1402 inc = spread_vsiten(ia, ip, x, f, fshift, pbc_null2, g);
1405 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
1406 ftype, __FILE__, __LINE__);
1408 clear_rvec(f[ia[1]]);
1410 /* Increment loop variables */
1418 void spread_vsite_f(gmx_vsite_t *vsite,
1419 rvec x[], rvec f[], rvec *fshift,
1420 gmx_bool VirCorr, matrix vir,
1421 t_nrnb *nrnb, t_idef *idef,
1422 int ePBC, gmx_bool bMolPBC, t_graph *g, matrix box,
1425 t_pbc pbc, *pbc_null;
1428 /* We only need to do pbc when we have inter-cg vsites */
1429 if ((DOMAINDECOMP(cr) || bMolPBC) && vsite->n_intercg_vsite)
1431 /* This is wasting some CPU time as we now do this multiple times
1432 * per MD step. But how often do we have vsites with full pbc?
1434 pbc_null = set_pbc_dd(&pbc, ePBC, cr->dd, FALSE, box);
1441 if (DOMAINDECOMP(cr))
1443 dd_clear_f_vsites(cr->dd, f);
1446 if (vsite->nthreads == 1)
1448 spread_vsite_f_thread(vsite,
1450 VirCorr, vsite->tdata[0].dxdf,
1451 idef->iparams, idef->il,
1456 /* First spread the vsites that might depend on other vsites */
1457 spread_vsite_f_thread(vsite,
1459 VirCorr, vsite->tdata[vsite->nthreads].dxdf,
1461 vsite->tdata[vsite->nthreads].ilist,
1464 #pragma omp parallel num_threads(vsite->nthreads)
1469 thread = gmx_omp_get_thread_num();
1471 if (thread == 0 || fshift == NULL)
1479 fshift_t = vsite->tdata[thread].fshift;
1481 for (i = 0; i < SHIFTS; i++)
1483 clear_rvec(fshift_t[i]);
1487 spread_vsite_f_thread(vsite,
1489 VirCorr, vsite->tdata[thread].dxdf,
1491 vsite->tdata[thread].ilist,
1499 for (th = 1; th < vsite->nthreads; th++)
1501 for (i = 0; i < SHIFTS; i++)
1503 rvec_inc(fshift[i], vsite->tdata[th].fshift[i]);
1513 for (th = 0; th < (vsite->nthreads == 1 ? 1 : vsite->nthreads+1); th++)
1515 for (i = 0; i < DIM; i++)
1517 for (j = 0; j < DIM; j++)
1519 vir[i][j] += -0.5*vsite->tdata[th].dxdf[i][j];
1525 if (DOMAINDECOMP(cr))
1527 dd_move_f_vsites(cr->dd, f, fshift);
1530 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(idef->il, F_VSITE2));
1531 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(idef->il, F_VSITE3));
1532 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(idef->il, F_VSITE3FD));
1533 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(idef->il, F_VSITE3FAD));
1534 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(idef->il, F_VSITE3OUT));
1535 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(idef->il, F_VSITE4FD));
1536 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(idef->il, F_VSITE4FDN));
1537 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(idef->il, F_VSITEN));
1540 static int *atom2cg(t_block *cgs)
1544 snew(a2cg, cgs->index[cgs->nr]);
1545 for (cg = 0; cg < cgs->nr; cg++)
1547 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
1556 static int count_intercg_vsite(gmx_mtop_t *mtop,
1557 gmx_bool *bHaveChargeGroups)
1559 int mb, mt, ftype, nral, i, cg, a;
1560 gmx_molblock_t *molb;
1561 gmx_moltype_t *molt;
1565 int n_intercg_vsite;
1567 *bHaveChargeGroups = FALSE;
1569 n_intercg_vsite = 0;
1570 for (mb = 0; mb < mtop->nmolblock; mb++)
1572 molb = &mtop->molblock[mb];
1573 molt = &mtop->moltype[molb->type];
1575 if (molt->cgs.nr < molt->atoms.nr)
1577 *bHaveChargeGroups = TRUE;
1580 a2cg = atom2cg(&molt->cgs);
1581 for (ftype = 0; ftype < F_NRE; ftype++)
1583 if (interaction_function[ftype].flags & IF_VSITE)
1586 il = &molt->ilist[ftype];
1588 for (i = 0; i < il->nr; i += 1+nral)
1591 for (a = 1; a < nral; a++)
1593 if (a2cg[ia[1+a]] != cg)
1595 n_intercg_vsite += molb->nmol;
1605 return n_intercg_vsite;
1608 static int **get_vsite_pbc(t_iparams *iparams, t_ilist *ilist,
1609 t_atom *atom, t_mdatoms *md,
1610 t_block *cgs, int *a2cg)
1612 int ftype, nral, i, j, vsi, vsite, cg_v, cg_c, a, nc3 = 0;
1615 int **vsite_pbc, *vsite_pbc_f;
1617 gmx_bool bViteOnlyCG_and_FirstAtom;
1619 /* Make an array that tells if the pbc of an atom is set */
1620 snew(pbc_set, cgs->index[cgs->nr]);
1621 /* PBC is set for all non vsites */
1622 for (a = 0; a < cgs->index[cgs->nr]; a++)
1624 if ((atom && atom[a].ptype != eptVSite) ||
1625 (md && md->ptype[a] != eptVSite))
1631 snew(vsite_pbc, F_VSITEN-F_VSITE2+1);
1633 for (ftype = 0; ftype < F_NRE; ftype++)
1635 if (interaction_function[ftype].flags & IF_VSITE)
1641 snew(vsite_pbc[ftype-F_VSITE2], il->nr/(1+nral));
1642 vsite_pbc_f = vsite_pbc[ftype-F_VSITE2];
1650 /* A value of -2 signals that this vsite and its contructing
1651 * atoms are all within the same cg, so no pbc is required.
1653 vsite_pbc_f[vsi] = -2;
1654 /* Check if constructing atoms are outside the vsite's cg */
1656 if (ftype == F_VSITEN)
1658 nc3 = 3*iparams[ia[i]].vsiten.n;
1659 for (j = 0; j < nc3; j += 3)
1661 if (a2cg[ia[i+j+2]] != cg_v)
1663 vsite_pbc_f[vsi] = -1;
1669 for (a = 1; a < nral; a++)
1671 if (a2cg[ia[i+1+a]] != cg_v)
1673 vsite_pbc_f[vsi] = -1;
1677 if (vsite_pbc_f[vsi] == -1)
1679 /* Check if this is the first processed atom of a vsite only cg */
1680 bViteOnlyCG_and_FirstAtom = TRUE;
1681 for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1683 /* Non-vsites already have pbc set, so simply check for pbc_set */
1686 bViteOnlyCG_and_FirstAtom = FALSE;
1690 if (bViteOnlyCG_and_FirstAtom)
1692 /* First processed atom of a vsite only charge group.
1693 * The pbc of the input coordinates to construct_vsites
1694 * should be preserved.
1696 vsite_pbc_f[vsi] = vsite;
1698 else if (cg_v != a2cg[ia[1+i+1]])
1700 /* This vsite has a different charge group index
1701 * than it's first constructing atom
1702 * and the charge group has more than one atom,
1703 * search for the first normal particle
1704 * or vsite that already had its pbc defined.
1705 * If nothing is found, use full pbc for this vsite.
1707 for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1709 if (a != vsite && pbc_set[a])
1711 vsite_pbc_f[vsi] = a;
1714 fprintf(debug, "vsite %d match pbc with atom %d\n",
1722 fprintf(debug, "vsite atom %d cg %d - %d pbc atom %d\n",
1723 vsite+1, cgs->index[cg_v]+1, cgs->index[cg_v+1],
1724 vsite_pbc_f[vsi]+1);
1728 if (ftype == F_VSITEN)
1730 /* The other entries in vsite_pbc_f are not used for center vsites */
1738 /* This vsite now has its pbc defined */
1750 gmx_vsite_t *init_vsite(gmx_mtop_t *mtop, t_commrec *cr,
1751 gmx_bool bSerial_NoPBC)
1757 gmx_moltype_t *molt;
1760 /* check if there are vsites */
1762 for (i = 0; i < F_NRE; i++)
1764 if (interaction_function[i].flags & IF_VSITE)
1766 nvsite += gmx_mtop_ftype_count(mtop, i);
1777 vsite->n_intercg_vsite = count_intercg_vsite(mtop,
1778 &vsite->bHaveChargeGroups);
1780 /* If we don't have charge groups, the vsite follows its own pbc */
1781 if (!bSerial_NoPBC &&
1782 vsite->bHaveChargeGroups &&
1783 vsite->n_intercg_vsite > 0 && DOMAINDECOMP(cr))
1785 vsite->nvsite_pbc_molt = mtop->nmoltype;
1786 snew(vsite->vsite_pbc_molt, vsite->nvsite_pbc_molt);
1787 for (mt = 0; mt < mtop->nmoltype; mt++)
1789 molt = &mtop->moltype[mt];
1790 /* Make an atom to charge group index */
1791 a2cg = atom2cg(&molt->cgs);
1792 vsite->vsite_pbc_molt[mt] = get_vsite_pbc(mtop->ffparams.iparams,
1794 molt->atoms.atom, NULL,
1799 snew(vsite->vsite_pbc_loc_nalloc, F_VSITEN-F_VSITE2+1);
1800 snew(vsite->vsite_pbc_loc, F_VSITEN-F_VSITE2+1);
1805 vsite->nthreads = 1;
1809 vsite->nthreads = gmx_omp_nthreads_get(emntVSITE);
1813 /* We need one extra thread data structure for the overlap vsites */
1814 snew(vsite->tdata, vsite->nthreads+1);
1817 vsite->th_ind = NULL;
1818 vsite->th_ind_nalloc = 0;
1823 static void prepare_vsite_thread(const t_ilist *ilist,
1824 gmx_vsite_thread_t *vsite_th)
1828 for (ftype = 0; ftype < F_NRE; ftype++)
1830 if (interaction_function[ftype].flags & IF_VSITE)
1832 if (ilist[ftype].nr > vsite_th->ilist[ftype].nalloc)
1834 vsite_th->ilist[ftype].nalloc = over_alloc_large(ilist[ftype].nr);
1835 srenew(vsite_th->ilist[ftype].iatoms, vsite_th->ilist[ftype].nalloc);
1838 vsite_th->ilist[ftype].nr = 0;
1843 void split_vsites_over_threads(const t_ilist *ilist,
1844 const t_mdatoms *mdatoms,
1845 gmx_bool bLimitRange,
1849 int vsite_atom_range, natperthread;
1854 int nral1, inc, i, j;
1856 if (vsite->nthreads == 1)
1862 #pragma omp parallel for num_threads(vsite->nthreads) schedule(static)
1863 for (th = 0; th < vsite->nthreads; th++)
1865 prepare_vsite_thread(ilist, &vsite->tdata[th]);
1867 /* Master threads does the (potential) overlap vsites */
1868 prepare_vsite_thread(ilist, &vsite->tdata[vsite->nthreads]);
1870 /* The current way of distributing the vsites over threads in primitive.
1871 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
1872 * without taking into account how the vsites are distributed.
1873 * Without domain decomposition we bLimitRange=TRUE and we at least
1874 * tighten the upper bound of the range (useful for common systems
1875 * such as a vsite-protein in 3-site water).
1879 vsite_atom_range = -1;
1880 for (ftype = 0; ftype < F_NRE; ftype++)
1882 if ((interaction_function[ftype].flags & IF_VSITE) &&
1885 nral1 = 1 + NRAL(ftype);
1886 iat = ilist[ftype].iatoms;
1887 for (i = 0; i < ilist[ftype].nr; i += nral1)
1889 for (j = i+1; j < i+nral1; j++)
1891 vsite_atom_range = max(vsite_atom_range, iat[j]);
1900 vsite_atom_range = mdatoms->homenr;
1902 natperthread = (vsite_atom_range + vsite->nthreads - 1)/vsite->nthreads;
1906 fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n", mdatoms->nr, vsite_atom_range, natperthread);
1909 /* To simplify the vsite assignment, we make an index which tells us
1910 * to which thread particles, both non-vsites and vsites, are assigned.
1912 if (mdatoms->nr > vsite->th_ind_nalloc)
1914 vsite->th_ind_nalloc = over_alloc_large(mdatoms->nr);
1915 srenew(vsite->th_ind, vsite->th_ind_nalloc);
1917 th_ind = vsite->th_ind;
1919 for (i = 0; i < mdatoms->nr; i++)
1921 if (mdatoms->ptype[i] == eptVSite)
1923 /* vsites are not assigned to a thread yet */
1928 /* assign non-vsite particles to thread th */
1931 if (i == (th + 1)*natperthread && th < vsite->nthreads)
1937 for (ftype = 0; ftype < F_NRE; ftype++)
1939 if ((interaction_function[ftype].flags & IF_VSITE) &&
1942 nral1 = 1 + NRAL(ftype);
1944 iat = ilist[ftype].iatoms;
1945 for (i = 0; i < ilist[ftype].nr; )
1947 th = iat[1+i]/natperthread;
1948 /* We would like to assign this vsite the thread th,
1949 * but it might depend on atoms outside the atom range of th
1950 * or on another vsite not assigned to thread th.
1952 if (ftype != F_VSITEN)
1954 for (j = i+2; j < i+nral1; j++)
1956 if (th_ind[iat[j]] != th)
1958 /* Some constructing atoms are not assigned to
1959 * thread th, move this vsite to a separate batch.
1961 th = vsite->nthreads;
1968 for (j = i+2; j < i+inc; j += 3)
1970 if (th_ind[iat[j]] != th)
1972 th = vsite->nthreads;
1976 /* Copy this vsite to the thread data struct of thread th */
1977 il_th = &vsite->tdata[th].ilist[ftype];
1978 for (j = i; j < i+inc; j++)
1980 il_th->iatoms[il_th->nr++] = iat[j];
1982 /* Update this vsite's thread index entry */
1983 th_ind[iat[1+i]] = th;
1992 for (ftype = 0; ftype < F_NRE; ftype++)
1994 if ((interaction_function[ftype].flags & IF_VSITE) &&
1995 ilist[ftype].nr > 0)
1997 fprintf(debug, "%-20s thread dist:",
1998 interaction_function[ftype].longname);
1999 for (th = 0; th < vsite->nthreads+1; th++)
2001 fprintf(debug, " %4d", vsite->tdata[th].ilist[ftype].nr);
2003 fprintf(debug, "\n");
2009 void set_vsite_top(gmx_vsite_t *vsite, gmx_localtop_t *top, t_mdatoms *md,
2014 if (vsite->n_intercg_vsite > 0)
2016 if (vsite->bHaveChargeGroups)
2018 /* Make an atom to charge group index */
2019 a2cg = atom2cg(&top->cgs);
2020 vsite->vsite_pbc_loc = get_vsite_pbc(top->idef.iparams,
2021 top->idef.il, NULL, md,
2028 if (vsite->nthreads > 1)
2030 if (vsite->bHaveChargeGroups)
2032 gmx_fatal(FARGS, "The combination of threading, virtual sites and charge groups is not implemented");
2035 split_vsites_over_threads(top->idef.il, md, !DOMAINDECOMP(cr), vsite);