1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
53 #include "mtop_util.h"
54 #include "gmx_omp_nthreads.h"
56 #include "gromacs/utility/gmxomp.h"
58 /* Routines to send/recieve coordinates and force
59 * of constructing atoms.
62 static void move_construct_x(t_comm_vsites *vsitecomm, rvec x[], t_commrec *cr)
68 sendbuf = vsitecomm->send_buf;
69 recvbuf = vsitecomm->recv_buf;
72 /* Prepare pulse left by copying to send buffer */
73 for (i = 0; i < vsitecomm->left_export_nconstruct; i++)
75 ia = vsitecomm->left_export_construct[i];
76 copy_rvec(x[ia], sendbuf[i]);
79 /* Pulse coordinates left */
80 gmx_tx_rx_real(cr, GMX_LEFT, (real *)sendbuf, 3*vsitecomm->left_export_nconstruct, GMX_RIGHT, (real *)recvbuf, 3*vsitecomm->right_import_nconstruct);
82 /* Copy from receive buffer to coordinate array */
83 for (i = 0; i < vsitecomm->right_import_nconstruct; i++)
85 ia = vsitecomm->right_import_construct[i];
86 copy_rvec(recvbuf[i], x[ia]);
89 /* Prepare pulse right by copying to send buffer */
90 for (i = 0; i < vsitecomm->right_export_nconstruct; i++)
92 ia = vsitecomm->right_export_construct[i];
93 copy_rvec(x[ia], sendbuf[i]);
96 /* Pulse coordinates right */
97 gmx_tx_rx_real(cr, GMX_RIGHT, (real *)sendbuf, 3*vsitecomm->right_export_nconstruct, GMX_LEFT, (real *)recvbuf, 3*vsitecomm->left_import_nconstruct);
99 /* Copy from receive buffer to coordinate array */
100 for (i = 0; i < vsitecomm->left_import_nconstruct; i++)
102 ia = vsitecomm->left_import_construct[i];
103 copy_rvec(recvbuf[i], x[ia]);
108 static void move_construct_f(t_comm_vsites *vsitecomm, rvec f[], t_commrec *cr)
114 sendbuf = vsitecomm->send_buf;
115 recvbuf = vsitecomm->recv_buf;
117 /* Prepare pulse right by copying to send buffer */
118 for (i = 0; i < vsitecomm->right_import_nconstruct; i++)
120 ia = vsitecomm->right_import_construct[i];
121 copy_rvec(f[ia], sendbuf[i]);
122 clear_rvec(f[ia]); /* Zero it here after moving, just to simplify debug book-keeping... */
125 /* Pulse forces right */
126 gmx_tx_rx_real(cr, GMX_RIGHT, (real *)sendbuf, 3*vsitecomm->right_import_nconstruct, GMX_LEFT, (real *)recvbuf, 3*vsitecomm->left_export_nconstruct);
128 /* Copy from receive buffer to coordinate array */
129 for (i = 0; i < vsitecomm->left_export_nconstruct; i++)
131 ia = vsitecomm->left_export_construct[i];
132 rvec_inc(f[ia], recvbuf[i]);
135 /* Prepare pulse left by copying to send buffer */
136 for (i = 0; i < vsitecomm->left_import_nconstruct; i++)
138 ia = vsitecomm->left_import_construct[i];
139 copy_rvec(f[ia], sendbuf[i]);
140 clear_rvec(f[ia]); /* Zero it here after moving, just to simplify debug book-keeping... */
143 /* Pulse coordinates left */
144 gmx_tx_rx_real(cr, GMX_LEFT, (real *)sendbuf, 3*vsitecomm->left_import_nconstruct, GMX_RIGHT, (real *)recvbuf, 3*vsitecomm->right_export_nconstruct);
146 /* Copy from receive buffer to coordinate array */
147 for (i = 0; i < vsitecomm->right_export_nconstruct; i++)
149 ia = vsitecomm->right_export_construct[i];
150 rvec_inc(f[ia], recvbuf[i]);
153 /* All forces are now on the home processors */
158 pd_clear_nonlocal_constructs(t_comm_vsites *vsitecomm, rvec f[])
162 for (i = 0; i < vsitecomm->left_import_nconstruct; i++)
164 ia = vsitecomm->left_import_construct[i];
167 for (i = 0; i < vsitecomm->right_import_nconstruct; i++)
169 ia = vsitecomm->right_import_construct[i];
176 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
180 return pbc_dx_aiuc(pbc, xi, xj, dx);
184 rvec_sub(xi, xj, dx);
189 /* Vsite construction routines */
191 static void constr_vsite2(rvec xi, rvec xj, rvec x, real a, t_pbc *pbc)
201 pbc_dx_aiuc(pbc, xj, xi, dx);
202 x[XX] = xi[XX] + a*dx[XX];
203 x[YY] = xi[YY] + a*dx[YY];
204 x[ZZ] = xi[ZZ] + a*dx[ZZ];
208 x[XX] = b*xi[XX] + a*xj[XX];
209 x[YY] = b*xi[YY] + a*xj[YY];
210 x[ZZ] = b*xi[ZZ] + a*xj[ZZ];
214 /* TOTAL: 10 flops */
217 static void constr_vsite3(rvec xi, rvec xj, rvec xk, rvec x, real a, real b,
228 pbc_dx_aiuc(pbc, xj, xi, dxj);
229 pbc_dx_aiuc(pbc, xk, xi, dxk);
230 x[XX] = xi[XX] + a*dxj[XX] + b*dxk[XX];
231 x[YY] = xi[YY] + a*dxj[YY] + b*dxk[YY];
232 x[ZZ] = xi[ZZ] + a*dxj[ZZ] + b*dxk[ZZ];
236 x[XX] = c*xi[XX] + a*xj[XX] + b*xk[XX];
237 x[YY] = c*xi[YY] + a*xj[YY] + b*xk[YY];
238 x[ZZ] = c*xi[ZZ] + a*xj[ZZ] + b*xk[ZZ];
242 /* TOTAL: 17 flops */
245 static void constr_vsite3FD(rvec xi, rvec xj, rvec xk, rvec x, real a, real b,
251 pbc_rvec_sub(pbc, xj, xi, xij);
252 pbc_rvec_sub(pbc, xk, xj, xjk);
255 /* temp goes from i to a point on the line jk */
256 temp[XX] = xij[XX] + a*xjk[XX];
257 temp[YY] = xij[YY] + a*xjk[YY];
258 temp[ZZ] = xij[ZZ] + a*xjk[ZZ];
261 c = b*gmx_invsqrt(iprod(temp, temp));
264 x[XX] = xi[XX] + c*temp[XX];
265 x[YY] = xi[YY] + c*temp[YY];
266 x[ZZ] = xi[ZZ] + c*temp[ZZ];
269 /* TOTAL: 34 flops */
272 static void constr_vsite3FAD(rvec xi, rvec xj, rvec xk, rvec x, real a, real b, t_pbc *pbc)
275 real a1, b1, c1, invdij;
277 pbc_rvec_sub(pbc, xj, xi, xij);
278 pbc_rvec_sub(pbc, xk, xj, xjk);
281 invdij = gmx_invsqrt(iprod(xij, xij));
282 c1 = invdij * invdij * iprod(xij, xjk);
283 xp[XX] = xjk[XX] - c1*xij[XX];
284 xp[YY] = xjk[YY] - c1*xij[YY];
285 xp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
287 b1 = b*gmx_invsqrt(iprod(xp, xp));
290 x[XX] = xi[XX] + a1*xij[XX] + b1*xp[XX];
291 x[YY] = xi[YY] + a1*xij[YY] + b1*xp[YY];
292 x[ZZ] = xi[ZZ] + a1*xij[ZZ] + b1*xp[ZZ];
295 /* TOTAL: 63 flops */
298 static void constr_vsite3OUT(rvec xi, rvec xj, rvec xk, rvec x,
299 real a, real b, real c, t_pbc *pbc)
303 pbc_rvec_sub(pbc, xj, xi, xij);
304 pbc_rvec_sub(pbc, xk, xi, xik);
305 cprod(xij, xik, temp);
308 x[XX] = xi[XX] + a*xij[XX] + b*xik[XX] + c*temp[XX];
309 x[YY] = xi[YY] + a*xij[YY] + b*xik[YY] + c*temp[YY];
310 x[ZZ] = xi[ZZ] + a*xij[ZZ] + b*xik[ZZ] + c*temp[ZZ];
313 /* TOTAL: 33 flops */
316 static void constr_vsite4FD(rvec xi, rvec xj, rvec xk, rvec xl, rvec x,
317 real a, real b, real c, t_pbc *pbc)
319 rvec xij, xjk, xjl, temp;
322 pbc_rvec_sub(pbc, xj, xi, xij);
323 pbc_rvec_sub(pbc, xk, xj, xjk);
324 pbc_rvec_sub(pbc, xl, xj, xjl);
327 /* temp goes from i to a point on the plane jkl */
328 temp[XX] = xij[XX] + a*xjk[XX] + b*xjl[XX];
329 temp[YY] = xij[YY] + a*xjk[YY] + b*xjl[YY];
330 temp[ZZ] = xij[ZZ] + a*xjk[ZZ] + b*xjl[ZZ];
333 d = c*gmx_invsqrt(iprod(temp, temp));
336 x[XX] = xi[XX] + d*temp[XX];
337 x[YY] = xi[YY] + d*temp[YY];
338 x[ZZ] = xi[ZZ] + d*temp[ZZ];
341 /* TOTAL: 43 flops */
345 static void constr_vsite4FDN(rvec xi, rvec xj, rvec xk, rvec xl, rvec x,
346 real a, real b, real c, t_pbc *pbc)
348 rvec xij, xik, xil, ra, rb, rja, rjb, rm;
351 pbc_rvec_sub(pbc, xj, xi, xij);
352 pbc_rvec_sub(pbc, xk, xi, xik);
353 pbc_rvec_sub(pbc, xl, xi, xil);
366 rvec_sub(ra, xij, rja);
367 rvec_sub(rb, xij, rjb);
373 d = c*gmx_invsqrt(norm2(rm));
376 x[XX] = xi[XX] + d*rm[XX];
377 x[YY] = xi[YY] + d*rm[YY];
378 x[ZZ] = xi[ZZ] + d*rm[ZZ];
381 /* TOTAL: 47 flops */
385 static int constr_vsiten(t_iatom *ia, t_iparams ip[],
393 n3 = 3*ip[ia[0]].vsiten.n;
396 copy_rvec(x[ai], x1);
398 for (i = 3; i < n3; i += 3)
401 a = ip[ia[i]].vsiten.a;
404 pbc_dx_aiuc(pbc, x[ai], x1, dx);
408 rvec_sub(x[ai], x1, dx);
410 dsum[XX] += a*dx[XX];
411 dsum[YY] += a*dx[YY];
412 dsum[ZZ] += a*dx[ZZ];
416 x[av][XX] = x1[XX] + dsum[XX];
417 x[av][YY] = x1[YY] + dsum[YY];
418 x[av][ZZ] = x1[ZZ] + dsum[ZZ];
424 void construct_vsites_thread(gmx_vsite_t *vsite,
427 t_iparams ip[], t_ilist ilist[],
431 rvec xpbc, xv, vv, dx;
432 real a1, b1, c1, inv_dt;
433 int i, inc, ii, nra, nr, tp, ftype;
434 t_iatom avsite, ai, aj, ak, al, pbc_atom;
437 int *vsite_pbc, ishift;
438 rvec reftmp, vtmp, rtmp;
449 bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
453 for (ftype = 0; (ftype < F_NRE); ftype++)
455 if ((interaction_function[ftype].flags & IF_VSITE) &&
458 nra = interaction_function[ftype].nratoms;
460 nr = ilist[ftype].nr;
461 ia = ilist[ftype].iatoms;
465 pbc_null2 = pbc_null;
467 else if (pbc_null != NULL)
469 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
472 for (i = 0; i < nr; )
476 /* The vsite and constructing atoms */
480 /* Constants for constructing vsites */
482 /* Check what kind of pbc we need to use */
485 /* No charge groups, vsite follows its own pbc */
487 copy_rvec(x[avsite], xpbc);
489 else if (vsite_pbc != NULL)
491 pbc_atom = vsite_pbc[i/(1+nra)];
496 /* We need to copy the coordinates here,
497 * single for single atom cg's pbc_atom
498 * is the vsite itself.
500 copy_rvec(x[pbc_atom], xpbc);
502 pbc_null2 = pbc_null;
513 /* Copy the old position */
514 copy_rvec(x[avsite], xv);
516 /* Construct the vsite depending on type */
521 constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
527 constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
533 constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
539 constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
546 constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
554 constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
563 constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
567 inc = constr_vsiten(ia, ip, x, pbc_null2);
570 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
571 ftype, __FILE__, __LINE__);
576 /* Match the pbc of this vsite to the rest of its charge group */
577 ishift = pbc_dx_aiuc(pbc_null, x[avsite], xpbc, dx);
578 if (ishift != CENTRAL)
580 rvec_add(xpbc, dx, x[avsite]);
585 /* Calculate velocity of vsite... */
586 rvec_sub(x[avsite], xv, vv);
587 svmul(inv_dt, vv, v[avsite]);
590 /* Increment loop variables */
598 void construct_vsites(gmx_vsite_t *vsite,
601 t_iparams ip[], t_ilist ilist[],
602 int ePBC, gmx_bool bMolPBC, t_graph *graph,
603 t_commrec *cr, matrix box)
605 t_pbc pbc, *pbc_null;
609 bDomDec = cr && DOMAINDECOMP(cr);
611 /* We only need to do pbc when we have inter-cg vsites */
612 if (ePBC != epbcNONE && (bDomDec || bMolPBC) && vsite->n_intercg_vsite)
614 /* This is wasting some CPU time as we now do this multiple times
615 * per MD step. But how often do we have vsites with full pbc?
617 pbc_null = set_pbc_dd(&pbc, ePBC, cr != NULL ? cr->dd : NULL, FALSE, box);
628 dd_move_x_vsites(cr->dd, box, x);
630 else if (vsite->bPDvsitecomm)
632 /* I'm not sure whether the periodicity and shift are guaranteed
633 * to be consistent between different nodes when running e.g. polymers
634 * in parallel. In this special case we thus unshift/shift,
635 * but only when necessary. This is to make sure the coordinates
636 * we move don't end up a box away...
640 unshift_self(graph, box, x);
643 move_construct_x(vsite->vsitecomm, x, cr);
647 shift_self(graph, box, x);
652 if (vsite->nthreads == 1)
654 construct_vsites_thread(vsite,
661 #pragma omp parallel num_threads(vsite->nthreads)
663 construct_vsites_thread(vsite,
665 ip, vsite->tdata[gmx_omp_get_thread_num()].ilist,
668 /* Now we can construct the vsites that might depend on other vsites */
669 construct_vsites_thread(vsite,
671 ip, vsite->tdata[vsite->nthreads].ilist,
676 static void spread_vsite2(t_iatom ia[], real a,
677 rvec x[], rvec f[], rvec fshift[],
678 t_pbc *pbc, t_graph *g)
690 svmul(1-a, f[av], fi);
691 svmul( a, f[av], fj);
700 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
702 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
707 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
708 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
716 if (fshift && (siv != CENTRAL || sij != CENTRAL))
718 rvec_inc(fshift[siv], f[av]);
719 rvec_dec(fshift[CENTRAL], fi);
720 rvec_dec(fshift[sij], fj);
723 /* TOTAL: 13 flops */
726 void construct_vsites_mtop(gmx_vsite_t *vsite,
727 gmx_mtop_t *mtop, rvec x[])
730 gmx_molblock_t *molb;
734 for (mb = 0; mb < mtop->nmolblock; mb++)
736 molb = &mtop->molblock[mb];
737 molt = &mtop->moltype[molb->type];
738 for (mol = 0; mol < molb->nmol; mol++)
740 construct_vsites(vsite, x+as, 0.0, NULL,
741 mtop->ffparams.iparams, molt->ilist,
742 epbcNONE, TRUE, NULL, NULL, NULL);
743 as += molt->atoms.nr;
748 static void spread_vsite3(t_iatom ia[], real a, real b,
749 rvec x[], rvec f[], rvec fshift[],
750 t_pbc *pbc, t_graph *g)
753 atom_id av, ai, aj, ak;
762 svmul(1-a-b, f[av], fi);
763 svmul( a, f[av], fj);
764 svmul( b, f[av], fk);
774 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ia[1]), di);
776 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
778 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), di);
783 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
784 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
785 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
794 if (fshift && (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL))
796 rvec_inc(fshift[siv], f[av]);
797 rvec_dec(fshift[CENTRAL], fi);
798 rvec_dec(fshift[sij], fj);
799 rvec_dec(fshift[sik], fk);
802 /* TOTAL: 20 flops */
805 static void spread_vsite3FD(t_iatom ia[], real a, real b,
806 rvec x[], rvec f[], rvec fshift[],
807 gmx_bool VirCorr, matrix dxdf,
808 t_pbc *pbc, t_graph *g)
810 real fx, fy, fz, c, invl, fproj, a1;
811 rvec xvi, xij, xjk, xix, fv, temp;
812 t_iatom av, ai, aj, ak;
813 int svi, sji, skj, d;
820 copy_rvec(f[av], fv);
822 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
823 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
826 /* xix goes from i to point x on the line jk */
827 xix[XX] = xij[XX]+a*xjk[XX];
828 xix[YY] = xij[YY]+a*xjk[YY];
829 xix[ZZ] = xij[ZZ]+a*xjk[ZZ];
832 invl = gmx_invsqrt(iprod(xix, xix));
836 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
838 temp[XX] = c*(fv[XX]-fproj*xix[XX]);
839 temp[YY] = c*(fv[YY]-fproj*xix[YY]);
840 temp[ZZ] = c*(fv[ZZ]-fproj*xix[ZZ]);
843 /* c is already calculated in constr_vsite3FD
844 storing c somewhere will save 26 flops! */
847 f[ai][XX] += fv[XX] - temp[XX];
848 f[ai][YY] += fv[YY] - temp[YY];
849 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
850 f[aj][XX] += a1*temp[XX];
851 f[aj][YY] += a1*temp[YY];
852 f[aj][ZZ] += a1*temp[ZZ];
853 f[ak][XX] += a*temp[XX];
854 f[ak][YY] += a*temp[YY];
855 f[ak][ZZ] += a*temp[ZZ];
860 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
862 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
864 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
869 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
876 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
878 rvec_dec(fshift[svi], fv);
879 fshift[CENTRAL][XX] += fv[XX] - (1 + a)*temp[XX];
880 fshift[CENTRAL][YY] += fv[YY] - (1 + a)*temp[YY];
881 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a)*temp[ZZ];
882 fshift[ sji][XX] += temp[XX];
883 fshift[ sji][YY] += temp[YY];
884 fshift[ sji][ZZ] += temp[ZZ];
885 fshift[ skj][XX] += a*temp[XX];
886 fshift[ skj][YY] += a*temp[YY];
887 fshift[ skj][ZZ] += a*temp[ZZ];
892 /* When VirCorr=TRUE, the virial for the current forces is not
893 * calculated from the redistributed forces. This means that
894 * the effect of non-linear virtual site constructions on the virial
895 * needs to be added separately. This contribution can be calculated
896 * in many ways, but the simplest and cheapest way is to use
897 * the first constructing atom ai as a reference position in space:
898 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
903 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
905 for (i = 0; i < DIM; i++)
907 for (j = 0; j < DIM; j++)
909 /* As xix is a linear combination of j and k, use that here */
910 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
915 /* TOTAL: 61 flops */
918 static void spread_vsite3FAD(t_iatom ia[], real a, real b,
919 rvec x[], rvec f[], rvec fshift[],
920 gmx_bool VirCorr, matrix dxdf,
921 t_pbc *pbc, t_graph *g)
923 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
924 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
925 t_iatom av, ai, aj, ak;
926 int svi, sji, skj, d;
933 copy_rvec(f[ia[1]], fv);
935 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
936 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
939 invdij = gmx_invsqrt(iprod(xij, xij));
940 invdij2 = invdij * invdij;
941 c1 = iprod(xij, xjk) * invdij2;
942 xperp[XX] = xjk[XX] - c1*xij[XX];
943 xperp[YY] = xjk[YY] - c1*xij[YY];
944 xperp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
945 /* xperp in plane ijk, perp. to ij */
946 invdp = gmx_invsqrt(iprod(xperp, xperp));
951 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
952 storing them somewhere will save 45 flops! */
954 fproj = iprod(xij, fv)*invdij2;
955 svmul(fproj, xij, Fpij); /* proj. f on xij */
956 svmul(iprod(xperp, fv)*invdp*invdp, xperp, Fppp); /* proj. f on xperp */
957 svmul(b1*fproj, xperp, f3);
960 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
961 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
962 for (d = 0; (d < DIM); d++)
970 f[ai][XX] += fv[XX] - f1[XX] + c1*f2[XX] + f3[XX];
971 f[ai][YY] += fv[YY] - f1[YY] + c1*f2[YY] + f3[YY];
972 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1*f2[ZZ] + f3[ZZ];
973 f[aj][XX] += f1[XX] - c2*f2[XX] - f3[XX];
974 f[aj][YY] += f1[YY] - c2*f2[YY] - f3[YY];
975 f[aj][ZZ] += f1[ZZ] - c2*f2[ZZ] - f3[ZZ];
983 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
985 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
987 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
992 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
999 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
1001 rvec_dec(fshift[svi], fv);
1002 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1-c1)*f2[XX] + f3[XX];
1003 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1-c1)*f2[YY] + f3[YY];
1004 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1-c1)*f2[ZZ] + f3[ZZ];
1005 fshift[ sji][XX] += f1[XX] - c1 *f2[XX] - f3[XX];
1006 fshift[ sji][YY] += f1[YY] - c1 *f2[YY] - f3[YY];
1007 fshift[ sji][ZZ] += f1[ZZ] - c1 *f2[ZZ] - f3[ZZ];
1008 fshift[ skj][XX] += f2[XX];
1009 fshift[ skj][YY] += f2[YY];
1010 fshift[ skj][ZZ] += f2[ZZ];
1018 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1020 for (i = 0; i < DIM; i++)
1022 for (j = 0; j < DIM; j++)
1024 /* Note that xik=xij+xjk, so we have to add xij*f2 */
1027 + xij[i]*(f1[j] + (1 - c2)*f2[j] - f3[j])
1033 /* TOTAL: 113 flops */
1036 static void spread_vsite3OUT(t_iatom ia[], real a, real b, real c,
1037 rvec x[], rvec f[], rvec fshift[],
1038 gmx_bool VirCorr, matrix dxdf,
1039 t_pbc *pbc, t_graph *g)
1041 rvec xvi, xij, xik, fv, fj, fk;
1043 atom_id av, ai, aj, ak;
1052 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1053 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1056 copy_rvec(f[av], fv);
1063 fj[XX] = a*fv[XX] - xik[ZZ]*cfy + xik[YY]*cfz;
1064 fj[YY] = xik[ZZ]*cfx + a*fv[YY] - xik[XX]*cfz;
1065 fj[ZZ] = -xik[YY]*cfx + xik[XX]*cfy + a*fv[ZZ];
1067 fk[XX] = b*fv[XX] + xij[ZZ]*cfy - xij[YY]*cfz;
1068 fk[YY] = -xij[ZZ]*cfx + b*fv[YY] + xij[XX]*cfz;
1069 fk[ZZ] = xij[YY]*cfx - xij[XX]*cfy + b*fv[ZZ];
1072 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
1073 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
1074 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1075 rvec_inc(f[aj], fj);
1076 rvec_inc(f[ak], fk);
1081 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1083 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1085 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1090 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1097 if (fshift && (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL))
1099 rvec_dec(fshift[svi], fv);
1100 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
1101 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
1102 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1103 rvec_inc(fshift[sji], fj);
1104 rvec_inc(fshift[ski], fk);
1112 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1114 for (i = 0; i < DIM; i++)
1116 for (j = 0; j < DIM; j++)
1118 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j];
1123 /* TOTAL: 54 flops */
1126 static void spread_vsite4FD(t_iatom ia[], real a, real b, real c,
1127 rvec x[], rvec f[], rvec fshift[],
1128 gmx_bool VirCorr, matrix dxdf,
1129 t_pbc *pbc, t_graph *g)
1131 real d, invl, fproj, a1;
1132 rvec xvi, xij, xjk, xjl, xix, fv, temp;
1133 atom_id av, ai, aj, ak, al;
1135 int svi, sji, skj, slj, m;
1143 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1144 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1145 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1148 /* xix goes from i to point x on the plane jkl */
1149 for (m = 0; m < DIM; m++)
1151 xix[m] = xij[m] + a*xjk[m] + b*xjl[m];
1155 invl = gmx_invsqrt(iprod(xix, xix));
1157 /* 4 + ?10? flops */
1159 copy_rvec(f[av], fv);
1161 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
1163 for (m = 0; m < DIM; m++)
1165 temp[m] = d*(fv[m] - fproj*xix[m]);
1169 /* c is already calculated in constr_vsite3FD
1170 storing c somewhere will save 35 flops! */
1173 for (m = 0; m < DIM; m++)
1175 f[ai][m] += fv[m] - temp[m];
1176 f[aj][m] += a1*temp[m];
1177 f[ak][m] += a*temp[m];
1178 f[al][m] += b*temp[m];
1184 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1186 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1188 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
1190 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, aj), di);
1195 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1203 (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL))
1205 rvec_dec(fshift[svi], fv);
1206 for (m = 0; m < DIM; m++)
1208 fshift[CENTRAL][m] += fv[m] - (1 + a + b)*temp[m];
1209 fshift[ sji][m] += temp[m];
1210 fshift[ skj][m] += a*temp[m];
1211 fshift[ slj][m] += b*temp[m];
1220 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1222 for (i = 0; i < DIM; i++)
1224 for (j = 0; j < DIM; j++)
1226 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
1231 /* TOTAL: 77 flops */
1235 static void spread_vsite4FDN(t_iatom ia[], real a, real b, real c,
1236 rvec x[], rvec f[], rvec fshift[],
1237 gmx_bool VirCorr, matrix dxdf,
1238 t_pbc *pbc, t_graph *g)
1240 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1241 rvec fv, fj, fk, fl;
1245 int av, ai, aj, ak, al;
1246 int svi, sij, sik, sil;
1248 /* DEBUG: check atom indices */
1255 copy_rvec(f[av], fv);
1257 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1258 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1259 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1272 rvec_sub(ra, xij, rja);
1273 rvec_sub(rb, xij, rjb);
1274 rvec_sub(rb, ra, rab);
1277 cprod(rja, rjb, rm);
1280 invrm = gmx_invsqrt(norm2(rm));
1281 denom = invrm*invrm;
1284 cfx = c*invrm*fv[XX];
1285 cfy = c*invrm*fv[YY];
1286 cfz = c*invrm*fv[ZZ];
1297 fj[XX] = ( -rm[XX]*rt[XX]) * cfx + ( rab[ZZ]-rm[YY]*rt[XX]) * cfy + (-rab[YY]-rm[ZZ]*rt[XX]) * cfz;
1298 fj[YY] = (-rab[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( rab[XX]-rm[ZZ]*rt[YY]) * cfz;
1299 fj[ZZ] = ( rab[YY]-rm[XX]*rt[ZZ]) * cfx + (-rab[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1310 fk[XX] = ( -rm[XX]*rt[XX]) * cfx + (-a*rjb[ZZ]-rm[YY]*rt[XX]) * cfy + ( a*rjb[YY]-rm[ZZ]*rt[XX]) * cfz;
1311 fk[YY] = ( a*rjb[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + (-a*rjb[XX]-rm[ZZ]*rt[YY]) * cfz;
1312 fk[ZZ] = (-a*rjb[YY]-rm[XX]*rt[ZZ]) * cfx + ( a*rjb[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1323 fl[XX] = ( -rm[XX]*rt[XX]) * cfx + ( b*rja[ZZ]-rm[YY]*rt[XX]) * cfy + (-b*rja[YY]-rm[ZZ]*rt[XX]) * cfz;
1324 fl[YY] = (-b*rja[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( b*rja[XX]-rm[ZZ]*rt[YY]) * cfz;
1325 fl[ZZ] = ( b*rja[YY]-rm[XX]*rt[ZZ]) * cfx + (-b*rja[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1328 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1329 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1330 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1331 rvec_inc(f[aj], fj);
1332 rvec_inc(f[ak], fk);
1333 rvec_inc(f[al], fl);
1338 ivec_sub(SHIFT_IVEC(g, av), SHIFT_IVEC(g, ai), di);
1340 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1342 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1344 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, ai), di);
1349 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1356 if (fshift && (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL))
1358 rvec_dec(fshift[svi], fv);
1359 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1360 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1361 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1362 rvec_inc(fshift[sij], fj);
1363 rvec_inc(fshift[sik], fk);
1364 rvec_inc(fshift[sil], fl);
1372 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1374 for (i = 0; i < DIM; i++)
1376 for (j = 0; j < DIM; j++)
1378 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j] + xil[i]*fl[j];
1383 /* Total: 207 flops (Yuck!) */
1387 static int spread_vsiten(t_iatom ia[], t_iparams ip[],
1388 rvec x[], rvec f[], rvec fshift[],
1389 t_pbc *pbc, t_graph *g)
1397 n3 = 3*ip[ia[0]].vsiten.n;
1399 copy_rvec(x[av], xv);
1401 for (i = 0; i < n3; i += 3)
1406 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
1411 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1417 a = ip[ia[i]].vsiten.a;
1418 svmul(a, f[av], fi);
1419 rvec_inc(f[ai], fi);
1420 if (fshift && siv != CENTRAL)
1422 rvec_inc(fshift[siv], fi);
1423 rvec_dec(fshift[CENTRAL], fi);
1432 static int vsite_count(const t_ilist *ilist, int ftype)
1434 if (ftype == F_VSITEN)
1436 return ilist[ftype].nr/3;
1440 return ilist[ftype].nr/(1 + interaction_function[ftype].nratoms);
1444 static void spread_vsite_f_thread(gmx_vsite_t *vsite,
1445 rvec x[], rvec f[], rvec *fshift,
1446 gmx_bool VirCorr, matrix dxdf,
1447 t_iparams ip[], t_ilist ilist[],
1448 t_graph *g, t_pbc *pbc_null)
1452 int i, inc, m, nra, nr, tp, ftype;
1462 bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
1464 /* this loop goes backwards to be able to build *
1465 * higher type vsites from lower types */
1468 for (ftype = F_NRE-1; (ftype >= 0); ftype--)
1470 if ((interaction_function[ftype].flags & IF_VSITE) &&
1471 ilist[ftype].nr > 0)
1473 nra = interaction_function[ftype].nratoms;
1475 nr = ilist[ftype].nr;
1476 ia = ilist[ftype].iatoms;
1480 pbc_null2 = pbc_null;
1482 else if (pbc_null != NULL)
1484 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
1487 for (i = 0; i < nr; )
1489 if (vsite_pbc != NULL)
1491 if (vsite_pbc[i/(1+nra)] > -2)
1493 pbc_null2 = pbc_null;
1503 /* Constants for constructing */
1504 a1 = ip[tp].vsite.a;
1505 /* Construct the vsite depending on type */
1509 spread_vsite2(ia, a1, x, f, fshift, pbc_null2, g);
1512 b1 = ip[tp].vsite.b;
1513 spread_vsite3(ia, a1, b1, x, f, fshift, pbc_null2, g);
1516 b1 = ip[tp].vsite.b;
1517 spread_vsite3FD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1520 b1 = ip[tp].vsite.b;
1521 spread_vsite3FAD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1524 b1 = ip[tp].vsite.b;
1525 c1 = ip[tp].vsite.c;
1526 spread_vsite3OUT(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1529 b1 = ip[tp].vsite.b;
1530 c1 = ip[tp].vsite.c;
1531 spread_vsite4FD(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1534 b1 = ip[tp].vsite.b;
1535 c1 = ip[tp].vsite.c;
1536 spread_vsite4FDN(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1539 inc = spread_vsiten(ia, ip, x, f, fshift, pbc_null2, g);
1542 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
1543 ftype, __FILE__, __LINE__);
1545 clear_rvec(f[ia[1]]);
1547 /* Increment loop variables */
1555 void spread_vsite_f(gmx_vsite_t *vsite,
1556 rvec x[], rvec f[], rvec *fshift,
1557 gmx_bool VirCorr, matrix vir,
1558 t_nrnb *nrnb, t_idef *idef,
1559 int ePBC, gmx_bool bMolPBC, t_graph *g, matrix box,
1562 t_pbc pbc, *pbc_null;
1565 /* We only need to do pbc when we have inter-cg vsites */
1566 if ((DOMAINDECOMP(cr) || bMolPBC) && vsite->n_intercg_vsite)
1568 /* This is wasting some CPU time as we now do this multiple times
1569 * per MD step. But how often do we have vsites with full pbc?
1571 pbc_null = set_pbc_dd(&pbc, ePBC, cr->dd, FALSE, box);
1578 if (DOMAINDECOMP(cr))
1580 dd_clear_f_vsites(cr->dd, f);
1582 else if (PARTDECOMP(cr) && vsite->vsitecomm != NULL)
1584 pd_clear_nonlocal_constructs(vsite->vsitecomm, f);
1587 if (vsite->nthreads == 1)
1589 spread_vsite_f_thread(vsite,
1591 VirCorr, vsite->tdata[0].dxdf,
1592 idef->iparams, idef->il,
1597 /* First spread the vsites that might depend on other vsites */
1598 spread_vsite_f_thread(vsite,
1600 VirCorr, vsite->tdata[vsite->nthreads].dxdf,
1602 vsite->tdata[vsite->nthreads].ilist,
1605 #pragma omp parallel num_threads(vsite->nthreads)
1610 thread = gmx_omp_get_thread_num();
1612 if (thread == 0 || fshift == NULL)
1620 fshift_t = vsite->tdata[thread].fshift;
1622 for (i = 0; i < SHIFTS; i++)
1624 clear_rvec(fshift_t[i]);
1628 spread_vsite_f_thread(vsite,
1630 VirCorr, vsite->tdata[thread].dxdf,
1632 vsite->tdata[thread].ilist,
1640 for (th = 1; th < vsite->nthreads; th++)
1642 for (i = 0; i < SHIFTS; i++)
1644 rvec_inc(fshift[i], vsite->tdata[th].fshift[i]);
1654 for (th = 0; th < (vsite->nthreads == 1 ? 1 : vsite->nthreads+1); th++)
1656 for (i = 0; i < DIM; i++)
1658 for (j = 0; j < DIM; j++)
1660 vir[i][j] += -0.5*vsite->tdata[th].dxdf[i][j];
1666 if (DOMAINDECOMP(cr))
1668 dd_move_f_vsites(cr->dd, f, fshift);
1670 else if (vsite->bPDvsitecomm)
1672 /* We only move forces here, and they are independent of shifts */
1673 move_construct_f(vsite->vsitecomm, f, cr);
1676 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(idef->il, F_VSITE2));
1677 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(idef->il, F_VSITE3));
1678 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(idef->il, F_VSITE3FD));
1679 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(idef->il, F_VSITE3FAD));
1680 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(idef->il, F_VSITE3OUT));
1681 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(idef->il, F_VSITE4FD));
1682 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(idef->il, F_VSITE4FDN));
1683 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(idef->il, F_VSITEN));
1686 static int *atom2cg(t_block *cgs)
1690 snew(a2cg, cgs->index[cgs->nr]);
1691 for (cg = 0; cg < cgs->nr; cg++)
1693 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
1702 static int count_intercg_vsite(gmx_mtop_t *mtop,
1703 gmx_bool *bHaveChargeGroups)
1705 int mb, mt, ftype, nral, i, cg, a;
1706 gmx_molblock_t *molb;
1707 gmx_moltype_t *molt;
1711 int n_intercg_vsite;
1713 *bHaveChargeGroups = FALSE;
1715 n_intercg_vsite = 0;
1716 for (mb = 0; mb < mtop->nmolblock; mb++)
1718 molb = &mtop->molblock[mb];
1719 molt = &mtop->moltype[molb->type];
1721 if (molt->cgs.nr < molt->atoms.nr)
1723 *bHaveChargeGroups = TRUE;
1726 a2cg = atom2cg(&molt->cgs);
1727 for (ftype = 0; ftype < F_NRE; ftype++)
1729 if (interaction_function[ftype].flags & IF_VSITE)
1732 il = &molt->ilist[ftype];
1734 for (i = 0; i < il->nr; i += 1+nral)
1737 for (a = 1; a < nral; a++)
1739 if (a2cg[ia[1+a]] != cg)
1741 n_intercg_vsite += molb->nmol;
1751 return n_intercg_vsite;
1754 static int **get_vsite_pbc(t_iparams *iparams, t_ilist *ilist,
1755 t_atom *atom, t_mdatoms *md,
1756 t_block *cgs, int *a2cg)
1758 int ftype, nral, i, j, vsi, vsite, cg_v, cg_c, a, nc3 = 0;
1761 int **vsite_pbc, *vsite_pbc_f;
1763 gmx_bool bViteOnlyCG_and_FirstAtom;
1765 /* Make an array that tells if the pbc of an atom is set */
1766 snew(pbc_set, cgs->index[cgs->nr]);
1767 /* PBC is set for all non vsites */
1768 for (a = 0; a < cgs->index[cgs->nr]; a++)
1770 if ((atom && atom[a].ptype != eptVSite) ||
1771 (md && md->ptype[a] != eptVSite))
1777 snew(vsite_pbc, F_VSITEN-F_VSITE2+1);
1779 for (ftype = 0; ftype < F_NRE; ftype++)
1781 if (interaction_function[ftype].flags & IF_VSITE)
1787 snew(vsite_pbc[ftype-F_VSITE2], il->nr/(1+nral));
1788 vsite_pbc_f = vsite_pbc[ftype-F_VSITE2];
1796 /* A value of -2 signals that this vsite and its contructing
1797 * atoms are all within the same cg, so no pbc is required.
1799 vsite_pbc_f[vsi] = -2;
1800 /* Check if constructing atoms are outside the vsite's cg */
1802 if (ftype == F_VSITEN)
1804 nc3 = 3*iparams[ia[i]].vsiten.n;
1805 for (j = 0; j < nc3; j += 3)
1807 if (a2cg[ia[i+j+2]] != cg_v)
1809 vsite_pbc_f[vsi] = -1;
1815 for (a = 1; a < nral; a++)
1817 if (a2cg[ia[i+1+a]] != cg_v)
1819 vsite_pbc_f[vsi] = -1;
1823 if (vsite_pbc_f[vsi] == -1)
1825 /* Check if this is the first processed atom of a vsite only cg */
1826 bViteOnlyCG_and_FirstAtom = TRUE;
1827 for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1829 /* Non-vsites already have pbc set, so simply check for pbc_set */
1832 bViteOnlyCG_and_FirstAtom = FALSE;
1836 if (bViteOnlyCG_and_FirstAtom)
1838 /* First processed atom of a vsite only charge group.
1839 * The pbc of the input coordinates to construct_vsites
1840 * should be preserved.
1842 vsite_pbc_f[vsi] = vsite;
1844 else if (cg_v != a2cg[ia[1+i+1]])
1846 /* This vsite has a different charge group index
1847 * than it's first constructing atom
1848 * and the charge group has more than one atom,
1849 * search for the first normal particle
1850 * or vsite that already had its pbc defined.
1851 * If nothing is found, use full pbc for this vsite.
1853 for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1855 if (a != vsite && pbc_set[a])
1857 vsite_pbc_f[vsi] = a;
1860 fprintf(debug, "vsite %d match pbc with atom %d\n",
1868 fprintf(debug, "vsite atom %d cg %d - %d pbc atom %d\n",
1869 vsite+1, cgs->index[cg_v]+1, cgs->index[cg_v+1],
1870 vsite_pbc_f[vsi]+1);
1874 if (ftype == F_VSITEN)
1876 /* The other entries in vsite_pbc_f are not used for center vsites */
1884 /* This vsite now has its pbc defined */
1896 gmx_vsite_t *init_vsite(gmx_mtop_t *mtop, t_commrec *cr,
1897 gmx_bool bSerial_NoPBC)
1903 gmx_moltype_t *molt;
1906 /* check if there are vsites */
1908 for (i = 0; i < F_NRE; i++)
1910 if (interaction_function[i].flags & IF_VSITE)
1912 nvsite += gmx_mtop_ftype_count(mtop, i);
1923 vsite->n_intercg_vsite = count_intercg_vsite(mtop,
1924 &vsite->bHaveChargeGroups);
1926 /* If we don't have charge groups, the vsite follows its own pbc */
1927 if (!bSerial_NoPBC &&
1928 vsite->bHaveChargeGroups &&
1929 vsite->n_intercg_vsite > 0 && DOMAINDECOMP(cr))
1931 vsite->nvsite_pbc_molt = mtop->nmoltype;
1932 snew(vsite->vsite_pbc_molt, vsite->nvsite_pbc_molt);
1933 for (mt = 0; mt < mtop->nmoltype; mt++)
1935 molt = &mtop->moltype[mt];
1936 /* Make an atom to charge group index */
1937 a2cg = atom2cg(&molt->cgs);
1938 vsite->vsite_pbc_molt[mt] = get_vsite_pbc(mtop->ffparams.iparams,
1940 molt->atoms.atom, NULL,
1945 snew(vsite->vsite_pbc_loc_nalloc, F_VSITEN-F_VSITE2+1);
1946 snew(vsite->vsite_pbc_loc, F_VSITEN-F_VSITE2+1);
1951 vsite->nthreads = 1;
1955 vsite->nthreads = gmx_omp_nthreads_get(emntVSITE);
1959 /* We need one extra thread data structure for the overlap vsites */
1960 snew(vsite->tdata, vsite->nthreads+1);
1963 vsite->th_ind = NULL;
1964 vsite->th_ind_nalloc = 0;
1969 static void prepare_vsite_thread(const t_ilist *ilist,
1970 gmx_vsite_thread_t *vsite_th)
1974 for (ftype = 0; ftype < F_NRE; ftype++)
1976 if (interaction_function[ftype].flags & IF_VSITE)
1978 if (ilist[ftype].nr > vsite_th->ilist[ftype].nalloc)
1980 vsite_th->ilist[ftype].nalloc = over_alloc_large(ilist[ftype].nr);
1981 srenew(vsite_th->ilist[ftype].iatoms, vsite_th->ilist[ftype].nalloc);
1984 vsite_th->ilist[ftype].nr = 0;
1989 void split_vsites_over_threads(const t_ilist *ilist,
1990 const t_mdatoms *mdatoms,
1991 gmx_bool bLimitRange,
1995 int vsite_atom_range, natperthread;
2000 int nral1, inc, i, j;
2002 if (vsite->nthreads == 1)
2008 #pragma omp parallel for num_threads(vsite->nthreads) schedule(static)
2009 for (th = 0; th < vsite->nthreads; th++)
2011 prepare_vsite_thread(ilist, &vsite->tdata[th]);
2013 /* Master threads does the (potential) overlap vsites */
2014 prepare_vsite_thread(ilist, &vsite->tdata[vsite->nthreads]);
2016 /* The current way of distributing the vsites over threads in primitive.
2017 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
2018 * without taking into account how the vsites are distributed.
2019 * Without domain decomposition we bLimitRange=TRUE and we at least
2020 * tighten the upper bound of the range (useful for common systems
2021 * such as a vsite-protein in 3-site water).
2025 vsite_atom_range = -1;
2026 for (ftype = 0; ftype < F_NRE; ftype++)
2028 if ((interaction_function[ftype].flags & IF_VSITE) &&
2031 nral1 = 1 + NRAL(ftype);
2032 iat = ilist[ftype].iatoms;
2033 for (i = 0; i < ilist[ftype].nr; i += nral1)
2035 for (j = i+1; j < i+nral1; j++)
2037 vsite_atom_range = max(vsite_atom_range, iat[j]);
2046 vsite_atom_range = mdatoms->homenr;
2048 natperthread = (vsite_atom_range + vsite->nthreads - 1)/vsite->nthreads;
2052 fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n", mdatoms->nr, vsite_atom_range, natperthread);
2055 /* To simplify the vsite assignment, we make an index which tells us
2056 * to which thread particles, both non-vsites and vsites, are assigned.
2058 if (mdatoms->nr > vsite->th_ind_nalloc)
2060 vsite->th_ind_nalloc = over_alloc_large(mdatoms->nr);
2061 srenew(vsite->th_ind, vsite->th_ind_nalloc);
2063 th_ind = vsite->th_ind;
2065 for (i = 0; i < mdatoms->nr; i++)
2067 if (mdatoms->ptype[i] == eptVSite)
2069 /* vsites are not assigned to a thread yet */
2074 /* assign non-vsite particles to thread th */
2077 if (i == (th + 1)*natperthread && th < vsite->nthreads)
2083 for (ftype = 0; ftype < F_NRE; ftype++)
2085 if ((interaction_function[ftype].flags & IF_VSITE) &&
2088 nral1 = 1 + NRAL(ftype);
2090 iat = ilist[ftype].iatoms;
2091 for (i = 0; i < ilist[ftype].nr; )
2093 th = iat[1+i]/natperthread;
2094 /* We would like to assign this vsite the thread th,
2095 * but it might depend on atoms outside the atom range of th
2096 * or on another vsite not assigned to thread th.
2098 if (ftype != F_VSITEN)
2100 for (j = i+2; j < i+nral1; j++)
2102 if (th_ind[iat[j]] != th)
2104 /* Some constructing atoms are not assigned to
2105 * thread th, move this vsite to a separate batch.
2107 th = vsite->nthreads;
2114 for (j = i+2; j < i+inc; j += 3)
2116 if (th_ind[iat[j]] != th)
2118 th = vsite->nthreads;
2122 /* Copy this vsite to the thread data struct of thread th */
2123 il_th = &vsite->tdata[th].ilist[ftype];
2124 for (j = i; j < i+inc; j++)
2126 il_th->iatoms[il_th->nr++] = iat[j];
2128 /* Update this vsite's thread index entry */
2129 th_ind[iat[1+i]] = th;
2138 for (ftype = 0; ftype < F_NRE; ftype++)
2140 if ((interaction_function[ftype].flags & IF_VSITE) &&
2141 ilist[ftype].nr > 0)
2143 fprintf(debug, "%-20s thread dist:",
2144 interaction_function[ftype].longname);
2145 for (th = 0; th < vsite->nthreads+1; th++)
2147 fprintf(debug, " %4d", vsite->tdata[th].ilist[ftype].nr);
2149 fprintf(debug, "\n");
2155 void set_vsite_top(gmx_vsite_t *vsite, gmx_localtop_t *top, t_mdatoms *md,
2160 if (vsite->n_intercg_vsite > 0)
2162 if (vsite->bHaveChargeGroups)
2164 /* Make an atom to charge group index */
2165 a2cg = atom2cg(&top->cgs);
2166 vsite->vsite_pbc_loc = get_vsite_pbc(top->idef.iparams,
2167 top->idef.il, NULL, md,
2174 snew(vsite->vsitecomm, 1);
2175 vsite->bPDvsitecomm =
2176 setup_parallel_vsites(&(top->idef), cr, vsite->vsitecomm);
2180 if (vsite->nthreads > 1)
2182 if (vsite->bHaveChargeGroups || PARTDECOMP(cr))
2184 gmx_incons("Can not use threads virtual sites combined with charge groups or particle decomposition");
2187 split_vsites_over_threads(top->idef.il, md, !DOMAINDECOMP(cr), vsite);