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 */
479 /* Constants for constructing vsites */
481 /* Check what kind of pbc we need to use */
484 /* No charge groups, vsite follows its own pbc */
486 copy_rvec(x[avsite], xpbc);
488 else if (vsite_pbc != NULL)
490 pbc_atom = vsite_pbc[i/(1+nra)];
495 /* We need to copy the coordinates here,
496 * single for single atom cg's pbc_atom
497 * is the vsite itself.
499 copy_rvec(x[pbc_atom], xpbc);
501 pbc_null2 = pbc_null;
512 /* Copy the old position */
513 copy_rvec(x[avsite], xv);
515 /* Construct the vsite depending on type */
520 constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
526 constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
532 constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
538 constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
545 constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
553 constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
562 constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
566 inc = constr_vsiten(ia, ip, x, pbc_null2);
569 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
570 ftype, __FILE__, __LINE__);
575 /* Match the pbc of this vsite to the rest of its charge group */
576 ishift = pbc_dx_aiuc(pbc_null, x[avsite], xpbc, dx);
577 if (ishift != CENTRAL)
579 rvec_add(xpbc, dx, x[avsite]);
584 /* Calculate velocity of vsite... */
585 rvec_sub(x[avsite], xv, vv);
586 svmul(inv_dt, vv, v[avsite]);
589 /* Increment loop variables */
597 void construct_vsites(FILE *log, gmx_vsite_t *vsite,
598 rvec x[], t_nrnb *nrnb,
600 t_iparams ip[], t_ilist ilist[],
601 int ePBC, gmx_bool bMolPBC, t_graph *graph,
602 t_commrec *cr, matrix box)
604 t_pbc pbc, *pbc_null;
608 bDomDec = cr && DOMAINDECOMP(cr);
610 /* We only need to do pbc when we have inter-cg vsites */
611 if (ePBC != epbcNONE && (bDomDec || bMolPBC) && vsite->n_intercg_vsite)
613 /* This is wasting some CPU time as we now do this multiple times
614 * per MD step. But how often do we have vsites with full pbc?
616 pbc_null = set_pbc_dd(&pbc, ePBC, cr != NULL ? cr->dd : NULL, FALSE, box);
627 dd_move_x_vsites(cr->dd, box, x);
629 else if (vsite->bPDvsitecomm)
631 /* I'm not sure whether the periodicity and shift are guaranteed
632 * to be consistent between different nodes when running e.g. polymers
633 * in parallel. In this special case we thus unshift/shift,
634 * but only when necessary. This is to make sure the coordinates
635 * we move don't end up a box away...
639 unshift_self(graph, box, x);
642 move_construct_x(vsite->vsitecomm, x, cr);
646 shift_self(graph, box, x);
651 if (vsite->nthreads == 1)
653 construct_vsites_thread(vsite,
660 #pragma omp parallel num_threads(vsite->nthreads)
662 construct_vsites_thread(vsite,
664 ip, vsite->tdata[gmx_omp_get_thread_num()].ilist,
667 /* Now we can construct the vsites that might depend on other vsites */
668 construct_vsites_thread(vsite,
670 ip, vsite->tdata[vsite->nthreads].ilist,
675 static void spread_vsite2(t_iatom ia[], real a,
676 rvec x[], rvec f[], rvec fshift[],
677 t_pbc *pbc, t_graph *g)
689 svmul(1-a, f[av], fi);
690 svmul( a, f[av], fj);
699 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
701 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
706 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
707 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
715 if (fshift && (siv != CENTRAL || sij != CENTRAL))
717 rvec_inc(fshift[siv], f[av]);
718 rvec_dec(fshift[CENTRAL], fi);
719 rvec_dec(fshift[sij], fj);
722 /* TOTAL: 13 flops */
725 void construct_vsites_mtop(FILE *log, gmx_vsite_t *vsite,
726 gmx_mtop_t *mtop, rvec x[])
729 gmx_molblock_t *molb;
733 for (mb = 0; mb < mtop->nmolblock; mb++)
735 molb = &mtop->molblock[mb];
736 molt = &mtop->moltype[molb->type];
737 for (mol = 0; mol < molb->nmol; mol++)
739 construct_vsites(log, vsite, x+as, NULL, 0.0, NULL,
740 mtop->ffparams.iparams, molt->ilist,
741 epbcNONE, TRUE, NULL, NULL, NULL);
742 as += molt->atoms.nr;
747 static void spread_vsite3(t_iatom ia[], real a, real b,
748 rvec x[], rvec f[], rvec fshift[],
749 t_pbc *pbc, t_graph *g)
752 atom_id av, ai, aj, ak;
761 svmul(1-a-b, f[av], fi);
762 svmul( a, f[av], fj);
763 svmul( b, f[av], fk);
773 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ia[1]), di);
775 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
777 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), di);
782 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
783 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
784 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
793 if (fshift && (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL))
795 rvec_inc(fshift[siv], f[av]);
796 rvec_dec(fshift[CENTRAL], fi);
797 rvec_dec(fshift[sij], fj);
798 rvec_dec(fshift[sik], fk);
801 /* TOTAL: 20 flops */
804 static void spread_vsite3FD(t_iatom ia[], real a, real b,
805 rvec x[], rvec f[], rvec fshift[],
806 gmx_bool VirCorr, matrix dxdf,
807 t_pbc *pbc, t_graph *g)
809 real fx, fy, fz, c, invl, fproj, a1;
810 rvec xvi, xij, xjk, xix, fv, temp;
811 t_iatom av, ai, aj, ak;
812 int svi, sji, skj, d;
819 copy_rvec(f[av], fv);
821 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
822 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
825 /* xix goes from i to point x on the line jk */
826 xix[XX] = xij[XX]+a*xjk[XX];
827 xix[YY] = xij[YY]+a*xjk[YY];
828 xix[ZZ] = xij[ZZ]+a*xjk[ZZ];
831 invl = gmx_invsqrt(iprod(xix, xix));
835 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
837 temp[XX] = c*(fv[XX]-fproj*xix[XX]);
838 temp[YY] = c*(fv[YY]-fproj*xix[YY]);
839 temp[ZZ] = c*(fv[ZZ]-fproj*xix[ZZ]);
842 /* c is already calculated in constr_vsite3FD
843 storing c somewhere will save 26 flops! */
846 f[ai][XX] += fv[XX] - temp[XX];
847 f[ai][YY] += fv[YY] - temp[YY];
848 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
849 f[aj][XX] += a1*temp[XX];
850 f[aj][YY] += a1*temp[YY];
851 f[aj][ZZ] += a1*temp[ZZ];
852 f[ak][XX] += a*temp[XX];
853 f[ak][YY] += a*temp[YY];
854 f[ak][ZZ] += a*temp[ZZ];
859 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
861 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
863 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
868 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
875 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
877 rvec_dec(fshift[svi], fv);
878 fshift[CENTRAL][XX] += fv[XX] - (1 + a)*temp[XX];
879 fshift[CENTRAL][YY] += fv[YY] - (1 + a)*temp[YY];
880 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a)*temp[ZZ];
881 fshift[ sji][XX] += temp[XX];
882 fshift[ sji][YY] += temp[YY];
883 fshift[ sji][ZZ] += temp[ZZ];
884 fshift[ skj][XX] += a*temp[XX];
885 fshift[ skj][YY] += a*temp[YY];
886 fshift[ skj][ZZ] += a*temp[ZZ];
891 /* When VirCorr=TRUE, the virial for the current forces is not
892 * calculated from the redistributed forces. This means that
893 * the effect of non-linear virtual site constructions on the virial
894 * needs to be added separately. This contribution can be calculated
895 * in many ways, but the simplest and cheapest way is to use
896 * the first constructing atom ai as a reference position in space:
897 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
902 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
904 for (i = 0; i < DIM; i++)
906 for (j = 0; j < DIM; j++)
908 /* As xix is a linear combination of j and k, use that here */
909 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
914 /* TOTAL: 61 flops */
917 static void spread_vsite3FAD(t_iatom ia[], real a, real b,
918 rvec x[], rvec f[], rvec fshift[],
919 gmx_bool VirCorr, matrix dxdf,
920 t_pbc *pbc, t_graph *g)
922 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
923 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
924 t_iatom av, ai, aj, ak;
925 int svi, sji, skj, d;
932 copy_rvec(f[ia[1]], fv);
934 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
935 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
938 invdij = gmx_invsqrt(iprod(xij, xij));
939 invdij2 = invdij * invdij;
940 c1 = iprod(xij, xjk) * invdij2;
941 xperp[XX] = xjk[XX] - c1*xij[XX];
942 xperp[YY] = xjk[YY] - c1*xij[YY];
943 xperp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
944 /* xperp in plane ijk, perp. to ij */
945 invdp = gmx_invsqrt(iprod(xperp, xperp));
950 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
951 storing them somewhere will save 45 flops! */
953 fproj = iprod(xij, fv)*invdij2;
954 svmul(fproj, xij, Fpij); /* proj. f on xij */
955 svmul(iprod(xperp, fv)*invdp*invdp, xperp, Fppp); /* proj. f on xperp */
956 svmul(b1*fproj, xperp, f3);
959 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
960 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
961 for (d = 0; (d < DIM); d++)
969 f[ai][XX] += fv[XX] - f1[XX] + c1*f2[XX] + f3[XX];
970 f[ai][YY] += fv[YY] - f1[YY] + c1*f2[YY] + f3[YY];
971 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1*f2[ZZ] + f3[ZZ];
972 f[aj][XX] += f1[XX] - c2*f2[XX] - f3[XX];
973 f[aj][YY] += f1[YY] - c2*f2[YY] - f3[YY];
974 f[aj][ZZ] += f1[ZZ] - c2*f2[ZZ] - f3[ZZ];
982 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
984 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
986 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
991 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
998 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
1000 rvec_dec(fshift[svi], fv);
1001 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1-c1)*f2[XX] + f3[XX];
1002 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1-c1)*f2[YY] + f3[YY];
1003 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1-c1)*f2[ZZ] + f3[ZZ];
1004 fshift[ sji][XX] += f1[XX] - c1 *f2[XX] - f3[XX];
1005 fshift[ sji][YY] += f1[YY] - c1 *f2[YY] - f3[YY];
1006 fshift[ sji][ZZ] += f1[ZZ] - c1 *f2[ZZ] - f3[ZZ];
1007 fshift[ skj][XX] += f2[XX];
1008 fshift[ skj][YY] += f2[YY];
1009 fshift[ skj][ZZ] += f2[ZZ];
1017 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1019 for (i = 0; i < DIM; i++)
1021 for (j = 0; j < DIM; j++)
1023 /* Note that xik=xij+xjk, so we have to add xij*f2 */
1026 + xij[i]*(f1[j] + (1 - c2)*f2[j] - f3[j])
1032 /* TOTAL: 113 flops */
1035 static void spread_vsite3OUT(t_iatom ia[], real a, real b, real c,
1036 rvec x[], rvec f[], rvec fshift[],
1037 gmx_bool VirCorr, matrix dxdf,
1038 t_pbc *pbc, t_graph *g)
1040 rvec xvi, xij, xik, fv, fj, fk;
1042 atom_id av, ai, aj, ak;
1051 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1052 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1055 copy_rvec(f[av], fv);
1062 fj[XX] = a*fv[XX] - xik[ZZ]*cfy + xik[YY]*cfz;
1063 fj[YY] = xik[ZZ]*cfx + a*fv[YY] - xik[XX]*cfz;
1064 fj[ZZ] = -xik[YY]*cfx + xik[XX]*cfy + a*fv[ZZ];
1066 fk[XX] = b*fv[XX] + xij[ZZ]*cfy - xij[YY]*cfz;
1067 fk[YY] = -xij[ZZ]*cfx + b*fv[YY] + xij[XX]*cfz;
1068 fk[ZZ] = xij[YY]*cfx - xij[XX]*cfy + b*fv[ZZ];
1071 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
1072 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
1073 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1074 rvec_inc(f[aj], fj);
1075 rvec_inc(f[ak], fk);
1080 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1082 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1084 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1089 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1096 if (fshift && (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL))
1098 rvec_dec(fshift[svi], fv);
1099 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
1100 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
1101 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
1102 rvec_inc(fshift[sji], fj);
1103 rvec_inc(fshift[ski], fk);
1111 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1113 for (i = 0; i < DIM; i++)
1115 for (j = 0; j < DIM; j++)
1117 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j];
1122 /* TOTAL: 54 flops */
1125 static void spread_vsite4FD(t_iatom ia[], real a, real b, real c,
1126 rvec x[], rvec f[], rvec fshift[],
1127 gmx_bool VirCorr, matrix dxdf,
1128 t_pbc *pbc, t_graph *g)
1130 real d, invl, fproj, a1;
1131 rvec xvi, xij, xjk, xjl, xix, fv, temp;
1132 atom_id av, ai, aj, ak, al;
1134 int svi, sji, skj, slj, m;
1142 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1143 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1144 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1147 /* xix goes from i to point x on the plane jkl */
1148 for (m = 0; m < DIM; m++)
1150 xix[m] = xij[m] + a*xjk[m] + b*xjl[m];
1154 invl = gmx_invsqrt(iprod(xix, xix));
1156 /* 4 + ?10? flops */
1158 copy_rvec(f[av], fv);
1160 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
1162 for (m = 0; m < DIM; m++)
1164 temp[m] = d*(fv[m] - fproj*xix[m]);
1168 /* c is already calculated in constr_vsite3FD
1169 storing c somewhere will save 35 flops! */
1172 for (m = 0; m < DIM; m++)
1174 f[ai][m] += fv[m] - temp[m];
1175 f[aj][m] += a1*temp[m];
1176 f[ak][m] += a*temp[m];
1177 f[al][m] += b*temp[m];
1183 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1185 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1187 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
1189 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, aj), di);
1194 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1202 (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL))
1204 rvec_dec(fshift[svi], fv);
1205 for (m = 0; m < DIM; m++)
1207 fshift[CENTRAL][m] += fv[m] - (1 + a + b)*temp[m];
1208 fshift[ sji][m] += temp[m];
1209 fshift[ skj][m] += a*temp[m];
1210 fshift[ slj][m] += b*temp[m];
1219 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1221 for (i = 0; i < DIM; i++)
1223 for (j = 0; j < DIM; j++)
1225 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
1230 /* TOTAL: 77 flops */
1234 static void spread_vsite4FDN(t_iatom ia[], real a, real b, real c,
1235 rvec x[], rvec f[], rvec fshift[],
1236 gmx_bool VirCorr, matrix dxdf,
1237 t_pbc *pbc, t_graph *g)
1239 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1240 rvec fv, fj, fk, fl;
1244 int av, ai, aj, ak, al;
1245 int svi, sij, sik, sil;
1247 /* DEBUG: check atom indices */
1254 copy_rvec(f[av], fv);
1256 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1257 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1258 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1271 rvec_sub(ra, xij, rja);
1272 rvec_sub(rb, xij, rjb);
1273 rvec_sub(rb, ra, rab);
1276 cprod(rja, rjb, rm);
1279 invrm = gmx_invsqrt(norm2(rm));
1280 denom = invrm*invrm;
1283 cfx = c*invrm*fv[XX];
1284 cfy = c*invrm*fv[YY];
1285 cfz = c*invrm*fv[ZZ];
1296 fj[XX] = ( -rm[XX]*rt[XX]) * cfx + ( rab[ZZ]-rm[YY]*rt[XX]) * cfy + (-rab[YY]-rm[ZZ]*rt[XX]) * cfz;
1297 fj[YY] = (-rab[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( rab[XX]-rm[ZZ]*rt[YY]) * cfz;
1298 fj[ZZ] = ( rab[YY]-rm[XX]*rt[ZZ]) * cfx + (-rab[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1309 fk[XX] = ( -rm[XX]*rt[XX]) * cfx + (-a*rjb[ZZ]-rm[YY]*rt[XX]) * cfy + ( a*rjb[YY]-rm[ZZ]*rt[XX]) * cfz;
1310 fk[YY] = ( a*rjb[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + (-a*rjb[XX]-rm[ZZ]*rt[YY]) * cfz;
1311 fk[ZZ] = (-a*rjb[YY]-rm[XX]*rt[ZZ]) * cfx + ( a*rjb[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1322 fl[XX] = ( -rm[XX]*rt[XX]) * cfx + ( b*rja[ZZ]-rm[YY]*rt[XX]) * cfy + (-b*rja[YY]-rm[ZZ]*rt[XX]) * cfz;
1323 fl[YY] = (-b*rja[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( b*rja[XX]-rm[ZZ]*rt[YY]) * cfz;
1324 fl[ZZ] = ( b*rja[YY]-rm[XX]*rt[ZZ]) * cfx + (-b*rja[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1327 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1328 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1329 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1330 rvec_inc(f[aj], fj);
1331 rvec_inc(f[ak], fk);
1332 rvec_inc(f[al], fl);
1337 ivec_sub(SHIFT_IVEC(g, av), SHIFT_IVEC(g, ai), di);
1339 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1341 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1343 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, ai), di);
1348 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1355 if (fshift && (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL))
1357 rvec_dec(fshift[svi], fv);
1358 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1359 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1360 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1361 rvec_inc(fshift[sij], fj);
1362 rvec_inc(fshift[sik], fk);
1363 rvec_inc(fshift[sil], fl);
1371 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1373 for (i = 0; i < DIM; i++)
1375 for (j = 0; j < DIM; j++)
1377 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j] + xil[i]*fl[j];
1382 /* Total: 207 flops (Yuck!) */
1386 static int spread_vsiten(t_iatom ia[], t_iparams ip[],
1387 rvec x[], rvec f[], rvec fshift[],
1388 t_pbc *pbc, t_graph *g)
1396 n3 = 3*ip[ia[0]].vsiten.n;
1398 copy_rvec(x[av], xv);
1400 for (i = 0; i < n3; i += 3)
1405 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
1410 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1416 a = ip[ia[i]].vsiten.a;
1417 svmul(a, f[av], fi);
1418 rvec_inc(f[ai], fi);
1419 if (fshift && siv != CENTRAL)
1421 rvec_inc(fshift[siv], fi);
1422 rvec_dec(fshift[CENTRAL], fi);
1431 static int vsite_count(const t_ilist *ilist, int ftype)
1433 if (ftype == F_VSITEN)
1435 return ilist[ftype].nr/3;
1439 return ilist[ftype].nr/(1 + interaction_function[ftype].nratoms);
1443 static void spread_vsite_f_thread(gmx_vsite_t *vsite,
1444 rvec x[], rvec f[], rvec *fshift,
1445 gmx_bool VirCorr, matrix dxdf,
1446 t_iparams ip[], t_ilist ilist[],
1447 t_graph *g, t_pbc *pbc_null)
1451 int i, inc, m, nra, nr, tp, ftype;
1461 bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
1463 /* this loop goes backwards to be able to build *
1464 * higher type vsites from lower types */
1467 for (ftype = F_NRE-1; (ftype >= 0); ftype--)
1469 if ((interaction_function[ftype].flags & IF_VSITE) &&
1470 ilist[ftype].nr > 0)
1472 nra = interaction_function[ftype].nratoms;
1474 nr = ilist[ftype].nr;
1475 ia = ilist[ftype].iatoms;
1479 pbc_null2 = pbc_null;
1481 else if (pbc_null != NULL)
1483 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
1486 for (i = 0; i < nr; )
1488 if (vsite_pbc != NULL)
1490 if (vsite_pbc[i/(1+nra)] > -2)
1492 pbc_null2 = pbc_null;
1502 /* Constants for constructing */
1503 a1 = ip[tp].vsite.a;
1504 /* Construct the vsite depending on type */
1508 spread_vsite2(ia, a1, x, f, fshift, pbc_null2, g);
1511 b1 = ip[tp].vsite.b;
1512 spread_vsite3(ia, a1, b1, x, f, fshift, pbc_null2, g);
1515 b1 = ip[tp].vsite.b;
1516 spread_vsite3FD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1519 b1 = ip[tp].vsite.b;
1520 spread_vsite3FAD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1523 b1 = ip[tp].vsite.b;
1524 c1 = ip[tp].vsite.c;
1525 spread_vsite3OUT(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1528 b1 = ip[tp].vsite.b;
1529 c1 = ip[tp].vsite.c;
1530 spread_vsite4FD(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1533 b1 = ip[tp].vsite.b;
1534 c1 = ip[tp].vsite.c;
1535 spread_vsite4FDN(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1538 inc = spread_vsiten(ia, ip, x, f, fshift, pbc_null2, g);
1541 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
1542 ftype, __FILE__, __LINE__);
1544 clear_rvec(f[ia[1]]);
1546 /* Increment loop variables */
1554 void spread_vsite_f(FILE *log, gmx_vsite_t *vsite,
1555 rvec x[], rvec f[], rvec *fshift,
1556 gmx_bool VirCorr, matrix vir,
1557 t_nrnb *nrnb, t_idef *idef,
1558 int ePBC, gmx_bool bMolPBC, t_graph *g, matrix box,
1561 t_pbc pbc, *pbc_null;
1564 /* We only need to do pbc when we have inter-cg vsites */
1565 if ((DOMAINDECOMP(cr) || bMolPBC) && vsite->n_intercg_vsite)
1567 /* This is wasting some CPU time as we now do this multiple times
1568 * per MD step. But how often do we have vsites with full pbc?
1570 pbc_null = set_pbc_dd(&pbc, ePBC, cr->dd, FALSE, box);
1577 if (DOMAINDECOMP(cr))
1579 dd_clear_f_vsites(cr->dd, f);
1581 else if (PARTDECOMP(cr) && vsite->vsitecomm != NULL)
1583 pd_clear_nonlocal_constructs(vsite->vsitecomm, f);
1586 if (vsite->nthreads == 1)
1588 spread_vsite_f_thread(vsite,
1590 VirCorr, vsite->tdata[0].dxdf,
1591 idef->iparams, idef->il,
1596 /* First spread the vsites that might depend on other vsites */
1597 spread_vsite_f_thread(vsite,
1599 VirCorr, vsite->tdata[vsite->nthreads].dxdf,
1601 vsite->tdata[vsite->nthreads].ilist,
1604 #pragma omp parallel num_threads(vsite->nthreads)
1609 thread = gmx_omp_get_thread_num();
1611 if (thread == 0 || fshift == NULL)
1619 fshift_t = vsite->tdata[thread].fshift;
1621 for (i = 0; i < SHIFTS; i++)
1623 clear_rvec(fshift_t[i]);
1627 spread_vsite_f_thread(vsite,
1629 VirCorr, vsite->tdata[thread].dxdf,
1631 vsite->tdata[thread].ilist,
1639 for (th = 1; th < vsite->nthreads; th++)
1641 for (i = 0; i < SHIFTS; i++)
1643 rvec_inc(fshift[i], vsite->tdata[th].fshift[i]);
1653 for (th = 0; th < (vsite->nthreads == 1 ? 1 : vsite->nthreads+1); th++)
1655 for (i = 0; i < DIM; i++)
1657 for (j = 0; j < DIM; j++)
1659 vir[i][j] += -0.5*vsite->tdata[th].dxdf[i][j];
1665 if (DOMAINDECOMP(cr))
1667 dd_move_f_vsites(cr->dd, f, fshift);
1669 else if (vsite->bPDvsitecomm)
1671 /* We only move forces here, and they are independent of shifts */
1672 move_construct_f(vsite->vsitecomm, f, cr);
1675 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(idef->il, F_VSITE2));
1676 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(idef->il, F_VSITE3));
1677 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(idef->il, F_VSITE3FD));
1678 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(idef->il, F_VSITE3FAD));
1679 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(idef->il, F_VSITE3OUT));
1680 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(idef->il, F_VSITE4FD));
1681 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(idef->il, F_VSITE4FDN));
1682 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(idef->il, F_VSITEN));
1685 static int *atom2cg(t_block *cgs)
1689 snew(a2cg, cgs->index[cgs->nr]);
1690 for (cg = 0; cg < cgs->nr; cg++)
1692 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
1701 static int count_intercg_vsite(gmx_mtop_t *mtop,
1702 gmx_bool *bHaveChargeGroups)
1704 int mb, mt, ftype, nral, i, cg, a;
1705 gmx_molblock_t *molb;
1706 gmx_moltype_t *molt;
1710 int n_intercg_vsite;
1712 *bHaveChargeGroups = FALSE;
1714 n_intercg_vsite = 0;
1715 for (mb = 0; mb < mtop->nmolblock; mb++)
1717 molb = &mtop->molblock[mb];
1718 molt = &mtop->moltype[molb->type];
1720 if (molt->cgs.nr < molt->atoms.nr)
1722 *bHaveChargeGroups = TRUE;
1725 a2cg = atom2cg(&molt->cgs);
1726 for (ftype = 0; ftype < F_NRE; ftype++)
1728 if (interaction_function[ftype].flags & IF_VSITE)
1731 il = &molt->ilist[ftype];
1733 for (i = 0; i < il->nr; i += 1+nral)
1736 for (a = 1; a < nral; a++)
1738 if (a2cg[ia[1+a]] != cg)
1740 n_intercg_vsite += molb->nmol;
1750 return n_intercg_vsite;
1753 static int **get_vsite_pbc(t_iparams *iparams, t_ilist *ilist,
1754 t_atom *atom, t_mdatoms *md,
1755 t_block *cgs, int *a2cg)
1757 int ftype, nral, i, j, vsi, vsite, cg_v, cg_c, a, nc3 = 0;
1760 int **vsite_pbc, *vsite_pbc_f;
1762 gmx_bool bViteOnlyCG_and_FirstAtom;
1764 /* Make an array that tells if the pbc of an atom is set */
1765 snew(pbc_set, cgs->index[cgs->nr]);
1766 /* PBC is set for all non vsites */
1767 for (a = 0; a < cgs->index[cgs->nr]; a++)
1769 if ((atom && atom[a].ptype != eptVSite) ||
1770 (md && md->ptype[a] != eptVSite))
1776 snew(vsite_pbc, F_VSITEN-F_VSITE2+1);
1778 for (ftype = 0; ftype < F_NRE; ftype++)
1780 if (interaction_function[ftype].flags & IF_VSITE)
1786 snew(vsite_pbc[ftype-F_VSITE2], il->nr/(1+nral));
1787 vsite_pbc_f = vsite_pbc[ftype-F_VSITE2];
1795 /* A value of -2 signals that this vsite and its contructing
1796 * atoms are all within the same cg, so no pbc is required.
1798 vsite_pbc_f[vsi] = -2;
1799 /* Check if constructing atoms are outside the vsite's cg */
1801 if (ftype == F_VSITEN)
1803 nc3 = 3*iparams[ia[i]].vsiten.n;
1804 for (j = 0; j < nc3; j += 3)
1806 if (a2cg[ia[i+j+2]] != cg_v)
1808 vsite_pbc_f[vsi] = -1;
1814 for (a = 1; a < nral; a++)
1816 if (a2cg[ia[i+1+a]] != cg_v)
1818 vsite_pbc_f[vsi] = -1;
1822 if (vsite_pbc_f[vsi] == -1)
1824 /* Check if this is the first processed atom of a vsite only cg */
1825 bViteOnlyCG_and_FirstAtom = TRUE;
1826 for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1828 /* Non-vsites already have pbc set, so simply check for pbc_set */
1831 bViteOnlyCG_and_FirstAtom = FALSE;
1835 if (bViteOnlyCG_and_FirstAtom)
1837 /* First processed atom of a vsite only charge group.
1838 * The pbc of the input coordinates to construct_vsites
1839 * should be preserved.
1841 vsite_pbc_f[vsi] = vsite;
1843 else if (cg_v != a2cg[ia[1+i+1]])
1845 /* This vsite has a different charge group index
1846 * than it's first constructing atom
1847 * and the charge group has more than one atom,
1848 * search for the first normal particle
1849 * or vsite that already had its pbc defined.
1850 * If nothing is found, use full pbc for this vsite.
1852 for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1854 if (a != vsite && pbc_set[a])
1856 vsite_pbc_f[vsi] = a;
1859 fprintf(debug, "vsite %d match pbc with atom %d\n",
1867 fprintf(debug, "vsite atom %d cg %d - %d pbc atom %d\n",
1868 vsite+1, cgs->index[cg_v]+1, cgs->index[cg_v+1],
1869 vsite_pbc_f[vsi]+1);
1873 if (ftype == F_VSITEN)
1875 /* The other entries in vsite_pbc_f are not used for center vsites */
1883 /* This vsite now has its pbc defined */
1895 gmx_vsite_t *init_vsite(gmx_mtop_t *mtop, t_commrec *cr,
1896 gmx_bool bSerial_NoPBC)
1902 gmx_moltype_t *molt;
1905 /* check if there are vsites */
1907 for (i = 0; i < F_NRE; i++)
1909 if (interaction_function[i].flags & IF_VSITE)
1911 nvsite += gmx_mtop_ftype_count(mtop, i);
1922 vsite->n_intercg_vsite = count_intercg_vsite(mtop,
1923 &vsite->bHaveChargeGroups);
1925 /* If we don't have charge groups, the vsite follows its own pbc */
1926 if (!bSerial_NoPBC &&
1927 vsite->bHaveChargeGroups &&
1928 vsite->n_intercg_vsite > 0 && DOMAINDECOMP(cr))
1930 vsite->nvsite_pbc_molt = mtop->nmoltype;
1931 snew(vsite->vsite_pbc_molt, vsite->nvsite_pbc_molt);
1932 for (mt = 0; mt < mtop->nmoltype; mt++)
1934 molt = &mtop->moltype[mt];
1935 /* Make an atom to charge group index */
1936 a2cg = atom2cg(&molt->cgs);
1937 vsite->vsite_pbc_molt[mt] = get_vsite_pbc(mtop->ffparams.iparams,
1939 molt->atoms.atom, NULL,
1944 snew(vsite->vsite_pbc_loc_nalloc, F_VSITEN-F_VSITE2+1);
1945 snew(vsite->vsite_pbc_loc, F_VSITEN-F_VSITE2+1);
1950 vsite->nthreads = 1;
1954 vsite->nthreads = gmx_omp_nthreads_get(emntVSITE);
1958 /* We need one extra thread data structure for the overlap vsites */
1959 snew(vsite->tdata, vsite->nthreads+1);
1962 vsite->th_ind = NULL;
1963 vsite->th_ind_nalloc = 0;
1968 static void prepare_vsite_thread(const t_ilist *ilist,
1969 gmx_vsite_thread_t *vsite_th)
1973 for (ftype = 0; ftype < F_NRE; ftype++)
1975 if (interaction_function[ftype].flags & IF_VSITE)
1977 if (ilist[ftype].nr > vsite_th->ilist[ftype].nalloc)
1979 vsite_th->ilist[ftype].nalloc = over_alloc_large(ilist[ftype].nr);
1980 srenew(vsite_th->ilist[ftype].iatoms, vsite_th->ilist[ftype].nalloc);
1983 vsite_th->ilist[ftype].nr = 0;
1988 void split_vsites_over_threads(const t_ilist *ilist,
1989 const t_mdatoms *mdatoms,
1990 gmx_bool bLimitRange,
1994 int vsite_atom_range, natperthread;
1999 int nral1, inc, i, j;
2001 if (vsite->nthreads == 1)
2007 #pragma omp parallel for num_threads(vsite->nthreads) schedule(static)
2008 for (th = 0; th < vsite->nthreads; th++)
2010 prepare_vsite_thread(ilist, &vsite->tdata[th]);
2012 /* Master threads does the (potential) overlap vsites */
2013 prepare_vsite_thread(ilist, &vsite->tdata[vsite->nthreads]);
2015 /* The current way of distributing the vsites over threads in primitive.
2016 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
2017 * without taking into account how the vsites are distributed.
2018 * Without domain decomposition we bLimitRange=TRUE and we at least
2019 * tighten the upper bound of the range (useful for common systems
2020 * such as a vsite-protein in 3-site water).
2024 vsite_atom_range = -1;
2025 for (ftype = 0; ftype < F_NRE; ftype++)
2027 if ((interaction_function[ftype].flags & IF_VSITE) &&
2030 nral1 = 1 + NRAL(ftype);
2031 iat = ilist[ftype].iatoms;
2032 for (i = 0; i < ilist[ftype].nr; i += nral1)
2034 for (j = i+1; j < i+nral1; j++)
2036 vsite_atom_range = max(vsite_atom_range, iat[j]);
2045 vsite_atom_range = mdatoms->homenr;
2047 natperthread = (vsite_atom_range + vsite->nthreads - 1)/vsite->nthreads;
2051 fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n", mdatoms->nr, vsite_atom_range, natperthread);
2054 /* To simplify the vsite assignment, we make an index which tells us
2055 * to which thread particles, both non-vsites and vsites, are assigned.
2057 if (mdatoms->nr > vsite->th_ind_nalloc)
2059 vsite->th_ind_nalloc = over_alloc_large(mdatoms->nr);
2060 srenew(vsite->th_ind, vsite->th_ind_nalloc);
2062 th_ind = vsite->th_ind;
2064 for (i = 0; i < mdatoms->nr; i++)
2066 if (mdatoms->ptype[i] == eptVSite)
2068 /* vsites are not assigned to a thread yet */
2073 /* assign non-vsite particles to thread th */
2076 if (i == (th + 1)*natperthread && th < vsite->nthreads)
2082 for (ftype = 0; ftype < F_NRE; ftype++)
2084 if ((interaction_function[ftype].flags & IF_VSITE) &&
2087 nral1 = 1 + NRAL(ftype);
2089 iat = ilist[ftype].iatoms;
2090 for (i = 0; i < ilist[ftype].nr; )
2092 th = iat[1+i]/natperthread;
2093 /* We would like to assign this vsite the thread th,
2094 * but it might depend on atoms outside the atom range of th
2095 * or on another vsite not assigned to thread th.
2097 if (ftype != F_VSITEN)
2099 for (j = i+2; j < i+nral1; j++)
2101 if (th_ind[iat[j]] != th)
2103 /* Some constructing atoms are not assigned to
2104 * thread th, move this vsite to a separate batch.
2106 th = vsite->nthreads;
2113 for (j = i+2; j < i+inc; j += 3)
2115 if (th_ind[iat[j]] != th)
2117 th = vsite->nthreads;
2121 /* Copy this vsite to the thread data struct of thread th */
2122 il_th = &vsite->tdata[th].ilist[ftype];
2123 for (j = i; j < i+inc; j++)
2125 il_th->iatoms[il_th->nr++] = iat[j];
2127 /* Update this vsite's thread index entry */
2128 th_ind[iat[1+i]] = th;
2137 for (ftype = 0; ftype < F_NRE; ftype++)
2139 if ((interaction_function[ftype].flags & IF_VSITE) &&
2140 ilist[ftype].nr > 0)
2142 fprintf(debug, "%-20s thread dist:",
2143 interaction_function[ftype].longname);
2144 for (th = 0; th < vsite->nthreads+1; th++)
2146 fprintf(debug, " %4d", vsite->tdata[th].ilist[ftype].nr);
2148 fprintf(debug, "\n");
2154 void set_vsite_top(gmx_vsite_t *vsite, gmx_localtop_t *top, t_mdatoms *md,
2159 if (vsite->n_intercg_vsite > 0)
2161 if (vsite->bHaveChargeGroups)
2163 /* Make an atom to charge group index */
2164 a2cg = atom2cg(&top->cgs);
2165 vsite->vsite_pbc_loc = get_vsite_pbc(top->idef.iparams,
2166 top->idef.il, NULL, md,
2173 snew(vsite->vsitecomm, 1);
2174 vsite->bPDvsitecomm =
2175 setup_parallel_vsites(&(top->idef), cr, vsite->vsitecomm);
2179 if (vsite->nthreads > 1)
2181 if (vsite->bHaveChargeGroups || PARTDECOMP(cr))
2183 gmx_incons("Can not use threads virtual sites combined with charge groups or particle decomposition");
2186 split_vsites_over_threads(top->idef.il, md, !DOMAINDECOMP(cr), vsite);