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"
57 /* Routines to send/recieve coordinates and force
58 * of constructing atoms.
61 static void move_construct_x(t_comm_vsites *vsitecomm, rvec x[], t_commrec *cr)
67 sendbuf = vsitecomm->send_buf;
68 recvbuf = vsitecomm->recv_buf;
71 /* Prepare pulse left by copying to send buffer */
72 for (i = 0; i < vsitecomm->left_export_nconstruct; i++)
74 ia = vsitecomm->left_export_construct[i];
75 copy_rvec(x[ia], sendbuf[i]);
78 /* Pulse coordinates left */
79 gmx_tx_rx_real(cr, GMX_LEFT, (real *)sendbuf, 3*vsitecomm->left_export_nconstruct, GMX_RIGHT, (real *)recvbuf, 3*vsitecomm->right_import_nconstruct);
81 /* Copy from receive buffer to coordinate array */
82 for (i = 0; i < vsitecomm->right_import_nconstruct; i++)
84 ia = vsitecomm->right_import_construct[i];
85 copy_rvec(recvbuf[i], x[ia]);
88 /* Prepare pulse right by copying to send buffer */
89 for (i = 0; i < vsitecomm->right_export_nconstruct; i++)
91 ia = vsitecomm->right_export_construct[i];
92 copy_rvec(x[ia], sendbuf[i]);
95 /* Pulse coordinates right */
96 gmx_tx_rx_real(cr, GMX_RIGHT, (real *)sendbuf, 3*vsitecomm->right_export_nconstruct, GMX_LEFT, (real *)recvbuf, 3*vsitecomm->left_import_nconstruct);
98 /* Copy from receive buffer to coordinate array */
99 for (i = 0; i < vsitecomm->left_import_nconstruct; i++)
101 ia = vsitecomm->left_import_construct[i];
102 copy_rvec(recvbuf[i], x[ia]);
107 static void move_construct_f(t_comm_vsites *vsitecomm, rvec f[], t_commrec *cr)
113 sendbuf = vsitecomm->send_buf;
114 recvbuf = vsitecomm->recv_buf;
116 /* Prepare pulse right by copying to send buffer */
117 for (i = 0; i < vsitecomm->right_import_nconstruct; i++)
119 ia = vsitecomm->right_import_construct[i];
120 copy_rvec(f[ia], sendbuf[i]);
121 clear_rvec(f[ia]); /* Zero it here after moving, just to simplify debug book-keeping... */
124 /* Pulse forces right */
125 gmx_tx_rx_real(cr, GMX_RIGHT, (real *)sendbuf, 3*vsitecomm->right_import_nconstruct, GMX_LEFT, (real *)recvbuf, 3*vsitecomm->left_export_nconstruct);
127 /* Copy from receive buffer to coordinate array */
128 for (i = 0; i < vsitecomm->left_export_nconstruct; i++)
130 ia = vsitecomm->left_export_construct[i];
131 rvec_inc(f[ia], recvbuf[i]);
134 /* Prepare pulse left by copying to send buffer */
135 for (i = 0; i < vsitecomm->left_import_nconstruct; i++)
137 ia = vsitecomm->left_import_construct[i];
138 copy_rvec(f[ia], sendbuf[i]);
139 clear_rvec(f[ia]); /* Zero it here after moving, just to simplify debug book-keeping... */
142 /* Pulse coordinates left */
143 gmx_tx_rx_real(cr, GMX_LEFT, (real *)sendbuf, 3*vsitecomm->left_import_nconstruct, GMX_RIGHT, (real *)recvbuf, 3*vsitecomm->right_export_nconstruct);
145 /* Copy from receive buffer to coordinate array */
146 for (i = 0; i < vsitecomm->right_export_nconstruct; i++)
148 ia = vsitecomm->right_export_construct[i];
149 rvec_inc(f[ia], recvbuf[i]);
152 /* All forces are now on the home processors */
157 pd_clear_nonlocal_constructs(t_comm_vsites *vsitecomm, rvec f[])
161 for (i = 0; i < vsitecomm->left_import_nconstruct; i++)
163 ia = vsitecomm->left_import_construct[i];
166 for (i = 0; i < vsitecomm->right_import_nconstruct; i++)
168 ia = vsitecomm->right_import_construct[i];
175 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
179 return pbc_dx_aiuc(pbc, xi, xj, dx);
183 rvec_sub(xi, xj, dx);
188 /* Vsite construction routines */
190 static void constr_vsite2(rvec xi, rvec xj, rvec x, real a, t_pbc *pbc)
200 pbc_dx_aiuc(pbc, xj, xi, dx);
201 x[XX] = xi[XX] + a*dx[XX];
202 x[YY] = xi[YY] + a*dx[YY];
203 x[ZZ] = xi[ZZ] + a*dx[ZZ];
207 x[XX] = b*xi[XX] + a*xj[XX];
208 x[YY] = b*xi[YY] + a*xj[YY];
209 x[ZZ] = b*xi[ZZ] + a*xj[ZZ];
213 /* TOTAL: 10 flops */
216 static void constr_vsite3(rvec xi, rvec xj, rvec xk, rvec x, real a, real b,
227 pbc_dx_aiuc(pbc, xj, xi, dxj);
228 pbc_dx_aiuc(pbc, xk, xi, dxk);
229 x[XX] = xi[XX] + a*dxj[XX] + b*dxk[XX];
230 x[YY] = xi[YY] + a*dxj[YY] + b*dxk[YY];
231 x[ZZ] = xi[ZZ] + a*dxj[ZZ] + b*dxk[ZZ];
235 x[XX] = c*xi[XX] + a*xj[XX] + b*xk[XX];
236 x[YY] = c*xi[YY] + a*xj[YY] + b*xk[YY];
237 x[ZZ] = c*xi[ZZ] + a*xj[ZZ] + b*xk[ZZ];
241 /* TOTAL: 17 flops */
244 static void constr_vsite3FD(rvec xi, rvec xj, rvec xk, rvec x, real a, real b,
250 pbc_rvec_sub(pbc, xj, xi, xij);
251 pbc_rvec_sub(pbc, xk, xj, xjk);
254 /* temp goes from i to a point on the line jk */
255 temp[XX] = xij[XX] + a*xjk[XX];
256 temp[YY] = xij[YY] + a*xjk[YY];
257 temp[ZZ] = xij[ZZ] + a*xjk[ZZ];
260 c = b*gmx_invsqrt(iprod(temp, temp));
263 x[XX] = xi[XX] + c*temp[XX];
264 x[YY] = xi[YY] + c*temp[YY];
265 x[ZZ] = xi[ZZ] + c*temp[ZZ];
268 /* TOTAL: 34 flops */
271 static void constr_vsite3FAD(rvec xi, rvec xj, rvec xk, rvec x, real a, real b, t_pbc *pbc)
274 real a1, b1, c1, invdij;
276 pbc_rvec_sub(pbc, xj, xi, xij);
277 pbc_rvec_sub(pbc, xk, xj, xjk);
280 invdij = gmx_invsqrt(iprod(xij, xij));
281 c1 = invdij * invdij * iprod(xij, xjk);
282 xp[XX] = xjk[XX] - c1*xij[XX];
283 xp[YY] = xjk[YY] - c1*xij[YY];
284 xp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
286 b1 = b*gmx_invsqrt(iprod(xp, xp));
289 x[XX] = xi[XX] + a1*xij[XX] + b1*xp[XX];
290 x[YY] = xi[YY] + a1*xij[YY] + b1*xp[YY];
291 x[ZZ] = xi[ZZ] + a1*xij[ZZ] + b1*xp[ZZ];
294 /* TOTAL: 63 flops */
297 static void constr_vsite3OUT(rvec xi, rvec xj, rvec xk, rvec x,
298 real a, real b, real c, t_pbc *pbc)
302 pbc_rvec_sub(pbc, xj, xi, xij);
303 pbc_rvec_sub(pbc, xk, xi, xik);
304 cprod(xij, xik, temp);
307 x[XX] = xi[XX] + a*xij[XX] + b*xik[XX] + c*temp[XX];
308 x[YY] = xi[YY] + a*xij[YY] + b*xik[YY] + c*temp[YY];
309 x[ZZ] = xi[ZZ] + a*xij[ZZ] + b*xik[ZZ] + c*temp[ZZ];
312 /* TOTAL: 33 flops */
315 static void constr_vsite4FD(rvec xi, rvec xj, rvec xk, rvec xl, rvec x,
316 real a, real b, real c, t_pbc *pbc)
318 rvec xij, xjk, xjl, temp;
321 pbc_rvec_sub(pbc, xj, xi, xij);
322 pbc_rvec_sub(pbc, xk, xj, xjk);
323 pbc_rvec_sub(pbc, xl, xj, xjl);
326 /* temp goes from i to a point on the plane jkl */
327 temp[XX] = xij[XX] + a*xjk[XX] + b*xjl[XX];
328 temp[YY] = xij[YY] + a*xjk[YY] + b*xjl[YY];
329 temp[ZZ] = xij[ZZ] + a*xjk[ZZ] + b*xjl[ZZ];
332 d = c*gmx_invsqrt(iprod(temp, temp));
335 x[XX] = xi[XX] + d*temp[XX];
336 x[YY] = xi[YY] + d*temp[YY];
337 x[ZZ] = xi[ZZ] + d*temp[ZZ];
340 /* TOTAL: 43 flops */
344 static void constr_vsite4FDN(rvec xi, rvec xj, rvec xk, rvec xl, rvec x,
345 real a, real b, real c, t_pbc *pbc)
347 rvec xij, xik, xil, ra, rb, rja, rjb, rm;
350 pbc_rvec_sub(pbc, xj, xi, xij);
351 pbc_rvec_sub(pbc, xk, xi, xik);
352 pbc_rvec_sub(pbc, xl, xi, xil);
365 rvec_sub(ra, xij, rja);
366 rvec_sub(rb, xij, rjb);
372 d = c*gmx_invsqrt(norm2(rm));
375 x[XX] = xi[XX] + d*rm[XX];
376 x[YY] = xi[YY] + d*rm[YY];
377 x[ZZ] = xi[ZZ] + d*rm[ZZ];
380 /* TOTAL: 47 flops */
384 static int constr_vsiten(t_iatom *ia, t_iparams ip[],
392 n3 = 3*ip[ia[0]].vsiten.n;
395 copy_rvec(x[ai], x1);
397 for (i = 3; i < n3; i += 3)
400 a = ip[ia[i]].vsiten.a;
403 pbc_dx_aiuc(pbc, x[ai], x1, dx);
407 rvec_sub(x[ai], x1, dx);
409 dsum[XX] += a*dx[XX];
410 dsum[YY] += a*dx[YY];
411 dsum[ZZ] += a*dx[ZZ];
415 x[av][XX] = x1[XX] + dsum[XX];
416 x[av][YY] = x1[YY] + dsum[YY];
417 x[av][ZZ] = x1[ZZ] + dsum[ZZ];
423 void construct_vsites_thread(gmx_vsite_t *vsite,
424 rvec x[], t_nrnb *nrnb,
426 t_iparams ip[], t_ilist ilist[],
430 rvec xpbc, xv, vv, dx;
431 real a1, b1, c1, inv_dt;
432 int i, inc, ii, nra, nr, tp, ftype;
433 t_iatom avsite, ai, aj, ak, al, pbc_atom;
436 int *vsite_pbc, ishift;
437 rvec reftmp, vtmp, rtmp;
448 bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
452 for (ftype = 0; (ftype < F_NRE); ftype++)
454 if ((interaction_function[ftype].flags & IF_VSITE) &&
457 nra = interaction_function[ftype].nratoms;
459 nr = ilist[ftype].nr;
460 ia = ilist[ftype].iatoms;
464 pbc_null2 = pbc_null;
466 else if (pbc_null != NULL)
468 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
471 for (i = 0; i < nr; )
475 /* 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 */
520 constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
525 constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
530 constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
535 constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
541 constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
548 constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
556 constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
560 inc = constr_vsiten(ia, ip, x, pbc_null2);
563 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
564 ftype, __FILE__, __LINE__);
569 /* Match the pbc of this vsite to the rest of its charge group */
570 ishift = pbc_dx_aiuc(pbc_null, x[avsite], xpbc, dx);
571 if (ishift != CENTRAL)
573 rvec_add(xpbc, dx, x[avsite]);
578 /* Calculate velocity of vsite... */
579 rvec_sub(x[avsite], xv, vv);
580 svmul(inv_dt, vv, v[avsite]);
583 /* Increment loop variables */
591 void construct_vsites(FILE *log, gmx_vsite_t *vsite,
592 rvec x[], t_nrnb *nrnb,
594 t_iparams ip[], t_ilist ilist[],
595 int ePBC, gmx_bool bMolPBC, t_graph *graph,
596 t_commrec *cr, matrix box)
598 t_pbc pbc, *pbc_null;
602 bDomDec = cr && DOMAINDECOMP(cr);
604 /* We only need to do pbc when we have inter-cg vsites */
605 if (ePBC != epbcNONE && (bDomDec || bMolPBC) && vsite->n_intercg_vsite)
607 /* This is wasting some CPU time as we now do this multiple times
608 * per MD step. But how often do we have vsites with full pbc?
610 pbc_null = set_pbc_dd(&pbc, ePBC, cr != NULL ? cr->dd : NULL, FALSE, box);
621 dd_move_x_vsites(cr->dd, box, x);
623 else if (vsite->bPDvsitecomm)
625 /* I'm not sure whether the periodicity and shift are guaranteed
626 * to be consistent between different nodes when running e.g. polymers
627 * in parallel. In this special case we thus unshift/shift,
628 * but only when necessary. This is to make sure the coordinates
629 * we move don't end up a box away...
633 unshift_self(graph, box, x);
636 move_construct_x(vsite->vsitecomm, x, cr);
640 shift_self(graph, box, x);
645 if (vsite->nthreads == 1)
647 construct_vsites_thread(vsite,
654 #pragma omp parallel num_threads(vsite->nthreads)
656 construct_vsites_thread(vsite,
658 ip, vsite->tdata[gmx_omp_get_thread_num()].ilist,
661 /* Now we can construct the vsites that might depend on other vsites */
662 construct_vsites_thread(vsite,
664 ip, vsite->tdata[vsite->nthreads].ilist,
669 static void spread_vsite2(t_iatom ia[], real a,
670 rvec x[], rvec f[], rvec fshift[],
671 t_pbc *pbc, t_graph *g)
683 svmul(1-a, f[av], fi);
684 svmul( a, f[av], fj);
693 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
695 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
700 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
701 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
709 if (fshift && (siv != CENTRAL || sij != CENTRAL))
711 rvec_inc(fshift[siv], f[av]);
712 rvec_dec(fshift[CENTRAL], fi);
713 rvec_dec(fshift[sij], fj);
716 /* TOTAL: 13 flops */
719 void construct_vsites_mtop(FILE *log, gmx_vsite_t *vsite,
720 gmx_mtop_t *mtop, rvec x[])
723 gmx_molblock_t *molb;
727 for (mb = 0; mb < mtop->nmolblock; mb++)
729 molb = &mtop->molblock[mb];
730 molt = &mtop->moltype[molb->type];
731 for (mol = 0; mol < molb->nmol; mol++)
733 construct_vsites(log, vsite, x+as, NULL, 0.0, NULL,
734 mtop->ffparams.iparams, molt->ilist,
735 epbcNONE, TRUE, NULL, NULL, NULL);
736 as += molt->atoms.nr;
741 static void spread_vsite3(t_iatom ia[], real a, real b,
742 rvec x[], rvec f[], rvec fshift[],
743 t_pbc *pbc, t_graph *g)
746 atom_id av, ai, aj, ak;
755 svmul(1-a-b, f[av], fi);
756 svmul( a, f[av], fj);
757 svmul( b, f[av], fk);
767 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ia[1]), di);
769 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
771 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), di);
776 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
777 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
778 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
787 if (fshift && (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL))
789 rvec_inc(fshift[siv], f[av]);
790 rvec_dec(fshift[CENTRAL], fi);
791 rvec_dec(fshift[sij], fj);
792 rvec_dec(fshift[sik], fk);
795 /* TOTAL: 20 flops */
798 static void spread_vsite3FD(t_iatom ia[], real a, real b,
799 rvec x[], rvec f[], rvec fshift[],
800 gmx_bool VirCorr, matrix dxdf,
801 t_pbc *pbc, t_graph *g)
803 real fx, fy, fz, c, invl, fproj, a1;
804 rvec xvi, xij, xjk, xix, fv, temp;
805 t_iatom av, ai, aj, ak;
806 int svi, sji, skj, d;
813 copy_rvec(f[av], fv);
815 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
816 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
819 /* xix goes from i to point x on the line jk */
820 xix[XX] = xij[XX]+a*xjk[XX];
821 xix[YY] = xij[YY]+a*xjk[YY];
822 xix[ZZ] = xij[ZZ]+a*xjk[ZZ];
825 invl = gmx_invsqrt(iprod(xix, xix));
829 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
831 temp[XX] = c*(fv[XX]-fproj*xix[XX]);
832 temp[YY] = c*(fv[YY]-fproj*xix[YY]);
833 temp[ZZ] = c*(fv[ZZ]-fproj*xix[ZZ]);
836 /* c is already calculated in constr_vsite3FD
837 storing c somewhere will save 26 flops! */
840 f[ai][XX] += fv[XX] - temp[XX];
841 f[ai][YY] += fv[YY] - temp[YY];
842 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
843 f[aj][XX] += a1*temp[XX];
844 f[aj][YY] += a1*temp[YY];
845 f[aj][ZZ] += a1*temp[ZZ];
846 f[ak][XX] += a*temp[XX];
847 f[ak][YY] += a*temp[YY];
848 f[ak][ZZ] += a*temp[ZZ];
853 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
855 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
857 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
862 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
869 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
871 rvec_dec(fshift[svi], fv);
872 fshift[CENTRAL][XX] += fv[XX] - (1 + a)*temp[XX];
873 fshift[CENTRAL][YY] += fv[YY] - (1 + a)*temp[YY];
874 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a)*temp[ZZ];
875 fshift[ sji][XX] += temp[XX];
876 fshift[ sji][YY] += temp[YY];
877 fshift[ sji][ZZ] += temp[ZZ];
878 fshift[ skj][XX] += a*temp[XX];
879 fshift[ skj][YY] += a*temp[YY];
880 fshift[ skj][ZZ] += a*temp[ZZ];
885 /* When VirCorr=TRUE, the virial for the current forces is not
886 * calculated from the redistributed forces. This means that
887 * the effect of non-linear virtual site constructions on the virial
888 * needs to be added separately. This contribution can be calculated
889 * in many ways, but the simplest and cheapest way is to use
890 * the first constructing atom ai as a reference position in space:
891 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
896 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
898 for (i = 0; i < DIM; i++)
900 for (j = 0; j < DIM; j++)
902 /* As xix is a linear combination of j and k, use that here */
903 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
908 /* TOTAL: 61 flops */
911 static void spread_vsite3FAD(t_iatom ia[], real a, real b,
912 rvec x[], rvec f[], rvec fshift[],
913 gmx_bool VirCorr, matrix dxdf,
914 t_pbc *pbc, t_graph *g)
916 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
917 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
918 t_iatom av, ai, aj, ak;
919 int svi, sji, skj, d;
926 copy_rvec(f[ia[1]], fv);
928 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
929 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
932 invdij = gmx_invsqrt(iprod(xij, xij));
933 invdij2 = invdij * invdij;
934 c1 = iprod(xij, xjk) * invdij2;
935 xperp[XX] = xjk[XX] - c1*xij[XX];
936 xperp[YY] = xjk[YY] - c1*xij[YY];
937 xperp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
938 /* xperp in plane ijk, perp. to ij */
939 invdp = gmx_invsqrt(iprod(xperp, xperp));
944 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
945 storing them somewhere will save 45 flops! */
947 fproj = iprod(xij, fv)*invdij2;
948 svmul(fproj, xij, Fpij); /* proj. f on xij */
949 svmul(iprod(xperp, fv)*invdp*invdp, xperp, Fppp); /* proj. f on xperp */
950 svmul(b1*fproj, xperp, f3);
953 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
954 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
955 for (d = 0; (d < DIM); d++)
963 f[ai][XX] += fv[XX] - f1[XX] + c1*f2[XX] + f3[XX];
964 f[ai][YY] += fv[YY] - f1[YY] + c1*f2[YY] + f3[YY];
965 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1*f2[ZZ] + f3[ZZ];
966 f[aj][XX] += f1[XX] - c2*f2[XX] - f3[XX];
967 f[aj][YY] += f1[YY] - c2*f2[YY] - f3[YY];
968 f[aj][ZZ] += f1[ZZ] - c2*f2[ZZ] - f3[ZZ];
976 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
978 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
980 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
985 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
992 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
994 rvec_dec(fshift[svi], fv);
995 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1-c1)*f2[XX] + f3[XX];
996 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1-c1)*f2[YY] + f3[YY];
997 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1-c1)*f2[ZZ] + f3[ZZ];
998 fshift[ sji][XX] += f1[XX] - c1 *f2[XX] - f3[XX];
999 fshift[ sji][YY] += f1[YY] - c1 *f2[YY] - f3[YY];
1000 fshift[ sji][ZZ] += f1[ZZ] - c1 *f2[ZZ] - f3[ZZ];
1001 fshift[ skj][XX] += f2[XX];
1002 fshift[ skj][YY] += f2[YY];
1003 fshift[ skj][ZZ] += f2[ZZ];
1011 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1013 for (i = 0; i < DIM; i++)
1015 for (j = 0; j < DIM; j++)
1017 /* Note that xik=xij+xjk, so we have to add xij*f2 */
1020 + xij[i]*(f1[j] + (1 - c2)*f2[j] - f3[j])
1026 /* TOTAL: 113 flops */
1029 static void spread_vsite3OUT(t_iatom ia[], real a, real b, real c,
1030 rvec x[], rvec f[], rvec fshift[],
1031 gmx_bool VirCorr, matrix dxdf,
1032 t_pbc *pbc, t_graph *g)
1034 rvec xvi, xij, xik, fv, fj, fk;
1036 atom_id av, ai, aj, ak;
1045 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1046 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1049 copy_rvec(f[av], fv);
1056 fj[XX] = a*fv[XX] - xik[ZZ]*cfy + xik[YY]*cfz;
1057 fj[YY] = xik[ZZ]*cfx + a*fv[YY] - xik[XX]*cfz;
1058 fj[ZZ] = -xik[YY]*cfx + xik[XX]*cfy + a*fv[ZZ];
1060 fk[XX] = b*fv[XX] + xij[ZZ]*cfy - xij[YY]*cfz;
1061 fk[YY] = -xij[ZZ]*cfx + b*fv[YY] + xij[XX]*cfz;
1062 fk[ZZ] = xij[YY]*cfx - xij[XX]*cfy + b*fv[ZZ];
1065 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
1066 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
1067 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1068 rvec_inc(f[aj], fj);
1069 rvec_inc(f[ak], fk);
1074 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1076 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1078 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1083 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1090 if (fshift && (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL))
1092 rvec_dec(fshift[svi], fv);
1093 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
1094 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
1095 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1096 rvec_inc(fshift[sji], fj);
1097 rvec_inc(fshift[ski], fk);
1105 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1107 for (i = 0; i < DIM; i++)
1109 for (j = 0; j < DIM; j++)
1111 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j];
1116 /* TOTAL: 54 flops */
1119 static void spread_vsite4FD(t_iatom ia[], real a, real b, real c,
1120 rvec x[], rvec f[], rvec fshift[],
1121 gmx_bool VirCorr, matrix dxdf,
1122 t_pbc *pbc, t_graph *g)
1124 real d, invl, fproj, a1;
1125 rvec xvi, xij, xjk, xjl, xix, fv, temp;
1126 atom_id av, ai, aj, ak, al;
1128 int svi, sji, skj, slj, m;
1136 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1137 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1138 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1141 /* xix goes from i to point x on the plane jkl */
1142 for (m = 0; m < DIM; m++)
1144 xix[m] = xij[m] + a*xjk[m] + b*xjl[m];
1148 invl = gmx_invsqrt(iprod(xix, xix));
1150 /* 4 + ?10? flops */
1152 copy_rvec(f[av], fv);
1154 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
1156 for (m = 0; m < DIM; m++)
1158 temp[m] = d*(fv[m] - fproj*xix[m]);
1162 /* c is already calculated in constr_vsite3FD
1163 storing c somewhere will save 35 flops! */
1166 for (m = 0; m < DIM; m++)
1168 f[ai][m] += fv[m] - temp[m];
1169 f[aj][m] += a1*temp[m];
1170 f[ak][m] += a*temp[m];
1171 f[al][m] += b*temp[m];
1177 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1179 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1181 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
1183 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, aj), di);
1188 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1196 (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL))
1198 rvec_dec(fshift[svi], fv);
1199 for (m = 0; m < DIM; m++)
1201 fshift[CENTRAL][m] += fv[m] - (1 + a + b)*temp[m];
1202 fshift[ sji][m] += temp[m];
1203 fshift[ skj][m] += a*temp[m];
1204 fshift[ slj][m] += b*temp[m];
1213 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1215 for (i = 0; i < DIM; i++)
1217 for (j = 0; j < DIM; j++)
1219 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
1224 /* TOTAL: 77 flops */
1228 static void spread_vsite4FDN(t_iatom ia[], real a, real b, real c,
1229 rvec x[], rvec f[], rvec fshift[],
1230 gmx_bool VirCorr, matrix dxdf,
1231 t_pbc *pbc, t_graph *g)
1233 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1234 rvec fv, fj, fk, fl;
1238 int av, ai, aj, ak, al;
1239 int svi, sij, sik, sil;
1241 /* DEBUG: check atom indices */
1248 copy_rvec(f[av], fv);
1250 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1251 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1252 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1265 rvec_sub(ra, xij, rja);
1266 rvec_sub(rb, xij, rjb);
1267 rvec_sub(rb, ra, rab);
1270 cprod(rja, rjb, rm);
1273 invrm = gmx_invsqrt(norm2(rm));
1274 denom = invrm*invrm;
1277 cfx = c*invrm*fv[XX];
1278 cfy = c*invrm*fv[YY];
1279 cfz = c*invrm*fv[ZZ];
1290 fj[XX] = ( -rm[XX]*rt[XX]) * cfx + ( rab[ZZ]-rm[YY]*rt[XX]) * cfy + (-rab[YY]-rm[ZZ]*rt[XX]) * cfz;
1291 fj[YY] = (-rab[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( rab[XX]-rm[ZZ]*rt[YY]) * cfz;
1292 fj[ZZ] = ( rab[YY]-rm[XX]*rt[ZZ]) * cfx + (-rab[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1303 fk[XX] = ( -rm[XX]*rt[XX]) * cfx + (-a*rjb[ZZ]-rm[YY]*rt[XX]) * cfy + ( a*rjb[YY]-rm[ZZ]*rt[XX]) * cfz;
1304 fk[YY] = ( a*rjb[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + (-a*rjb[XX]-rm[ZZ]*rt[YY]) * cfz;
1305 fk[ZZ] = (-a*rjb[YY]-rm[XX]*rt[ZZ]) * cfx + ( a*rjb[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1316 fl[XX] = ( -rm[XX]*rt[XX]) * cfx + ( b*rja[ZZ]-rm[YY]*rt[XX]) * cfy + (-b*rja[YY]-rm[ZZ]*rt[XX]) * cfz;
1317 fl[YY] = (-b*rja[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( b*rja[XX]-rm[ZZ]*rt[YY]) * cfz;
1318 fl[ZZ] = ( b*rja[YY]-rm[XX]*rt[ZZ]) * cfx + (-b*rja[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1321 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1322 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1323 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1324 rvec_inc(f[aj], fj);
1325 rvec_inc(f[ak], fk);
1326 rvec_inc(f[al], fl);
1331 ivec_sub(SHIFT_IVEC(g, av), SHIFT_IVEC(g, ai), di);
1333 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1335 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1337 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, ai), di);
1342 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1349 if (fshift && (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL))
1351 rvec_dec(fshift[svi], fv);
1352 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1353 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1354 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1355 rvec_inc(fshift[sij], fj);
1356 rvec_inc(fshift[sik], fk);
1357 rvec_inc(fshift[sil], fl);
1365 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1367 for (i = 0; i < DIM; i++)
1369 for (j = 0; j < DIM; j++)
1371 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j] + xil[i]*fl[j];
1376 /* Total: 207 flops (Yuck!) */
1380 static int spread_vsiten(t_iatom ia[], t_iparams ip[],
1381 rvec x[], rvec f[], rvec fshift[],
1382 t_pbc *pbc, t_graph *g)
1390 n3 = 3*ip[ia[0]].vsiten.n;
1392 copy_rvec(x[av], xv);
1394 for (i = 0; i < n3; i += 3)
1399 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
1404 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1410 a = ip[ia[i]].vsiten.a;
1411 svmul(a, f[av], fi);
1412 rvec_inc(f[ai], fi);
1413 if (fshift && siv != CENTRAL)
1415 rvec_inc(fshift[siv], fi);
1416 rvec_dec(fshift[CENTRAL], fi);
1425 static int vsite_count(const t_ilist *ilist, int ftype)
1427 if (ftype == F_VSITEN)
1429 return ilist[ftype].nr/3;
1433 return ilist[ftype].nr/(1 + interaction_function[ftype].nratoms);
1437 static void spread_vsite_f_thread(gmx_vsite_t *vsite,
1438 rvec x[], rvec f[], rvec *fshift,
1439 gmx_bool VirCorr, matrix dxdf,
1440 t_iparams ip[], t_ilist ilist[],
1441 t_graph *g, t_pbc *pbc_null)
1445 int i, inc, m, nra, nr, tp, ftype;
1455 bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
1457 /* this loop goes backwards to be able to build *
1458 * higher type vsites from lower types */
1461 for (ftype = F_NRE-1; (ftype >= 0); ftype--)
1463 if ((interaction_function[ftype].flags & IF_VSITE) &&
1464 ilist[ftype].nr > 0)
1466 nra = interaction_function[ftype].nratoms;
1468 nr = ilist[ftype].nr;
1469 ia = ilist[ftype].iatoms;
1473 pbc_null2 = pbc_null;
1475 else if (pbc_null != NULL)
1477 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
1480 for (i = 0; i < nr; )
1482 if (vsite_pbc != NULL)
1484 if (vsite_pbc[i/(1+nra)] > -2)
1486 pbc_null2 = pbc_null;
1496 /* Constants for constructing */
1497 a1 = ip[tp].vsite.a;
1498 /* Construct the vsite depending on type */
1502 spread_vsite2(ia, a1, x, f, fshift, pbc_null2, g);
1505 b1 = ip[tp].vsite.b;
1506 spread_vsite3(ia, a1, b1, x, f, fshift, pbc_null2, g);
1509 b1 = ip[tp].vsite.b;
1510 spread_vsite3FD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1513 b1 = ip[tp].vsite.b;
1514 spread_vsite3FAD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1517 b1 = ip[tp].vsite.b;
1518 c1 = ip[tp].vsite.c;
1519 spread_vsite3OUT(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1522 b1 = ip[tp].vsite.b;
1523 c1 = ip[tp].vsite.c;
1524 spread_vsite4FD(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1527 b1 = ip[tp].vsite.b;
1528 c1 = ip[tp].vsite.c;
1529 spread_vsite4FDN(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1532 inc = spread_vsiten(ia, ip, x, f, fshift, pbc_null2, g);
1535 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
1536 ftype, __FILE__, __LINE__);
1538 clear_rvec(f[ia[1]]);
1540 /* Increment loop variables */
1548 void spread_vsite_f(FILE *log, gmx_vsite_t *vsite,
1549 rvec x[], rvec f[], rvec *fshift,
1550 gmx_bool VirCorr, matrix vir,
1551 t_nrnb *nrnb, t_idef *idef,
1552 int ePBC, gmx_bool bMolPBC, t_graph *g, matrix box,
1555 t_pbc pbc, *pbc_null;
1558 /* We only need to do pbc when we have inter-cg vsites */
1559 if ((DOMAINDECOMP(cr) || bMolPBC) && vsite->n_intercg_vsite)
1561 /* This is wasting some CPU time as we now do this multiple times
1562 * per MD step. But how often do we have vsites with full pbc?
1564 pbc_null = set_pbc_dd(&pbc, ePBC, cr->dd, FALSE, box);
1571 if (DOMAINDECOMP(cr))
1573 dd_clear_f_vsites(cr->dd, f);
1575 else if (PARTDECOMP(cr) && vsite->vsitecomm != NULL)
1577 pd_clear_nonlocal_constructs(vsite->vsitecomm, f);
1580 if (vsite->nthreads == 1)
1582 spread_vsite_f_thread(vsite,
1584 VirCorr, vsite->tdata[0].dxdf,
1585 idef->iparams, idef->il,
1590 /* First spread the vsites that might depend on other vsites */
1591 spread_vsite_f_thread(vsite,
1593 VirCorr, vsite->tdata[vsite->nthreads].dxdf,
1595 vsite->tdata[vsite->nthreads].ilist,
1598 #pragma omp parallel num_threads(vsite->nthreads)
1603 thread = gmx_omp_get_thread_num();
1605 if (thread == 0 || fshift == NULL)
1613 fshift_t = vsite->tdata[thread].fshift;
1615 for (i = 0; i < SHIFTS; i++)
1617 clear_rvec(fshift_t[i]);
1621 spread_vsite_f_thread(vsite,
1623 VirCorr, vsite->tdata[thread].dxdf,
1625 vsite->tdata[thread].ilist,
1633 for (th = 1; th < vsite->nthreads; th++)
1635 for (i = 0; i < SHIFTS; i++)
1637 rvec_inc(fshift[i], vsite->tdata[th].fshift[i]);
1647 for (th = 0; th < (vsite->nthreads == 1 ? 1 : vsite->nthreads+1); th++)
1649 for (i = 0; i < DIM; i++)
1651 for (j = 0; j < DIM; j++)
1653 vir[i][j] += -0.5*vsite->tdata[th].dxdf[i][j];
1659 if (DOMAINDECOMP(cr))
1661 dd_move_f_vsites(cr->dd, f, fshift);
1663 else if (vsite->bPDvsitecomm)
1665 /* We only move forces here, and they are independent of shifts */
1666 move_construct_f(vsite->vsitecomm, f, cr);
1669 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(idef->il, F_VSITE2));
1670 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(idef->il, F_VSITE3));
1671 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(idef->il, F_VSITE3FD));
1672 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(idef->il, F_VSITE3FAD));
1673 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(idef->il, F_VSITE3OUT));
1674 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(idef->il, F_VSITE4FD));
1675 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(idef->il, F_VSITE4FDN));
1676 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(idef->il, F_VSITEN));
1679 static int *atom2cg(t_block *cgs)
1683 snew(a2cg, cgs->index[cgs->nr]);
1684 for (cg = 0; cg < cgs->nr; cg++)
1686 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
1695 static int count_intercg_vsite(gmx_mtop_t *mtop,
1696 gmx_bool *bHaveChargeGroups)
1698 int mb, mt, ftype, nral, i, cg, a;
1699 gmx_molblock_t *molb;
1700 gmx_moltype_t *molt;
1704 int n_intercg_vsite;
1706 *bHaveChargeGroups = FALSE;
1708 n_intercg_vsite = 0;
1709 for (mb = 0; mb < mtop->nmolblock; mb++)
1711 molb = &mtop->molblock[mb];
1712 molt = &mtop->moltype[molb->type];
1714 if (molt->cgs.nr < molt->atoms.nr)
1716 *bHaveChargeGroups = TRUE;
1719 a2cg = atom2cg(&molt->cgs);
1720 for (ftype = 0; ftype < F_NRE; ftype++)
1722 if (interaction_function[ftype].flags & IF_VSITE)
1725 il = &molt->ilist[ftype];
1727 for (i = 0; i < il->nr; i += 1+nral)
1730 for (a = 1; a < nral; a++)
1732 if (a2cg[ia[1+a]] != cg)
1734 n_intercg_vsite += molb->nmol;
1744 return n_intercg_vsite;
1747 static int **get_vsite_pbc(t_iparams *iparams, t_ilist *ilist,
1748 t_atom *atom, t_mdatoms *md,
1749 t_block *cgs, int *a2cg)
1751 int ftype, nral, i, j, vsi, vsite, cg_v, cg_c, a, nc3 = 0;
1754 int **vsite_pbc, *vsite_pbc_f;
1756 gmx_bool bViteOnlyCG_and_FirstAtom;
1758 /* Make an array that tells if the pbc of an atom is set */
1759 snew(pbc_set, cgs->index[cgs->nr]);
1760 /* PBC is set for all non vsites */
1761 for (a = 0; a < cgs->index[cgs->nr]; a++)
1763 if ((atom && atom[a].ptype != eptVSite) ||
1764 (md && md->ptype[a] != eptVSite))
1770 snew(vsite_pbc, F_VSITEN-F_VSITE2+1);
1772 for (ftype = 0; ftype < F_NRE; ftype++)
1774 if (interaction_function[ftype].flags & IF_VSITE)
1780 snew(vsite_pbc[ftype-F_VSITE2], il->nr/(1+nral));
1781 vsite_pbc_f = vsite_pbc[ftype-F_VSITE2];
1789 /* A value of -2 signals that this vsite and its contructing
1790 * atoms are all within the same cg, so no pbc is required.
1792 vsite_pbc_f[vsi] = -2;
1793 /* Check if constructing atoms are outside the vsite's cg */
1795 if (ftype == F_VSITEN)
1797 nc3 = 3*iparams[ia[i]].vsiten.n;
1798 for (j = 0; j < nc3; j += 3)
1800 if (a2cg[ia[i+j+2]] != cg_v)
1802 vsite_pbc_f[vsi] = -1;
1808 for (a = 1; a < nral; a++)
1810 if (a2cg[ia[i+1+a]] != cg_v)
1812 vsite_pbc_f[vsi] = -1;
1816 if (vsite_pbc_f[vsi] == -1)
1818 /* Check if this is the first processed atom of a vsite only cg */
1819 bViteOnlyCG_and_FirstAtom = TRUE;
1820 for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1822 /* Non-vsites already have pbc set, so simply check for pbc_set */
1825 bViteOnlyCG_and_FirstAtom = FALSE;
1829 if (bViteOnlyCG_and_FirstAtom)
1831 /* First processed atom of a vsite only charge group.
1832 * The pbc of the input coordinates to construct_vsites
1833 * should be preserved.
1835 vsite_pbc_f[vsi] = vsite;
1837 else if (cg_v != a2cg[ia[1+i+1]])
1839 /* This vsite has a different charge group index
1840 * than it's first constructing atom
1841 * and the charge group has more than one atom,
1842 * search for the first normal particle
1843 * or vsite that already had its pbc defined.
1844 * If nothing is found, use full pbc for this vsite.
1846 for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1848 if (a != vsite && pbc_set[a])
1850 vsite_pbc_f[vsi] = a;
1853 fprintf(debug, "vsite %d match pbc with atom %d\n",
1861 fprintf(debug, "vsite atom %d cg %d - %d pbc atom %d\n",
1862 vsite+1, cgs->index[cg_v]+1, cgs->index[cg_v+1],
1863 vsite_pbc_f[vsi]+1);
1867 if (ftype == F_VSITEN)
1869 /* The other entries in vsite_pbc_f are not used for center vsites */
1877 /* This vsite now has its pbc defined */
1889 gmx_vsite_t *init_vsite(gmx_mtop_t *mtop, t_commrec *cr,
1890 gmx_bool bSerial_NoPBC)
1896 gmx_moltype_t *molt;
1899 /* check if there are vsites */
1901 for (i = 0; i < F_NRE; i++)
1903 if (interaction_function[i].flags & IF_VSITE)
1905 nvsite += gmx_mtop_ftype_count(mtop, i);
1916 vsite->n_intercg_vsite = count_intercg_vsite(mtop,
1917 &vsite->bHaveChargeGroups);
1919 /* If we don't have charge groups, the vsite follows its own pbc */
1920 if (!bSerial_NoPBC &&
1921 vsite->bHaveChargeGroups &&
1922 vsite->n_intercg_vsite > 0 && DOMAINDECOMP(cr))
1924 vsite->nvsite_pbc_molt = mtop->nmoltype;
1925 snew(vsite->vsite_pbc_molt, vsite->nvsite_pbc_molt);
1926 for (mt = 0; mt < mtop->nmoltype; mt++)
1928 molt = &mtop->moltype[mt];
1929 /* Make an atom to charge group index */
1930 a2cg = atom2cg(&molt->cgs);
1931 vsite->vsite_pbc_molt[mt] = get_vsite_pbc(mtop->ffparams.iparams,
1933 molt->atoms.atom, NULL,
1938 snew(vsite->vsite_pbc_loc_nalloc, F_VSITEN-F_VSITE2+1);
1939 snew(vsite->vsite_pbc_loc, F_VSITEN-F_VSITE2+1);
1944 vsite->nthreads = 1;
1948 vsite->nthreads = gmx_omp_nthreads_get(emntVSITE);
1952 /* We need one extra thread data structure for the overlap vsites */
1953 snew(vsite->tdata, vsite->nthreads+1);
1956 vsite->th_ind = NULL;
1957 vsite->th_ind_nalloc = 0;
1962 static void prepare_vsite_thread(const t_ilist *ilist,
1963 gmx_vsite_thread_t *vsite_th)
1967 for (ftype = 0; ftype < F_NRE; ftype++)
1969 if (interaction_function[ftype].flags & IF_VSITE)
1971 if (ilist[ftype].nr > vsite_th->ilist[ftype].nalloc)
1973 vsite_th->ilist[ftype].nalloc = over_alloc_large(ilist[ftype].nr);
1974 srenew(vsite_th->ilist[ftype].iatoms, vsite_th->ilist[ftype].nalloc);
1977 vsite_th->ilist[ftype].nr = 0;
1982 void split_vsites_over_threads(const t_ilist *ilist,
1983 const t_mdatoms *mdatoms,
1984 gmx_bool bLimitRange,
1988 int vsite_atom_range, natperthread;
1993 int nral1, inc, i, j;
1995 if (vsite->nthreads == 1)
2001 #pragma omp parallel for num_threads(vsite->nthreads) schedule(static)
2002 for (th = 0; th < vsite->nthreads; th++)
2004 prepare_vsite_thread(ilist, &vsite->tdata[th]);
2006 /* Master threads does the (potential) overlap vsites */
2007 prepare_vsite_thread(ilist, &vsite->tdata[vsite->nthreads]);
2009 /* The current way of distributing the vsites over threads in primitive.
2010 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
2011 * without taking into account how the vsites are distributed.
2012 * Without domain decomposition we bLimitRange=TRUE and we at least
2013 * tighten the upper bound of the range (useful for common systems
2014 * such as a vsite-protein in 3-site water).
2018 vsite_atom_range = -1;
2019 for (ftype = 0; ftype < F_NRE; ftype++)
2021 if ((interaction_function[ftype].flags & IF_VSITE) &&
2024 nral1 = 1 + NRAL(ftype);
2025 iat = ilist[ftype].iatoms;
2026 for (i = 0; i < ilist[ftype].nr; i += nral1)
2028 for (j = i+1; j < i+nral1; j++)
2030 vsite_atom_range = max(vsite_atom_range, iat[j]);
2039 vsite_atom_range = mdatoms->homenr;
2041 natperthread = (vsite_atom_range + vsite->nthreads - 1)/vsite->nthreads;
2045 fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n", mdatoms->nr, vsite_atom_range, natperthread);
2048 /* To simplify the vsite assignment, we make an index which tells us
2049 * to which thread particles, both non-vsites and vsites, are assigned.
2051 if (mdatoms->nr > vsite->th_ind_nalloc)
2053 vsite->th_ind_nalloc = over_alloc_large(mdatoms->nr);
2054 srenew(vsite->th_ind, vsite->th_ind_nalloc);
2056 th_ind = vsite->th_ind;
2058 for (i = 0; i < mdatoms->nr; i++)
2060 if (mdatoms->ptype[i] == eptVSite)
2062 /* vsites are not assigned to a thread yet */
2067 /* assign non-vsite particles to thread th */
2070 if (i == (th + 1)*natperthread && th < vsite->nthreads)
2076 for (ftype = 0; ftype < F_NRE; ftype++)
2078 if ((interaction_function[ftype].flags & IF_VSITE) &&
2081 nral1 = 1 + NRAL(ftype);
2083 iat = ilist[ftype].iatoms;
2084 for (i = 0; i < ilist[ftype].nr; )
2086 th = iat[1+i]/natperthread;
2087 /* We would like to assign this vsite the thread th,
2088 * but it might depend on atoms outside the atom range of th
2089 * or on another vsite not assigned to thread th.
2091 if (ftype != F_VSITEN)
2093 for (j = i+2; j < i+nral1; j++)
2095 if (th_ind[iat[j]] != th)
2097 /* Some constructing atoms are not assigned to
2098 * thread th, move this vsite to a separate batch.
2100 th = vsite->nthreads;
2107 for (j = i+2; j < i+inc; j += 3)
2109 if (th_ind[iat[j]] != th)
2111 th = vsite->nthreads;
2115 /* Copy this vsite to the thread data struct of thread th */
2116 il_th = &vsite->tdata[th].ilist[ftype];
2117 for (j = i; j < i+inc; j++)
2119 il_th->iatoms[il_th->nr++] = iat[j];
2121 /* Update this vsite's thread index entry */
2122 th_ind[iat[1+i]] = th;
2131 for (ftype = 0; ftype < F_NRE; ftype++)
2133 if ((interaction_function[ftype].flags & IF_VSITE) &&
2134 ilist[ftype].nr > 0)
2136 fprintf(debug, "%-20s thread dist:",
2137 interaction_function[ftype].longname);
2138 for (th = 0; th < vsite->nthreads+1; th++)
2140 fprintf(debug, " %4d", vsite->tdata[th].ilist[ftype].nr);
2142 fprintf(debug, "\n");
2148 void set_vsite_top(gmx_vsite_t *vsite, gmx_localtop_t *top, t_mdatoms *md,
2153 if (vsite->n_intercg_vsite > 0)
2155 if (vsite->bHaveChargeGroups)
2157 /* Make an atom to charge group index */
2158 a2cg = atom2cg(&top->cgs);
2159 vsite->vsite_pbc_loc = get_vsite_pbc(top->idef.iparams,
2160 top->idef.il, NULL, md,
2167 snew(vsite->vsitecomm, 1);
2168 vsite->bPDvsitecomm =
2169 setup_parallel_vsites(&(top->idef), cr, vsite->vsitecomm);
2173 if (vsite->nthreads > 1)
2175 if (vsite->bHaveChargeGroups || PARTDECOMP(cr))
2177 gmx_incons("Can not use threads virtual sites combined with charge groups or particle decomposition");
2180 split_vsites_over_threads(top->idef.il, md, !DOMAINDECOMP(cr), vsite);