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