Merge branch release-4-6 into release-5-0
[alexxy/gromacs.git] / src / gromacs / mdlib / vsite.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <stdio.h>
42 #include "typedefs.h"
43 #include "vsite.h"
44 #include "macros.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "nrnb.h"
47 #include "vec.h"
48 #include "network.h"
49 #include "mshift.h"
50 #include "pbc.h"
51 #include "domdec.h"
52 #include "mtop_util.h"
53 #include "gmx_omp_nthreads.h"
54
55 #include "gromacs/utility/gmxomp.h"
56
57 /* Routines to send/recieve coordinates and force
58  * of constructing atoms.
59  */
60
61 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
62 {
63     if (pbc)
64     {
65         return pbc_dx_aiuc(pbc, xi, xj, dx);
66     }
67     else
68     {
69         rvec_sub(xi, xj, dx);
70         return CENTRAL;
71     }
72 }
73
74 /* Vsite construction routines */
75
76 static void constr_vsite2(rvec xi, rvec xj, rvec x, real a, t_pbc *pbc)
77 {
78     real b;
79     rvec dx;
80
81     b = 1.0-a;
82     /* 1 flop */
83
84     if (pbc)
85     {
86         pbc_dx_aiuc(pbc, xj, xi, dx);
87         x[XX] = xi[XX] + a*dx[XX];
88         x[YY] = xi[YY] + a*dx[YY];
89         x[ZZ] = xi[ZZ] + a*dx[ZZ];
90     }
91     else
92     {
93         x[XX] = b*xi[XX] + a*xj[XX];
94         x[YY] = b*xi[YY] + a*xj[YY];
95         x[ZZ] = b*xi[ZZ] + a*xj[ZZ];
96         /* 9 Flops */
97     }
98
99     /* TOTAL: 10 flops */
100 }
101
102 static void constr_vsite3(rvec xi, rvec xj, rvec xk, rvec x, real a, real b,
103                           t_pbc *pbc)
104 {
105     real c;
106     rvec dxj, dxk;
107
108     c = 1.0-a-b;
109     /* 2 flops */
110
111     if (pbc)
112     {
113         pbc_dx_aiuc(pbc, xj, xi, dxj);
114         pbc_dx_aiuc(pbc, xk, xi, dxk);
115         x[XX] = xi[XX] + a*dxj[XX] + b*dxk[XX];
116         x[YY] = xi[YY] + a*dxj[YY] + b*dxk[YY];
117         x[ZZ] = xi[ZZ] + a*dxj[ZZ] + b*dxk[ZZ];
118     }
119     else
120     {
121         x[XX] = c*xi[XX] + a*xj[XX] + b*xk[XX];
122         x[YY] = c*xi[YY] + a*xj[YY] + b*xk[YY];
123         x[ZZ] = c*xi[ZZ] + a*xj[ZZ] + b*xk[ZZ];
124         /* 15 Flops */
125     }
126
127     /* TOTAL: 17 flops */
128 }
129
130 static void constr_vsite3FD(rvec xi, rvec xj, rvec xk, rvec x, real a, real b,
131                             t_pbc *pbc)
132 {
133     rvec xij, xjk, temp;
134     real c;
135
136     pbc_rvec_sub(pbc, xj, xi, xij);
137     pbc_rvec_sub(pbc, xk, xj, xjk);
138     /* 6 flops */
139
140     /* temp goes from i to a point on the line jk */
141     temp[XX] = xij[XX] + a*xjk[XX];
142     temp[YY] = xij[YY] + a*xjk[YY];
143     temp[ZZ] = xij[ZZ] + a*xjk[ZZ];
144     /* 6 flops */
145
146     c = b*gmx_invsqrt(iprod(temp, temp));
147     /* 6 + 10 flops */
148
149     x[XX] = xi[XX] + c*temp[XX];
150     x[YY] = xi[YY] + c*temp[YY];
151     x[ZZ] = xi[ZZ] + c*temp[ZZ];
152     /* 6 Flops */
153
154     /* TOTAL: 34 flops */
155 }
156
157 static void constr_vsite3FAD(rvec xi, rvec xj, rvec xk, rvec x, real a, real b, t_pbc *pbc)
158 {
159     rvec xij, xjk, xp;
160     real a1, b1, c1, invdij;
161
162     pbc_rvec_sub(pbc, xj, xi, xij);
163     pbc_rvec_sub(pbc, xk, xj, xjk);
164     /* 6 flops */
165
166     invdij = gmx_invsqrt(iprod(xij, xij));
167     c1     = invdij * invdij * iprod(xij, xjk);
168     xp[XX] = xjk[XX] - c1*xij[XX];
169     xp[YY] = xjk[YY] - c1*xij[YY];
170     xp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
171     a1     = a*invdij;
172     b1     = b*gmx_invsqrt(iprod(xp, xp));
173     /* 45 */
174
175     x[XX] = xi[XX] + a1*xij[XX] + b1*xp[XX];
176     x[YY] = xi[YY] + a1*xij[YY] + b1*xp[YY];
177     x[ZZ] = xi[ZZ] + a1*xij[ZZ] + b1*xp[ZZ];
178     /* 12 Flops */
179
180     /* TOTAL: 63 flops */
181 }
182
183 static void constr_vsite3OUT(rvec xi, rvec xj, rvec xk, rvec x,
184                              real a, real b, real c, t_pbc *pbc)
185 {
186     rvec xij, xik, temp;
187
188     pbc_rvec_sub(pbc, xj, xi, xij);
189     pbc_rvec_sub(pbc, xk, xi, xik);
190     cprod(xij, xik, temp);
191     /* 15 Flops */
192
193     x[XX] = xi[XX] + a*xij[XX] + b*xik[XX] + c*temp[XX];
194     x[YY] = xi[YY] + a*xij[YY] + b*xik[YY] + c*temp[YY];
195     x[ZZ] = xi[ZZ] + a*xij[ZZ] + b*xik[ZZ] + c*temp[ZZ];
196     /* 18 Flops */
197
198     /* TOTAL: 33 flops */
199 }
200
201 static void constr_vsite4FD(rvec xi, rvec xj, rvec xk, rvec xl, rvec x,
202                             real a, real b, real c, t_pbc *pbc)
203 {
204     rvec xij, xjk, xjl, temp;
205     real d;
206
207     pbc_rvec_sub(pbc, xj, xi, xij);
208     pbc_rvec_sub(pbc, xk, xj, xjk);
209     pbc_rvec_sub(pbc, xl, xj, xjl);
210     /* 9 flops */
211
212     /* temp goes from i to a point on the plane jkl */
213     temp[XX] = xij[XX] + a*xjk[XX] + b*xjl[XX];
214     temp[YY] = xij[YY] + a*xjk[YY] + b*xjl[YY];
215     temp[ZZ] = xij[ZZ] + a*xjk[ZZ] + b*xjl[ZZ];
216     /* 12 flops */
217
218     d = c*gmx_invsqrt(iprod(temp, temp));
219     /* 6 + 10 flops */
220
221     x[XX] = xi[XX] + d*temp[XX];
222     x[YY] = xi[YY] + d*temp[YY];
223     x[ZZ] = xi[ZZ] + d*temp[ZZ];
224     /* 6 Flops */
225
226     /* TOTAL: 43 flops */
227 }
228
229
230 static void constr_vsite4FDN(rvec xi, rvec xj, rvec xk, rvec xl, rvec x,
231                              real a, real b, real c, t_pbc *pbc)
232 {
233     rvec xij, xik, xil, ra, rb, rja, rjb, rm;
234     real d;
235
236     pbc_rvec_sub(pbc, xj, xi, xij);
237     pbc_rvec_sub(pbc, xk, xi, xik);
238     pbc_rvec_sub(pbc, xl, xi, xil);
239     /* 9 flops */
240
241     ra[XX] = a*xik[XX];
242     ra[YY] = a*xik[YY];
243     ra[ZZ] = a*xik[ZZ];
244
245     rb[XX] = b*xil[XX];
246     rb[YY] = b*xil[YY];
247     rb[ZZ] = b*xil[ZZ];
248
249     /* 6 flops */
250
251     rvec_sub(ra, xij, rja);
252     rvec_sub(rb, xij, rjb);
253     /* 6 flops */
254
255     cprod(rja, rjb, rm);
256     /* 9 flops */
257
258     d = c*gmx_invsqrt(norm2(rm));
259     /* 5+5+1 flops */
260
261     x[XX] = xi[XX] + d*rm[XX];
262     x[YY] = xi[YY] + d*rm[YY];
263     x[ZZ] = xi[ZZ] + d*rm[ZZ];
264     /* 6 Flops */
265
266     /* TOTAL: 47 flops */
267 }
268
269
270 static int constr_vsiten(t_iatom *ia, t_iparams ip[],
271                          rvec *x, t_pbc *pbc)
272 {
273     rvec xs, x1, dx;
274     dvec dsum;
275     int  n3, av, ai, i;
276     real a;
277
278     n3 = 3*ip[ia[0]].vsiten.n;
279     av = ia[1];
280     ai = ia[2];
281     copy_rvec(x[ai], x1);
282     clear_dvec(dsum);
283     for (i = 3; i < n3; i += 3)
284     {
285         ai = ia[i+2];
286         a  = ip[ia[i]].vsiten.a;
287         if (pbc)
288         {
289             pbc_dx_aiuc(pbc, x[ai], x1, dx);
290         }
291         else
292         {
293             rvec_sub(x[ai], x1, dx);
294         }
295         dsum[XX] += a*dx[XX];
296         dsum[YY] += a*dx[YY];
297         dsum[ZZ] += a*dx[ZZ];
298         /* 9 Flops */
299     }
300
301     x[av][XX] = x1[XX] + dsum[XX];
302     x[av][YY] = x1[YY] + dsum[YY];
303     x[av][ZZ] = x1[ZZ] + dsum[ZZ];
304
305     return n3;
306 }
307
308
309 void construct_vsites_thread(gmx_vsite_t *vsite,
310                              rvec x[],
311                              real dt, rvec *v,
312                              t_iparams ip[], t_ilist ilist[],
313                              t_pbc *pbc_null)
314 {
315     gmx_bool   bPBCAll;
316     rvec       xpbc, xv, vv, dx;
317     real       a1, b1, c1, inv_dt;
318     int        i, inc, ii, nra, nr, tp, ftype;
319     t_iatom    avsite, ai, aj, ak, al, pbc_atom;
320     t_iatom   *ia;
321     t_pbc     *pbc_null2;
322     int       *vsite_pbc, ishift;
323     rvec       reftmp, vtmp, rtmp;
324
325     if (v != NULL)
326     {
327         inv_dt = 1.0/dt;
328     }
329     else
330     {
331         inv_dt = 1.0;
332     }
333
334     bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
335
336     pbc_null2 = NULL;
337     vsite_pbc = NULL;
338     for (ftype = 0; (ftype < F_NRE); ftype++)
339     {
340         if ((interaction_function[ftype].flags & IF_VSITE) &&
341             ilist[ftype].nr > 0)
342         {
343             nra    = interaction_function[ftype].nratoms;
344             inc    = 1 + nra;
345             nr     = ilist[ftype].nr;
346             ia     = ilist[ftype].iatoms;
347
348             if (bPBCAll)
349             {
350                 pbc_null2 = pbc_null;
351             }
352             else if (pbc_null != NULL)
353             {
354                 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
355             }
356
357             for (i = 0; i < nr; )
358             {
359                 tp   = ia[0];
360
361                 /* The vsite and constructing atoms */
362                 avsite = ia[1];
363                 ai     = ia[2];
364
365                 /* Constants for constructing vsites */
366                 a1   = ip[tp].vsite.a;
367                 /* Check what kind of pbc we need to use */
368                 if (bPBCAll)
369                 {
370                     /* No charge groups, vsite follows its own pbc */
371                     pbc_atom = avsite;
372                     copy_rvec(x[avsite], xpbc);
373                 }
374                 else if (vsite_pbc != NULL)
375                 {
376                     pbc_atom = vsite_pbc[i/(1+nra)];
377                     if (pbc_atom > -2)
378                     {
379                         if (pbc_atom >= 0)
380                         {
381                             /* We need to copy the coordinates here,
382                              * single for single atom cg's pbc_atom
383                              * is the vsite itself.
384                              */
385                             copy_rvec(x[pbc_atom], xpbc);
386                         }
387                         pbc_null2 = pbc_null;
388                     }
389                     else
390                     {
391                         pbc_null2 = NULL;
392                     }
393                 }
394                 else
395                 {
396                     pbc_atom = -2;
397                 }
398                 /* Copy the old position */
399                 copy_rvec(x[avsite], xv);
400
401                 /* Construct the vsite depending on type */
402                 switch (ftype)
403                 {
404                     case F_VSITE2:
405                         aj = ia[3];
406                         constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
407                         break;
408                     case F_VSITE3:
409                         aj = ia[3];
410                         ak = ia[4];
411                         b1 = ip[tp].vsite.b;
412                         constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
413                         break;
414                     case F_VSITE3FD:
415                         aj = ia[3];
416                         ak = ia[4];
417                         b1 = ip[tp].vsite.b;
418                         constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
419                         break;
420                     case F_VSITE3FAD:
421                         aj = ia[3];
422                         ak = ia[4];
423                         b1 = ip[tp].vsite.b;
424                         constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
425                         break;
426                     case F_VSITE3OUT:
427                         aj = ia[3];
428                         ak = ia[4];
429                         b1 = ip[tp].vsite.b;
430                         c1 = ip[tp].vsite.c;
431                         constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
432                         break;
433                     case F_VSITE4FD:
434                         aj = ia[3];
435                         ak = ia[4];
436                         al = ia[5];
437                         b1 = ip[tp].vsite.b;
438                         c1 = ip[tp].vsite.c;
439                         constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
440                                         pbc_null2);
441                         break;
442                     case F_VSITE4FDN:
443                         aj = ia[3];
444                         ak = ia[4];
445                         al = ia[5];
446                         b1 = ip[tp].vsite.b;
447                         c1 = ip[tp].vsite.c;
448                         constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
449                                          pbc_null2);
450                         break;
451                     case F_VSITEN:
452                         inc = constr_vsiten(ia, ip, x, pbc_null2);
453                         break;
454                     default:
455                         gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
456                                   ftype, __FILE__, __LINE__);
457                 }
458
459                 if (pbc_atom >= 0)
460                 {
461                     /* Match the pbc of this vsite to the rest of its charge group */
462                     ishift = pbc_dx_aiuc(pbc_null, x[avsite], xpbc, dx);
463                     if (ishift != CENTRAL)
464                     {
465                         rvec_add(xpbc, dx, x[avsite]);
466                     }
467                 }
468                 if (v != NULL)
469                 {
470                     /* Calculate velocity of vsite... */
471                     rvec_sub(x[avsite], xv, vv);
472                     svmul(inv_dt, vv, v[avsite]);
473                 }
474
475                 /* Increment loop variables */
476                 i  += inc;
477                 ia += inc;
478             }
479         }
480     }
481 }
482
483 void construct_vsites(gmx_vsite_t *vsite,
484                       rvec x[],
485                       real dt, rvec *v,
486                       t_iparams ip[], t_ilist ilist[],
487                       int ePBC, gmx_bool bMolPBC,
488                       t_commrec *cr, matrix box)
489 {
490     t_pbc     pbc, *pbc_null;
491     gmx_bool  bDomDec;
492     int       nthreads;
493
494     bDomDec = cr && DOMAINDECOMP(cr);
495
496     /* We only need to do pbc when we have inter-cg vsites */
497     if (ePBC != epbcNONE && (bDomDec || bMolPBC) && vsite->n_intercg_vsite)
498     {
499         /* This is wasting some CPU time as we now do this multiple times
500          * per MD step. But how often do we have vsites with full pbc?
501          */
502         pbc_null = set_pbc_dd(&pbc, ePBC, cr != NULL ? cr->dd : NULL, FALSE, box);
503     }
504     else
505     {
506         pbc_null = NULL;
507     }
508
509     if (bDomDec)
510     {
511         dd_move_x_vsites(cr->dd, box, x);
512     }
513
514     if (vsite->nthreads == 1)
515     {
516         construct_vsites_thread(vsite,
517                                 x, dt, v,
518                                 ip, ilist,
519                                 pbc_null);
520     }
521     else
522     {
523 #pragma omp parallel num_threads(vsite->nthreads)
524         {
525             construct_vsites_thread(vsite,
526                                     x, dt, v,
527                                     ip, vsite->tdata[gmx_omp_get_thread_num()].ilist,
528                                     pbc_null);
529         }
530         /* Now we can construct the vsites that might depend on other vsites */
531         construct_vsites_thread(vsite,
532                                 x, dt, v,
533                                 ip, vsite->tdata[vsite->nthreads].ilist,
534                                 pbc_null);
535     }
536 }
537
538 static void spread_vsite2(t_iatom ia[], real a,
539                           rvec x[], rvec f[], rvec fshift[],
540                           t_pbc *pbc, t_graph *g)
541 {
542     rvec    fi, fj, dx;
543     t_iatom av, ai, aj;
544     ivec    di;
545     real    b;
546     int     siv, sij;
547
548     av = ia[1];
549     ai = ia[2];
550     aj = ia[3];
551
552     svmul(1-a, f[av], fi);
553     svmul(  a, f[av], fj);
554     /* 7 flop */
555
556     rvec_inc(f[ai], fi);
557     rvec_inc(f[aj], fj);
558     /* 6 Flops */
559
560     if (g)
561     {
562         ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
563         siv = IVEC2IS(di);
564         ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
565         sij = IVEC2IS(di);
566     }
567     else if (pbc)
568     {
569         siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
570         sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
571     }
572     else
573     {
574         siv = CENTRAL;
575         sij = CENTRAL;
576     }
577
578     if (fshift && (siv != CENTRAL || sij != CENTRAL))
579     {
580         rvec_inc(fshift[siv], f[av]);
581         rvec_dec(fshift[CENTRAL], fi);
582         rvec_dec(fshift[sij], fj);
583     }
584
585     /* TOTAL: 13 flops */
586 }
587
588 void construct_vsites_mtop(gmx_vsite_t *vsite,
589                            gmx_mtop_t *mtop, rvec x[])
590 {
591     int             as, mb, mol;
592     gmx_molblock_t *molb;
593     gmx_moltype_t  *molt;
594
595     as = 0;
596     for (mb = 0; mb < mtop->nmolblock; mb++)
597     {
598         molb = &mtop->molblock[mb];
599         molt = &mtop->moltype[molb->type];
600         for (mol = 0; mol < molb->nmol; mol++)
601         {
602             construct_vsites(vsite, x+as, 0.0, NULL,
603                              mtop->ffparams.iparams, molt->ilist,
604                              epbcNONE, TRUE, NULL, NULL);
605             as += molt->atoms.nr;
606         }
607     }
608 }
609
610 static void spread_vsite3(t_iatom ia[], real a, real b,
611                           rvec x[], rvec f[], rvec fshift[],
612                           t_pbc *pbc, t_graph *g)
613 {
614     rvec    fi, fj, fk, dx;
615     atom_id av, ai, aj, ak;
616     ivec    di;
617     int     siv, sij, sik;
618
619     av = ia[1];
620     ai = ia[2];
621     aj = ia[3];
622     ak = ia[4];
623
624     svmul(1-a-b, f[av], fi);
625     svmul(    a, f[av], fj);
626     svmul(    b, f[av], fk);
627     /* 11 flops */
628
629     rvec_inc(f[ai], fi);
630     rvec_inc(f[aj], fj);
631     rvec_inc(f[ak], fk);
632     /* 9 Flops */
633
634     if (g)
635     {
636         ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ia[1]), di);
637         siv = IVEC2IS(di);
638         ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
639         sij = IVEC2IS(di);
640         ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), di);
641         sik = IVEC2IS(di);
642     }
643     else if (pbc)
644     {
645         siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
646         sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
647         sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
648     }
649     else
650     {
651         siv = CENTRAL;
652         sij = CENTRAL;
653         sik = CENTRAL;
654     }
655
656     if (fshift && (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL))
657     {
658         rvec_inc(fshift[siv], f[av]);
659         rvec_dec(fshift[CENTRAL], fi);
660         rvec_dec(fshift[sij], fj);
661         rvec_dec(fshift[sik], fk);
662     }
663
664     /* TOTAL: 20 flops */
665 }
666
667 static void spread_vsite3FD(t_iatom ia[], real a, real b,
668                             rvec x[], rvec f[], rvec fshift[],
669                             gmx_bool VirCorr, matrix dxdf,
670                             t_pbc *pbc, t_graph *g)
671 {
672     real    fx, fy, fz, c, invl, fproj, a1;
673     rvec    xvi, xij, xjk, xix, fv, temp;
674     t_iatom av, ai, aj, ak;
675     int     svi, sji, skj, d;
676     ivec    di;
677
678     av = ia[1];
679     ai = ia[2];
680     aj = ia[3];
681     ak = ia[4];
682     copy_rvec(f[av], fv);
683
684     sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
685     skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
686     /* 6 flops */
687
688     /* xix goes from i to point x on the line jk */
689     xix[XX] = xij[XX]+a*xjk[XX];
690     xix[YY] = xij[YY]+a*xjk[YY];
691     xix[ZZ] = xij[ZZ]+a*xjk[ZZ];
692     /* 6 flops */
693
694     invl = gmx_invsqrt(iprod(xix, xix));
695     c    = b*invl;
696     /* 4 + ?10? flops */
697
698     fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
699
700     temp[XX] = c*(fv[XX]-fproj*xix[XX]);
701     temp[YY] = c*(fv[YY]-fproj*xix[YY]);
702     temp[ZZ] = c*(fv[ZZ]-fproj*xix[ZZ]);
703     /* 16 */
704
705     /* c is already calculated in constr_vsite3FD
706        storing c somewhere will save 26 flops!     */
707
708     a1         = 1-a;
709     f[ai][XX] += fv[XX] - temp[XX];
710     f[ai][YY] += fv[YY] - temp[YY];
711     f[ai][ZZ] += fv[ZZ] - temp[ZZ];
712     f[aj][XX] += a1*temp[XX];
713     f[aj][YY] += a1*temp[YY];
714     f[aj][ZZ] += a1*temp[ZZ];
715     f[ak][XX] += a*temp[XX];
716     f[ak][YY] += a*temp[YY];
717     f[ak][ZZ] += a*temp[ZZ];
718     /* 19 Flops */
719
720     if (g)
721     {
722         ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
723         svi = IVEC2IS(di);
724         ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
725         sji = IVEC2IS(di);
726         ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
727         skj = IVEC2IS(di);
728     }
729     else if (pbc)
730     {
731         svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
732     }
733     else
734     {
735         svi = CENTRAL;
736     }
737
738     if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
739     {
740         rvec_dec(fshift[svi], fv);
741         fshift[CENTRAL][XX] += fv[XX] - (1 + a)*temp[XX];
742         fshift[CENTRAL][YY] += fv[YY] - (1 + a)*temp[YY];
743         fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a)*temp[ZZ];
744         fshift[    sji][XX] += temp[XX];
745         fshift[    sji][YY] += temp[YY];
746         fshift[    sji][ZZ] += temp[ZZ];
747         fshift[    skj][XX] += a*temp[XX];
748         fshift[    skj][YY] += a*temp[YY];
749         fshift[    skj][ZZ] += a*temp[ZZ];
750     }
751
752     if (VirCorr)
753     {
754         /* When VirCorr=TRUE, the virial for the current forces is not
755          * calculated from the redistributed forces. This means that
756          * the effect of non-linear virtual site constructions on the virial
757          * needs to be added separately. This contribution can be calculated
758          * in many ways, but the simplest and cheapest way is to use
759          * the first constructing atom ai as a reference position in space:
760          * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
761          */
762         rvec xiv;
763         int  i, j;
764
765         pbc_rvec_sub(pbc, x[av], x[ai], xiv);
766
767         for (i = 0; i < DIM; i++)
768         {
769             for (j = 0; j < DIM; j++)
770             {
771                 /* As xix is a linear combination of j and k, use that here */
772                 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
773             }
774         }
775     }
776
777     /* TOTAL: 61 flops */
778 }
779
780 static void spread_vsite3FAD(t_iatom ia[], real a, real b,
781                              rvec x[], rvec f[], rvec fshift[],
782                              gmx_bool VirCorr, matrix dxdf,
783                              t_pbc *pbc, t_graph *g)
784 {
785     rvec    xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
786     real    a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
787     t_iatom av, ai, aj, ak;
788     int     svi, sji, skj, d;
789     ivec    di;
790
791     av = ia[1];
792     ai = ia[2];
793     aj = ia[3];
794     ak = ia[4];
795     copy_rvec(f[ia[1]], fv);
796
797     sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
798     skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
799     /* 6 flops */
800
801     invdij    = gmx_invsqrt(iprod(xij, xij));
802     invdij2   = invdij * invdij;
803     c1        = iprod(xij, xjk) * invdij2;
804     xperp[XX] = xjk[XX] - c1*xij[XX];
805     xperp[YY] = xjk[YY] - c1*xij[YY];
806     xperp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
807     /* xperp in plane ijk, perp. to ij */
808     invdp = gmx_invsqrt(iprod(xperp, xperp));
809     a1    = a*invdij;
810     b1    = b*invdp;
811     /* 45 flops */
812
813     /* a1, b1 and c1 are already calculated in constr_vsite3FAD
814        storing them somewhere will save 45 flops!     */
815
816     fproj = iprod(xij, fv)*invdij2;
817     svmul(fproj,                      xij,  Fpij);    /* proj. f on xij */
818     svmul(iprod(xperp, fv)*invdp*invdp, xperp, Fppp); /* proj. f on xperp */
819     svmul(b1*fproj,                   xperp, f3);
820     /* 23 flops */
821
822     rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
823     rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
824     for (d = 0; (d < DIM); d++)
825     {
826         f1[d] *= a1;
827         f2[d] *= b1;
828     }
829     /* 12 flops */
830
831     c2         = 1+c1;
832     f[ai][XX] += fv[XX] - f1[XX] + c1*f2[XX] + f3[XX];
833     f[ai][YY] += fv[YY] - f1[YY] + c1*f2[YY] + f3[YY];
834     f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1*f2[ZZ] + f3[ZZ];
835     f[aj][XX] +=          f1[XX] - c2*f2[XX] - f3[XX];
836     f[aj][YY] +=          f1[YY] - c2*f2[YY] - f3[YY];
837     f[aj][ZZ] +=          f1[ZZ] - c2*f2[ZZ] - f3[ZZ];
838     f[ak][XX] +=                      f2[XX];
839     f[ak][YY] +=                      f2[YY];
840     f[ak][ZZ] +=                      f2[ZZ];
841     /* 30 Flops */
842
843     if (g)
844     {
845         ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
846         svi = IVEC2IS(di);
847         ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
848         sji = IVEC2IS(di);
849         ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
850         skj = IVEC2IS(di);
851     }
852     else if (pbc)
853     {
854         svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
855     }
856     else
857     {
858         svi = CENTRAL;
859     }
860
861     if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
862     {
863         rvec_dec(fshift[svi], fv);
864         fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1-c1)*f2[XX] + f3[XX];
865         fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1-c1)*f2[YY] + f3[YY];
866         fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1-c1)*f2[ZZ] + f3[ZZ];
867         fshift[    sji][XX] +=          f1[XX] -    c1 *f2[XX] - f3[XX];
868         fshift[    sji][YY] +=          f1[YY] -    c1 *f2[YY] - f3[YY];
869         fshift[    sji][ZZ] +=          f1[ZZ] -    c1 *f2[ZZ] - f3[ZZ];
870         fshift[    skj][XX] +=                          f2[XX];
871         fshift[    skj][YY] +=                          f2[YY];
872         fshift[    skj][ZZ] +=                          f2[ZZ];
873     }
874
875     if (VirCorr)
876     {
877         rvec xiv;
878         int  i, j;
879
880         pbc_rvec_sub(pbc, x[av], x[ai], xiv);
881
882         for (i = 0; i < DIM; i++)
883         {
884             for (j = 0; j < DIM; j++)
885             {
886                 /* Note that xik=xij+xjk, so we have to add xij*f2 */
887                 dxdf[i][j] +=
888                     -xiv[i]*fv[j]
889                     + xij[i]*(f1[j] + (1 - c2)*f2[j] - f3[j])
890                     + xjk[i]*f2[j];
891             }
892         }
893     }
894
895     /* TOTAL: 113 flops */
896 }
897
898 static void spread_vsite3OUT(t_iatom ia[], real a, real b, real c,
899                              rvec x[], rvec f[], rvec fshift[],
900                              gmx_bool VirCorr, matrix dxdf,
901                              t_pbc *pbc, t_graph *g)
902 {
903     rvec    xvi, xij, xik, fv, fj, fk;
904     real    cfx, cfy, cfz;
905     atom_id av, ai, aj, ak;
906     ivec    di;
907     int     svi, sji, ski;
908
909     av = ia[1];
910     ai = ia[2];
911     aj = ia[3];
912     ak = ia[4];
913
914     sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
915     ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
916     /* 6 Flops */
917
918     copy_rvec(f[av], fv);
919
920     cfx = c*fv[XX];
921     cfy = c*fv[YY];
922     cfz = c*fv[ZZ];
923     /* 3 Flops */
924
925     fj[XX] = a*fv[XX]     -  xik[ZZ]*cfy +  xik[YY]*cfz;
926     fj[YY] =  xik[ZZ]*cfx + a*fv[YY]     -  xik[XX]*cfz;
927     fj[ZZ] = -xik[YY]*cfx +  xik[XX]*cfy + a*fv[ZZ];
928
929     fk[XX] = b*fv[XX]     +  xij[ZZ]*cfy -  xij[YY]*cfz;
930     fk[YY] = -xij[ZZ]*cfx + b*fv[YY]     +  xij[XX]*cfz;
931     fk[ZZ] =  xij[YY]*cfx -  xij[XX]*cfy + b*fv[ZZ];
932     /* 30 Flops */
933
934     f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
935     f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
936     f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
937     rvec_inc(f[aj], fj);
938     rvec_inc(f[ak], fk);
939     /* 15 Flops */
940
941     if (g)
942     {
943         ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
944         svi = IVEC2IS(di);
945         ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
946         sji = IVEC2IS(di);
947         ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
948         ski = IVEC2IS(di);
949     }
950     else if (pbc)
951     {
952         svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
953     }
954     else
955     {
956         svi = CENTRAL;
957     }
958
959     if (fshift && (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL))
960     {
961         rvec_dec(fshift[svi], fv);
962         fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
963         fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
964         fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
965         rvec_inc(fshift[sji], fj);
966         rvec_inc(fshift[ski], fk);
967     }
968
969     if (VirCorr)
970     {
971         rvec xiv;
972         int  i, j;
973
974         pbc_rvec_sub(pbc, x[av], x[ai], xiv);
975
976         for (i = 0; i < DIM; i++)
977         {
978             for (j = 0; j < DIM; j++)
979             {
980                 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j];
981             }
982         }
983     }
984
985     /* TOTAL: 54 flops */
986 }
987
988 static void spread_vsite4FD(t_iatom ia[], real a, real b, real c,
989                             rvec x[], rvec f[], rvec fshift[],
990                             gmx_bool VirCorr, matrix dxdf,
991                             t_pbc *pbc, t_graph *g)
992 {
993     real    d, invl, fproj, a1;
994     rvec    xvi, xij, xjk, xjl, xix, fv, temp;
995     atom_id av, ai, aj, ak, al;
996     ivec    di;
997     int     svi, sji, skj, slj, m;
998
999     av = ia[1];
1000     ai = ia[2];
1001     aj = ia[3];
1002     ak = ia[4];
1003     al = ia[5];
1004
1005     sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1006     skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1007     slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1008     /* 9 flops */
1009
1010     /* xix goes from i to point x on the plane jkl */
1011     for (m = 0; m < DIM; m++)
1012     {
1013         xix[m] = xij[m] + a*xjk[m] + b*xjl[m];
1014     }
1015     /* 12 flops */
1016
1017     invl = gmx_invsqrt(iprod(xix, xix));
1018     d    = c*invl;
1019     /* 4 + ?10? flops */
1020
1021     copy_rvec(f[av], fv);
1022
1023     fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
1024
1025     for (m = 0; m < DIM; m++)
1026     {
1027         temp[m] = d*(fv[m] - fproj*xix[m]);
1028     }
1029     /* 16 */
1030
1031     /* c is already calculated in constr_vsite3FD
1032        storing c somewhere will save 35 flops!     */
1033
1034     a1 = 1 - a - b;
1035     for (m = 0; m < DIM; m++)
1036     {
1037         f[ai][m] += fv[m] - temp[m];
1038         f[aj][m] += a1*temp[m];
1039         f[ak][m] += a*temp[m];
1040         f[al][m] += b*temp[m];
1041     }
1042     /* 26 Flops */
1043
1044     if (g)
1045     {
1046         ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1047         svi = IVEC2IS(di);
1048         ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1049         sji = IVEC2IS(di);
1050         ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
1051         skj = IVEC2IS(di);
1052         ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, aj), di);
1053         slj = IVEC2IS(di);
1054     }
1055     else if (pbc)
1056     {
1057         svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1058     }
1059     else
1060     {
1061         svi = CENTRAL;
1062     }
1063
1064     if (fshift &&
1065         (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL))
1066     {
1067         rvec_dec(fshift[svi], fv);
1068         for (m = 0; m < DIM; m++)
1069         {
1070             fshift[CENTRAL][m] += fv[m] - (1 + a + b)*temp[m];
1071             fshift[    sji][m] += temp[m];
1072             fshift[    skj][m] += a*temp[m];
1073             fshift[    slj][m] += b*temp[m];
1074         }
1075     }
1076
1077     if (VirCorr)
1078     {
1079         rvec xiv;
1080         int  i, j;
1081
1082         pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1083
1084         for (i = 0; i < DIM; i++)
1085         {
1086             for (j = 0; j < DIM; j++)
1087             {
1088                 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
1089             }
1090         }
1091     }
1092
1093     /* TOTAL: 77 flops */
1094 }
1095
1096
1097 static void spread_vsite4FDN(t_iatom ia[], real a, real b, real c,
1098                              rvec x[], rvec f[], rvec fshift[],
1099                              gmx_bool VirCorr, matrix dxdf,
1100                              t_pbc *pbc, t_graph *g)
1101 {
1102     rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1103     rvec fv, fj, fk, fl;
1104     real invrm, denom;
1105     real cfx, cfy, cfz;
1106     ivec di;
1107     int  av, ai, aj, ak, al;
1108     int  svi, sij, sik, sil;
1109
1110     /* DEBUG: check atom indices */
1111     av = ia[1];
1112     ai = ia[2];
1113     aj = ia[3];
1114     ak = ia[4];
1115     al = ia[5];
1116
1117     copy_rvec(f[av], fv);
1118
1119     sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1120     sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1121     sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1122     /* 9 flops */
1123
1124     ra[XX] = a*xik[XX];
1125     ra[YY] = a*xik[YY];
1126     ra[ZZ] = a*xik[ZZ];
1127
1128     rb[XX] = b*xil[XX];
1129     rb[YY] = b*xil[YY];
1130     rb[ZZ] = b*xil[ZZ];
1131
1132     /* 6 flops */
1133
1134     rvec_sub(ra, xij, rja);
1135     rvec_sub(rb, xij, rjb);
1136     rvec_sub(rb, ra, rab);
1137     /* 9 flops */
1138
1139     cprod(rja, rjb, rm);
1140     /* 9 flops */
1141
1142     invrm = gmx_invsqrt(norm2(rm));
1143     denom = invrm*invrm;
1144     /* 5+5+2 flops */
1145
1146     cfx = c*invrm*fv[XX];
1147     cfy = c*invrm*fv[YY];
1148     cfz = c*invrm*fv[ZZ];
1149     /* 6 Flops */
1150
1151     cprod(rm, rab, rt);
1152     /* 9 flops */
1153
1154     rt[XX] *= denom;
1155     rt[YY] *= denom;
1156     rt[ZZ] *= denom;
1157     /* 3flops */
1158
1159     fj[XX] = (        -rm[XX]*rt[XX]) * cfx + ( rab[ZZ]-rm[YY]*rt[XX]) * cfy + (-rab[YY]-rm[ZZ]*rt[XX]) * cfz;
1160     fj[YY] = (-rab[ZZ]-rm[XX]*rt[YY]) * cfx + (        -rm[YY]*rt[YY]) * cfy + ( rab[XX]-rm[ZZ]*rt[YY]) * cfz;
1161     fj[ZZ] = ( rab[YY]-rm[XX]*rt[ZZ]) * cfx + (-rab[XX]-rm[YY]*rt[ZZ]) * cfy + (        -rm[ZZ]*rt[ZZ]) * cfz;
1162     /* 30 flops */
1163
1164     cprod(rjb, rm, rt);
1165     /* 9 flops */
1166
1167     rt[XX] *= denom*a;
1168     rt[YY] *= denom*a;
1169     rt[ZZ] *= denom*a;
1170     /* 3flops */
1171
1172     fk[XX] = (          -rm[XX]*rt[XX]) * cfx + (-a*rjb[ZZ]-rm[YY]*rt[XX]) * cfy + ( a*rjb[YY]-rm[ZZ]*rt[XX]) * cfz;
1173     fk[YY] = ( a*rjb[ZZ]-rm[XX]*rt[YY]) * cfx + (          -rm[YY]*rt[YY]) * cfy + (-a*rjb[XX]-rm[ZZ]*rt[YY]) * cfz;
1174     fk[ZZ] = (-a*rjb[YY]-rm[XX]*rt[ZZ]) * cfx + ( a*rjb[XX]-rm[YY]*rt[ZZ]) * cfy + (          -rm[ZZ]*rt[ZZ]) * cfz;
1175     /* 36 flops */
1176
1177     cprod(rm, rja, rt);
1178     /* 9 flops */
1179
1180     rt[XX] *= denom*b;
1181     rt[YY] *= denom*b;
1182     rt[ZZ] *= denom*b;
1183     /* 3flops */
1184
1185     fl[XX] = (          -rm[XX]*rt[XX]) * cfx + ( b*rja[ZZ]-rm[YY]*rt[XX]) * cfy + (-b*rja[YY]-rm[ZZ]*rt[XX]) * cfz;
1186     fl[YY] = (-b*rja[ZZ]-rm[XX]*rt[YY]) * cfx + (          -rm[YY]*rt[YY]) * cfy + ( b*rja[XX]-rm[ZZ]*rt[YY]) * cfz;
1187     fl[ZZ] = ( b*rja[YY]-rm[XX]*rt[ZZ]) * cfx + (-b*rja[XX]-rm[YY]*rt[ZZ]) * cfy + (          -rm[ZZ]*rt[ZZ]) * cfz;
1188     /* 36 flops */
1189
1190     f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1191     f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1192     f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1193     rvec_inc(f[aj], fj);
1194     rvec_inc(f[ak], fk);
1195     rvec_inc(f[al], fl);
1196     /* 21 flops */
1197
1198     if (g)
1199     {
1200         ivec_sub(SHIFT_IVEC(g, av), SHIFT_IVEC(g, ai), di);
1201         svi = IVEC2IS(di);
1202         ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1203         sij = IVEC2IS(di);
1204         ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1205         sik = IVEC2IS(di);
1206         ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, ai), di);
1207         sil = IVEC2IS(di);
1208     }
1209     else if (pbc)
1210     {
1211         svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1212     }
1213     else
1214     {
1215         svi = CENTRAL;
1216     }
1217
1218     if (fshift && (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL))
1219     {
1220         rvec_dec(fshift[svi], fv);
1221         fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1222         fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1223         fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1224         rvec_inc(fshift[sij], fj);
1225         rvec_inc(fshift[sik], fk);
1226         rvec_inc(fshift[sil], fl);
1227     }
1228
1229     if (VirCorr)
1230     {
1231         rvec xiv;
1232         int  i, j;
1233
1234         pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1235
1236         for (i = 0; i < DIM; i++)
1237         {
1238             for (j = 0; j < DIM; j++)
1239             {
1240                 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j] + xil[i]*fl[j];
1241             }
1242         }
1243     }
1244
1245     /* Total: 207 flops (Yuck!) */
1246 }
1247
1248
1249 static int spread_vsiten(t_iatom ia[], t_iparams ip[],
1250                          rvec x[], rvec f[], rvec fshift[],
1251                          t_pbc *pbc, t_graph *g)
1252 {
1253     rvec xv, dx, fi;
1254     int  n3, av, i, ai;
1255     real a;
1256     ivec di;
1257     int  siv;
1258
1259     n3 = 3*ip[ia[0]].vsiten.n;
1260     av = ia[1];
1261     copy_rvec(x[av], xv);
1262
1263     for (i = 0; i < n3; i += 3)
1264     {
1265         ai = ia[i+2];
1266         if (g)
1267         {
1268             ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
1269             siv = IVEC2IS(di);
1270         }
1271         else if (pbc)
1272         {
1273             siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1274         }
1275         else
1276         {
1277             siv = CENTRAL;
1278         }
1279         a = ip[ia[i]].vsiten.a;
1280         svmul(a, f[av], fi);
1281         rvec_inc(f[ai], fi);
1282         if (fshift && siv != CENTRAL)
1283         {
1284             rvec_inc(fshift[siv], fi);
1285             rvec_dec(fshift[CENTRAL], fi);
1286         }
1287         /* 6 Flops */
1288     }
1289
1290     return n3;
1291 }
1292
1293
1294 static int vsite_count(const t_ilist *ilist, int ftype)
1295 {
1296     if (ftype == F_VSITEN)
1297     {
1298         return ilist[ftype].nr/3;
1299     }
1300     else
1301     {
1302         return ilist[ftype].nr/(1 + interaction_function[ftype].nratoms);
1303     }
1304 }
1305
1306 static void spread_vsite_f_thread(gmx_vsite_t *vsite,
1307                                   rvec x[], rvec f[], rvec *fshift,
1308                                   gmx_bool VirCorr, matrix dxdf,
1309                                   t_iparams ip[], t_ilist ilist[],
1310                                   t_graph *g, t_pbc *pbc_null)
1311 {
1312     gmx_bool   bPBCAll;
1313     real       a1, b1, c1;
1314     int        i, inc, m, nra, nr, tp, ftype;
1315     t_iatom   *ia;
1316     t_pbc     *pbc_null2;
1317     int       *vsite_pbc;
1318
1319     if (VirCorr)
1320     {
1321         clear_mat(dxdf);
1322     }
1323
1324     bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
1325
1326     /* this loop goes backwards to be able to build *
1327      * higher type vsites from lower types         */
1328     pbc_null2 = NULL;
1329     vsite_pbc = NULL;
1330     for (ftype = F_NRE-1; (ftype >= 0); ftype--)
1331     {
1332         if ((interaction_function[ftype].flags & IF_VSITE) &&
1333             ilist[ftype].nr > 0)
1334         {
1335             nra    = interaction_function[ftype].nratoms;
1336             inc    = 1 + nra;
1337             nr     = ilist[ftype].nr;
1338             ia     = ilist[ftype].iatoms;
1339
1340             if (bPBCAll)
1341             {
1342                 pbc_null2 = pbc_null;
1343             }
1344             else if (pbc_null != NULL)
1345             {
1346                 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
1347             }
1348
1349             for (i = 0; i < nr; )
1350             {
1351                 if (vsite_pbc != NULL)
1352                 {
1353                     if (vsite_pbc[i/(1+nra)] > -2)
1354                     {
1355                         pbc_null2 = pbc_null;
1356                     }
1357                     else
1358                     {
1359                         pbc_null2 = NULL;
1360                     }
1361                 }
1362
1363                 tp   = ia[0];
1364
1365                 /* Constants for constructing */
1366                 a1   = ip[tp].vsite.a;
1367                 /* Construct the vsite depending on type */
1368                 switch (ftype)
1369                 {
1370                     case F_VSITE2:
1371                         spread_vsite2(ia, a1, x, f, fshift, pbc_null2, g);
1372                         break;
1373                     case F_VSITE3:
1374                         b1 = ip[tp].vsite.b;
1375                         spread_vsite3(ia, a1, b1, x, f, fshift, pbc_null2, g);
1376                         break;
1377                     case F_VSITE3FD:
1378                         b1 = ip[tp].vsite.b;
1379                         spread_vsite3FD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1380                         break;
1381                     case F_VSITE3FAD:
1382                         b1 = ip[tp].vsite.b;
1383                         spread_vsite3FAD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1384                         break;
1385                     case F_VSITE3OUT:
1386                         b1 = ip[tp].vsite.b;
1387                         c1 = ip[tp].vsite.c;
1388                         spread_vsite3OUT(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1389                         break;
1390                     case F_VSITE4FD:
1391                         b1 = ip[tp].vsite.b;
1392                         c1 = ip[tp].vsite.c;
1393                         spread_vsite4FD(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1394                         break;
1395                     case F_VSITE4FDN:
1396                         b1 = ip[tp].vsite.b;
1397                         c1 = ip[tp].vsite.c;
1398                         spread_vsite4FDN(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1399                         break;
1400                     case F_VSITEN:
1401                         inc = spread_vsiten(ia, ip, x, f, fshift, pbc_null2, g);
1402                         break;
1403                     default:
1404                         gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
1405                                   ftype, __FILE__, __LINE__);
1406                 }
1407                 clear_rvec(f[ia[1]]);
1408
1409                 /* Increment loop variables */
1410                 i  += inc;
1411                 ia += inc;
1412             }
1413         }
1414     }
1415 }
1416
1417 void spread_vsite_f(gmx_vsite_t *vsite,
1418                     rvec x[], rvec f[], rvec *fshift,
1419                     gmx_bool VirCorr, matrix vir,
1420                     t_nrnb *nrnb, t_idef *idef,
1421                     int ePBC, gmx_bool bMolPBC, t_graph *g, matrix box,
1422                     t_commrec *cr)
1423 {
1424     t_pbc pbc, *pbc_null;
1425     int   th;
1426
1427     /* We only need to do pbc when we have inter-cg vsites */
1428     if ((DOMAINDECOMP(cr) || bMolPBC) && vsite->n_intercg_vsite)
1429     {
1430         /* This is wasting some CPU time as we now do this multiple times
1431          * per MD step. But how often do we have vsites with full pbc?
1432          */
1433         pbc_null = set_pbc_dd(&pbc, ePBC, cr->dd, FALSE, box);
1434     }
1435     else
1436     {
1437         pbc_null = NULL;
1438     }
1439
1440     if (DOMAINDECOMP(cr))
1441     {
1442         dd_clear_f_vsites(cr->dd, f);
1443     }
1444
1445     if (vsite->nthreads == 1)
1446     {
1447         spread_vsite_f_thread(vsite,
1448                               x, f, fshift,
1449                               VirCorr, vsite->tdata[0].dxdf,
1450                               idef->iparams, idef->il,
1451                               g, pbc_null);
1452     }
1453     else
1454     {
1455         /* First spread the vsites that might depend on other vsites */
1456         spread_vsite_f_thread(vsite,
1457                               x, f, fshift,
1458                               VirCorr, vsite->tdata[vsite->nthreads].dxdf,
1459                               idef->iparams,
1460                               vsite->tdata[vsite->nthreads].ilist,
1461                               g, pbc_null);
1462
1463 #pragma omp parallel num_threads(vsite->nthreads)
1464         {
1465             int   thread;
1466             rvec *fshift_t;
1467
1468             thread = gmx_omp_get_thread_num();
1469
1470             if (thread == 0 || fshift == NULL)
1471             {
1472                 fshift_t = fshift;
1473             }
1474             else
1475             {
1476                 int i;
1477
1478                 fshift_t = vsite->tdata[thread].fshift;
1479
1480                 for (i = 0; i < SHIFTS; i++)
1481                 {
1482                     clear_rvec(fshift_t[i]);
1483                 }
1484             }
1485
1486             spread_vsite_f_thread(vsite,
1487                                   x, f, fshift_t,
1488                                   VirCorr, vsite->tdata[thread].dxdf,
1489                                   idef->iparams,
1490                                   vsite->tdata[thread].ilist,
1491                                   g, pbc_null);
1492         }
1493
1494         if (fshift != NULL)
1495         {
1496             int i;
1497
1498             for (th = 1; th < vsite->nthreads; th++)
1499             {
1500                 for (i = 0; i < SHIFTS; i++)
1501                 {
1502                     rvec_inc(fshift[i], vsite->tdata[th].fshift[i]);
1503                 }
1504             }
1505         }
1506     }
1507
1508     if (VirCorr)
1509     {
1510         int i, j;
1511
1512         for (th = 0; th < (vsite->nthreads == 1 ? 1 : vsite->nthreads+1); th++)
1513         {
1514             for (i = 0; i < DIM; i++)
1515             {
1516                 for (j = 0; j < DIM; j++)
1517                 {
1518                     vir[i][j] += -0.5*vsite->tdata[th].dxdf[i][j];
1519                 }
1520             }
1521         }
1522     }
1523
1524     if (DOMAINDECOMP(cr))
1525     {
1526         dd_move_f_vsites(cr->dd, f, fshift);
1527     }
1528
1529     inc_nrnb(nrnb, eNR_VSITE2,   vsite_count(idef->il, F_VSITE2));
1530     inc_nrnb(nrnb, eNR_VSITE3,   vsite_count(idef->il, F_VSITE3));
1531     inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(idef->il, F_VSITE3FD));
1532     inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(idef->il, F_VSITE3FAD));
1533     inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(idef->il, F_VSITE3OUT));
1534     inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(idef->il, F_VSITE4FD));
1535     inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(idef->il, F_VSITE4FDN));
1536     inc_nrnb(nrnb, eNR_VSITEN,   vsite_count(idef->il, F_VSITEN));
1537 }
1538
1539 static int *atom2cg(t_block *cgs)
1540 {
1541     int *a2cg, cg, i;
1542
1543     snew(a2cg, cgs->index[cgs->nr]);
1544     for (cg = 0; cg < cgs->nr; cg++)
1545     {
1546         for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
1547         {
1548             a2cg[i] = cg;
1549         }
1550     }
1551
1552     return a2cg;
1553 }
1554
1555 static int count_intercg_vsite(gmx_mtop_t *mtop,
1556                                gmx_bool   *bHaveChargeGroups)
1557 {
1558     int             mb, mt, ftype, nral, i, cg, a;
1559     gmx_molblock_t *molb;
1560     gmx_moltype_t  *molt;
1561     int            *a2cg;
1562     t_ilist        *il;
1563     t_iatom        *ia;
1564     int             n_intercg_vsite;
1565
1566     *bHaveChargeGroups = FALSE;
1567
1568     n_intercg_vsite = 0;
1569     for (mb = 0; mb < mtop->nmolblock; mb++)
1570     {
1571         molb = &mtop->molblock[mb];
1572         molt = &mtop->moltype[molb->type];
1573
1574         if (molt->cgs.nr < molt->atoms.nr)
1575         {
1576             *bHaveChargeGroups = TRUE;
1577         }
1578
1579         a2cg = atom2cg(&molt->cgs);
1580         for (ftype = 0; ftype < F_NRE; ftype++)
1581         {
1582             if (interaction_function[ftype].flags & IF_VSITE)
1583             {
1584                 nral = NRAL(ftype);
1585                 il   = &molt->ilist[ftype];
1586                 ia   = il->iatoms;
1587                 for (i = 0; i < il->nr; i += 1+nral)
1588                 {
1589                     cg = a2cg[ia[1+i]];
1590                     for (a = 1; a < nral; a++)
1591                     {
1592                         if (a2cg[ia[1+a]] != cg)
1593                         {
1594                             n_intercg_vsite += molb->nmol;
1595                             break;
1596                         }
1597                     }
1598                 }
1599             }
1600         }
1601         sfree(a2cg);
1602     }
1603
1604     return n_intercg_vsite;
1605 }
1606
1607 static int **get_vsite_pbc(t_iparams *iparams, t_ilist *ilist,
1608                            t_atom *atom, t_mdatoms *md,
1609                            t_block *cgs, int *a2cg)
1610 {
1611     int      ftype, nral, i, j, vsi, vsite, cg_v, cg_c, a, nc3 = 0;
1612     t_ilist *il;
1613     t_iatom *ia;
1614     int    **vsite_pbc, *vsite_pbc_f;
1615     char    *pbc_set;
1616     gmx_bool bViteOnlyCG_and_FirstAtom;
1617
1618     /* Make an array that tells if the pbc of an atom is set */
1619     snew(pbc_set, cgs->index[cgs->nr]);
1620     /* PBC is set for all non vsites */
1621     for (a = 0; a < cgs->index[cgs->nr]; a++)
1622     {
1623         if ((atom && atom[a].ptype != eptVSite) ||
1624             (md   && md->ptype[a]  != eptVSite))
1625         {
1626             pbc_set[a] = 1;
1627         }
1628     }
1629
1630     snew(vsite_pbc, F_VSITEN-F_VSITE2+1);
1631
1632     for (ftype = 0; ftype < F_NRE; ftype++)
1633     {
1634         if (interaction_function[ftype].flags & IF_VSITE)
1635         {
1636             nral = NRAL(ftype);
1637             il   = &ilist[ftype];
1638             ia   = il->iatoms;
1639
1640             snew(vsite_pbc[ftype-F_VSITE2], il->nr/(1+nral));
1641             vsite_pbc_f = vsite_pbc[ftype-F_VSITE2];
1642
1643             i = 0;
1644             while (i < il->nr)
1645             {
1646                 vsi   = i/(1+nral);
1647                 vsite = ia[i+1];
1648                 cg_v  = a2cg[vsite];
1649                 /* A value of -2 signals that this vsite and its contructing
1650                  * atoms are all within the same cg, so no pbc is required.
1651                  */
1652                 vsite_pbc_f[vsi] = -2;
1653                 /* Check if constructing atoms are outside the vsite's cg */
1654                 nc3 = 0;
1655                 if (ftype == F_VSITEN)
1656                 {
1657                     nc3 = 3*iparams[ia[i]].vsiten.n;
1658                     for (j = 0; j < nc3; j += 3)
1659                     {
1660                         if (a2cg[ia[i+j+2]] != cg_v)
1661                         {
1662                             vsite_pbc_f[vsi] = -1;
1663                         }
1664                     }
1665                 }
1666                 else
1667                 {
1668                     for (a = 1; a < nral; a++)
1669                     {
1670                         if (a2cg[ia[i+1+a]] != cg_v)
1671                         {
1672                             vsite_pbc_f[vsi] = -1;
1673                         }
1674                     }
1675                 }
1676                 if (vsite_pbc_f[vsi] == -1)
1677                 {
1678                     /* Check if this is the first processed atom of a vsite only cg */
1679                     bViteOnlyCG_and_FirstAtom = TRUE;
1680                     for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1681                     {
1682                         /* Non-vsites already have pbc set, so simply check for pbc_set */
1683                         if (pbc_set[a])
1684                         {
1685                             bViteOnlyCG_and_FirstAtom = FALSE;
1686                             break;
1687                         }
1688                     }
1689                     if (bViteOnlyCG_and_FirstAtom)
1690                     {
1691                         /* First processed atom of a vsite only charge group.
1692                          * The pbc of the input coordinates to construct_vsites
1693                          * should be preserved.
1694                          */
1695                         vsite_pbc_f[vsi] = vsite;
1696                     }
1697                     else if (cg_v != a2cg[ia[1+i+1]])
1698                     {
1699                         /* This vsite has a different charge group index
1700                          * than it's first constructing atom
1701                          * and the charge group has more than one atom,
1702                          * search for the first normal particle
1703                          * or vsite that already had its pbc defined.
1704                          * If nothing is found, use full pbc for this vsite.
1705                          */
1706                         for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1707                         {
1708                             if (a != vsite && pbc_set[a])
1709                             {
1710                                 vsite_pbc_f[vsi] = a;
1711                                 if (gmx_debug_at)
1712                                 {
1713                                     fprintf(debug, "vsite %d match pbc with atom %d\n",
1714                                             vsite+1, a+1);
1715                                 }
1716                                 break;
1717                             }
1718                         }
1719                         if (gmx_debug_at)
1720                         {
1721                             fprintf(debug, "vsite atom %d  cg %d - %d pbc atom %d\n",
1722                                     vsite+1, cgs->index[cg_v]+1, cgs->index[cg_v+1],
1723                                     vsite_pbc_f[vsi]+1);
1724                         }
1725                     }
1726                 }
1727                 if (ftype == F_VSITEN)
1728                 {
1729                     /* The other entries in vsite_pbc_f are not used for center vsites */
1730                     i += nc3;
1731                 }
1732                 else
1733                 {
1734                     i += 1+nral;
1735                 }
1736
1737                 /* This vsite now has its pbc defined */
1738                 pbc_set[vsite] = 1;
1739             }
1740         }
1741     }
1742
1743     sfree(pbc_set);
1744
1745     return vsite_pbc;
1746 }
1747
1748
1749 gmx_vsite_t *init_vsite(gmx_mtop_t *mtop, t_commrec *cr,
1750                         gmx_bool bSerial_NoPBC)
1751 {
1752     int            nvsite, i;
1753     int           *a2cg, cg;
1754     gmx_vsite_t   *vsite;
1755     int            mt;
1756     gmx_moltype_t *molt;
1757     int            nthreads;
1758
1759     /* check if there are vsites */
1760     nvsite = 0;
1761     for (i = 0; i < F_NRE; i++)
1762     {
1763         if (interaction_function[i].flags & IF_VSITE)
1764         {
1765             nvsite += gmx_mtop_ftype_count(mtop, i);
1766         }
1767     }
1768
1769     if (nvsite == 0)
1770     {
1771         return NULL;
1772     }
1773
1774     snew(vsite, 1);
1775
1776     vsite->n_intercg_vsite = count_intercg_vsite(mtop,
1777                                                  &vsite->bHaveChargeGroups);
1778
1779     /* If we don't have charge groups, the vsite follows its own pbc */
1780     if (!bSerial_NoPBC &&
1781         vsite->bHaveChargeGroups &&
1782         vsite->n_intercg_vsite > 0 && DOMAINDECOMP(cr))
1783     {
1784         vsite->nvsite_pbc_molt = mtop->nmoltype;
1785         snew(vsite->vsite_pbc_molt, vsite->nvsite_pbc_molt);
1786         for (mt = 0; mt < mtop->nmoltype; mt++)
1787         {
1788             molt = &mtop->moltype[mt];
1789             /* Make an atom to charge group index */
1790             a2cg = atom2cg(&molt->cgs);
1791             vsite->vsite_pbc_molt[mt] = get_vsite_pbc(mtop->ffparams.iparams,
1792                                                       molt->ilist,
1793                                                       molt->atoms.atom, NULL,
1794                                                       &molt->cgs, a2cg);
1795             sfree(a2cg);
1796         }
1797
1798         snew(vsite->vsite_pbc_loc_nalloc, F_VSITEN-F_VSITE2+1);
1799         snew(vsite->vsite_pbc_loc, F_VSITEN-F_VSITE2+1);
1800     }
1801
1802     if (bSerial_NoPBC)
1803     {
1804         vsite->nthreads = 1;
1805     }
1806     else
1807     {
1808         vsite->nthreads = gmx_omp_nthreads_get(emntVSITE);
1809     }
1810     if (!bSerial_NoPBC)
1811     {
1812         /* We need one extra thread data structure for the overlap vsites */
1813         snew(vsite->tdata, vsite->nthreads+1);
1814     }
1815
1816     vsite->th_ind        = NULL;
1817     vsite->th_ind_nalloc = 0;
1818
1819     return vsite;
1820 }
1821
1822 static void prepare_vsite_thread(const t_ilist      *ilist,
1823                                  gmx_vsite_thread_t *vsite_th)
1824 {
1825     int ftype;
1826
1827     for (ftype = 0; ftype < F_NRE; ftype++)
1828     {
1829         if (interaction_function[ftype].flags & IF_VSITE)
1830         {
1831             if (ilist[ftype].nr > vsite_th->ilist[ftype].nalloc)
1832             {
1833                 vsite_th->ilist[ftype].nalloc = over_alloc_large(ilist[ftype].nr);
1834                 srenew(vsite_th->ilist[ftype].iatoms, vsite_th->ilist[ftype].nalloc);
1835             }
1836
1837             vsite_th->ilist[ftype].nr = 0;
1838         }
1839     }
1840 }
1841
1842 void split_vsites_over_threads(const t_ilist   *ilist,
1843                                const t_mdatoms *mdatoms,
1844                                gmx_bool         bLimitRange,
1845                                gmx_vsite_t     *vsite)
1846 {
1847     int      th;
1848     int      vsite_atom_range, natperthread;
1849     int     *th_ind;
1850     int      ftype;
1851     t_iatom *iat;
1852     t_ilist *il_th;
1853     int      nral1, inc, i, j;
1854
1855     if (vsite->nthreads == 1)
1856     {
1857         /* Nothing to do */
1858         return;
1859     }
1860
1861 #pragma omp parallel for num_threads(vsite->nthreads) schedule(static)
1862     for (th = 0; th < vsite->nthreads; th++)
1863     {
1864         prepare_vsite_thread(ilist, &vsite->tdata[th]);
1865     }
1866     /* Master threads does the (potential) overlap vsites */
1867     prepare_vsite_thread(ilist, &vsite->tdata[vsite->nthreads]);
1868
1869     /* The current way of distributing the vsites over threads in primitive.
1870      * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
1871      * without taking into account how the vsites are distributed.
1872      * Without domain decomposition we bLimitRange=TRUE and we at least
1873      * tighten the upper bound of the range (useful for common systems
1874      * such as a vsite-protein in 3-site water).
1875      */
1876     if (bLimitRange)
1877     {
1878         vsite_atom_range = -1;
1879         for (ftype = 0; ftype < F_NRE; ftype++)
1880         {
1881             if ((interaction_function[ftype].flags & IF_VSITE) &&
1882                 ftype != F_VSITEN)
1883             {
1884                 nral1 = 1 + NRAL(ftype);
1885                 iat   = ilist[ftype].iatoms;
1886                 for (i = 0; i < ilist[ftype].nr; i += nral1)
1887                 {
1888                     for (j = i+1; j < i+nral1; j++)
1889                     {
1890                         vsite_atom_range = max(vsite_atom_range, iat[j]);
1891                     }
1892                 }
1893             }
1894         }
1895         vsite_atom_range++;
1896     }
1897     else
1898     {
1899         vsite_atom_range = mdatoms->homenr;
1900     }
1901     natperthread = (vsite_atom_range + vsite->nthreads - 1)/vsite->nthreads;
1902
1903     if (debug)
1904     {
1905         fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n", mdatoms->nr, vsite_atom_range, natperthread);
1906     }
1907
1908     /* To simplify the vsite assignment, we make an index which tells us
1909      * to which thread particles, both non-vsites and vsites, are assigned.
1910      */
1911     if (mdatoms->nr > vsite->th_ind_nalloc)
1912     {
1913         vsite->th_ind_nalloc = over_alloc_large(mdatoms->nr);
1914         srenew(vsite->th_ind, vsite->th_ind_nalloc);
1915     }
1916     th_ind = vsite->th_ind;
1917     th     = 0;
1918     for (i = 0; i < mdatoms->nr; i++)
1919     {
1920         if (mdatoms->ptype[i] == eptVSite)
1921         {
1922             /* vsites are not assigned to a thread yet */
1923             th_ind[i] = -1;
1924         }
1925         else
1926         {
1927             /* assign non-vsite particles to thread th */
1928             th_ind[i] = th;
1929         }
1930         if (i == (th + 1)*natperthread && th < vsite->nthreads)
1931         {
1932             th++;
1933         }
1934     }
1935
1936     for (ftype = 0; ftype < F_NRE; ftype++)
1937     {
1938         if ((interaction_function[ftype].flags & IF_VSITE) &&
1939             ftype != F_VSITEN)
1940         {
1941             nral1 = 1 + NRAL(ftype);
1942             inc   = nral1;
1943             iat   = ilist[ftype].iatoms;
1944             for (i = 0; i < ilist[ftype].nr; )
1945             {
1946                 th = iat[1+i]/natperthread;
1947                 /* We would like to assign this vsite the thread th,
1948                  * but it might depend on atoms outside the atom range of th
1949                  * or on another vsite not assigned to thread th.
1950                  */
1951                 if (ftype != F_VSITEN)
1952                 {
1953                     for (j = i+2; j < i+nral1; j++)
1954                     {
1955                         if (th_ind[iat[j]] != th)
1956                         {
1957                             /* Some constructing atoms are not assigned to
1958                              * thread th, move this vsite to a separate batch.
1959                              */
1960                             th = vsite->nthreads;
1961                         }
1962                     }
1963                 }
1964                 else
1965                 {
1966                     inc = iat[i];
1967                     for (j = i+2; j < i+inc; j += 3)
1968                     {
1969                         if (th_ind[iat[j]] != th)
1970                         {
1971                             th = vsite->nthreads;
1972                         }
1973                     }
1974                 }
1975                 /* Copy this vsite to the thread data struct of thread th */
1976                 il_th = &vsite->tdata[th].ilist[ftype];
1977                 for (j = i; j < i+inc; j++)
1978                 {
1979                     il_th->iatoms[il_th->nr++] = iat[j];
1980                 }
1981                 /* Update this vsite's thread index entry */
1982                 th_ind[iat[1+i]] = th;
1983
1984                 i += inc;
1985             }
1986         }
1987     }
1988
1989     if (debug)
1990     {
1991         for (ftype = 0; ftype < F_NRE; ftype++)
1992         {
1993             if ((interaction_function[ftype].flags & IF_VSITE) &&
1994                 ilist[ftype].nr > 0)
1995             {
1996                 fprintf(debug, "%-20s thread dist:",
1997                         interaction_function[ftype].longname);
1998                 for (th = 0; th < vsite->nthreads+1; th++)
1999                 {
2000                     fprintf(debug, " %4d", vsite->tdata[th].ilist[ftype].nr);
2001                 }
2002                 fprintf(debug, "\n");
2003             }
2004         }
2005     }
2006 }
2007
2008 void set_vsite_top(gmx_vsite_t *vsite, gmx_localtop_t *top, t_mdatoms *md,
2009                    t_commrec *cr)
2010 {
2011     int *a2cg;
2012
2013     if (vsite->n_intercg_vsite > 0)
2014     {
2015         if (vsite->bHaveChargeGroups)
2016         {
2017             /* Make an atom to charge group index */
2018             a2cg                 = atom2cg(&top->cgs);
2019             vsite->vsite_pbc_loc = get_vsite_pbc(top->idef.iparams,
2020                                                  top->idef.il, NULL, md,
2021                                                  &top->cgs, a2cg);
2022             sfree(a2cg);
2023         }
2024
2025     }
2026
2027     if (vsite->nthreads > 1)
2028     {
2029         if (vsite->bHaveChargeGroups)
2030         {
2031             gmx_fatal(FARGS, "The combination of threading, virtual sites and charge groups is not implemented");
2032         }
2033
2034         split_vsites_over_threads(top->idef.il, md, !DOMAINDECOMP(cr), vsite);
2035     }
2036 }