Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / clincs.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 /* This file is completely threadsafe - keep it that way! */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <math.h>
43 #include <stdlib.h>
44
45 #include "types/commrec.h"
46 #include "constr.h"
47 #include "copyrite.h"
48 #include "gromacs/math/units.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/pbcutil/pbc.h"
51 #include "mdrun.h"
52 #include "nrnb.h"
53 #include "domdec.h"
54 #include "gromacs/topology/mtop_util.h"
55 #include "gmx_omp_nthreads.h"
56
57 #include "gromacs/fileio/gmxfio.h"
58 #include "gromacs/topology/block.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/gmxomp.h"
61 #include "gromacs/utility/smalloc.h"
62
63 typedef struct {
64     int    b0;         /* first constraint for this thread */
65     int    b1;         /* b1-1 is the last constraint for this thread */
66     int    nind;       /* number of indices */
67     int   *ind;        /* constraint index for updating atom data */
68     int    nind_r;     /* number of indices */
69     int   *ind_r;      /* constraint index for updating atom data */
70     int    ind_nalloc; /* allocation size of ind and ind_r */
71     tensor vir_r_m_dr; /* temporary variable for virial calculation */
72 } lincs_thread_t;
73
74 typedef struct gmx_lincsdata {
75     int             ncg;          /* the global number of constraints */
76     int             ncg_flex;     /* the global number of flexible constraints */
77     int             ncg_triangle; /* the global number of constraints in triangles */
78     int             nIter;        /* the number of iterations */
79     int             nOrder;       /* the order of the matrix expansion */
80     int             nc;           /* the number of constraints */
81     int             nc_alloc;     /* the number we allocated memory for */
82     int             ncc;          /* the number of constraint connections */
83     int             ncc_alloc;    /* the number we allocated memory for */
84     real            matlam;       /* the FE lambda value used for filling blc and blmf */
85     real           *bllen0;       /* the reference distance in topology A */
86     real           *ddist;        /* the reference distance in top B - the r.d. in top A */
87     int            *bla;          /* the atom pairs involved in the constraints */
88     real           *blc;          /* 1/sqrt(invmass1 + invmass2) */
89     real           *blc1;         /* as blc, but with all masses 1 */
90     int            *blnr;         /* index into blbnb and blmf */
91     int            *blbnb;        /* list of constraint connections */
92     int             ntriangle;    /* the local number of constraints in triangles */
93     int            *triangle;     /* the list of triangle constraints */
94     int            *tri_bits;     /* the bits tell if the matrix element should be used */
95     int             ncc_triangle; /* the number of constraint connections in triangles */
96     gmx_bool        bCommIter;    /* communicate before each LINCS interation */
97     real           *blmf;         /* matrix of mass factors for constraint connections */
98     real           *blmf1;        /* as blmf, but with all masses 1 */
99     real           *bllen;        /* the reference bond length */
100     int             nth;          /* The number of threads doing LINCS */
101     lincs_thread_t *th;           /* LINCS thread division */
102     unsigned       *atf;          /* atom flags for thread parallelization */
103     int             atf_nalloc;   /* allocation size of atf */
104     /* arrays for temporary storage in the LINCS algorithm */
105     rvec           *tmpv;
106     real           *tmpncc;
107     real           *tmp1;
108     real           *tmp2;
109     real           *tmp3;
110     real           *tmp4;
111     real           *mlambda; /* the Lagrange multipliers * -1 */
112     /* storage for the constraint RMS relative deviation output */
113     real            rmsd_data[3];
114 } t_gmx_lincsdata;
115
116 real *lincs_rmsd_data(struct gmx_lincsdata *lincsd)
117 {
118     return lincsd->rmsd_data;
119 }
120
121 real lincs_rmsd(struct gmx_lincsdata *lincsd, gmx_bool bSD2)
122 {
123     if (lincsd->rmsd_data[0] > 0)
124     {
125         return sqrt(lincsd->rmsd_data[bSD2 ? 2 : 1]/lincsd->rmsd_data[0]);
126     }
127     else
128     {
129         return 0;
130     }
131 }
132
133 /* Do a set of nrec LINCS matrix multiplications.
134  * This function will return with up to date thread-local
135  * constraint data, without an OpenMP barrier.
136  */
137 static void lincs_matrix_expand(const struct gmx_lincsdata *lincsd,
138                                 int b0, int b1,
139                                 const real *blcc,
140                                 real *rhs1, real *rhs2, real *sol)
141 {
142     int        nrec, rec, b, j, n, nr0, nr1;
143     real       mvb, *swap;
144     int        ntriangle, tb, bits;
145     const int *blnr     = lincsd->blnr, *blbnb = lincsd->blbnb;
146     const int *triangle = lincsd->triangle, *tri_bits = lincsd->tri_bits;
147
148     ntriangle = lincsd->ntriangle;
149     nrec      = lincsd->nOrder;
150
151     for (rec = 0; rec < nrec; rec++)
152     {
153 #pragma omp barrier
154         for (b = b0; b < b1; b++)
155         {
156             mvb = 0;
157             for (n = blnr[b]; n < blnr[b+1]; n++)
158             {
159                 j   = blbnb[n];
160                 mvb = mvb + blcc[n]*rhs1[j];
161             }
162             rhs2[b] = mvb;
163             sol[b]  = sol[b] + mvb;
164         }
165         swap = rhs1;
166         rhs1 = rhs2;
167         rhs2 = swap;
168     } /* nrec*(ncons+2*nrtot) flops */
169
170     if (ntriangle > 0)
171     {
172         /* Perform an extra nrec recursions for only the constraints
173          * involved in rigid triangles.
174          * In this way their accuracy should come close to those of the other
175          * constraints, since traingles of constraints can produce eigenvalues
176          * around 0.7, while the effective eigenvalue for bond constraints
177          * is around 0.4 (and 0.7*0.7=0.5).
178          */
179         /* We need to copy the temporary array, since only the elements
180          * for constraints involved in triangles are updated and then
181          * the pointers are swapped. This saving copying the whole arrary.
182          * We need barrier as other threads might still be reading from rhs2.
183          */
184 #pragma omp barrier
185         for (b = b0; b < b1; b++)
186         {
187             rhs2[b] = rhs1[b];
188         }
189 #pragma omp barrier
190 #pragma omp master
191         {
192             for (rec = 0; rec < nrec; rec++)
193             {
194                 for (tb = 0; tb < ntriangle; tb++)
195                 {
196                     b    = triangle[tb];
197                     bits = tri_bits[tb];
198                     mvb  = 0;
199                     nr0  = blnr[b];
200                     nr1  = blnr[b+1];
201                     for (n = nr0; n < nr1; n++)
202                     {
203                         if (bits & (1<<(n-nr0)))
204                         {
205                             j   = blbnb[n];
206                             mvb = mvb + blcc[n]*rhs1[j];
207                         }
208                     }
209                     rhs2[b] = mvb;
210                     sol[b]  = sol[b] + mvb;
211                 }
212                 swap = rhs1;
213                 rhs1 = rhs2;
214                 rhs2 = swap;
215             }
216         } /* flops count is missing here */
217
218         /* We need a barrier here as the calling routine will continue
219          * to operate on the thread-local constraints without barrier.
220          */
221 #pragma omp barrier
222     }
223 }
224
225 static void lincs_update_atoms_noind(int ncons, const int *bla,
226                                      real prefac,
227                                      const real *fac, rvec *r,
228                                      const real *invmass,
229                                      rvec *x)
230 {
231     int  b, i, j;
232     real mvb, im1, im2, tmp0, tmp1, tmp2;
233
234     if (invmass != NULL)
235     {
236         for (b = 0; b < ncons; b++)
237         {
238             i        = bla[2*b];
239             j        = bla[2*b+1];
240             mvb      = prefac*fac[b];
241             im1      = invmass[i];
242             im2      = invmass[j];
243             tmp0     = r[b][0]*mvb;
244             tmp1     = r[b][1]*mvb;
245             tmp2     = r[b][2]*mvb;
246             x[i][0] -= tmp0*im1;
247             x[i][1] -= tmp1*im1;
248             x[i][2] -= tmp2*im1;
249             x[j][0] += tmp0*im2;
250             x[j][1] += tmp1*im2;
251             x[j][2] += tmp2*im2;
252         } /* 16 ncons flops */
253     }
254     else
255     {
256         for (b = 0; b < ncons; b++)
257         {
258             i        = bla[2*b];
259             j        = bla[2*b+1];
260             mvb      = prefac*fac[b];
261             tmp0     = r[b][0]*mvb;
262             tmp1     = r[b][1]*mvb;
263             tmp2     = r[b][2]*mvb;
264             x[i][0] -= tmp0;
265             x[i][1] -= tmp1;
266             x[i][2] -= tmp2;
267             x[j][0] += tmp0;
268             x[j][1] += tmp1;
269             x[j][2] += tmp2;
270         }
271     }
272 }
273
274 static void lincs_update_atoms_ind(int ncons, const int *ind, const int *bla,
275                                    real prefac,
276                                    const real *fac, rvec *r,
277                                    const real *invmass,
278                                    rvec *x)
279 {
280     int  bi, b, i, j;
281     real mvb, im1, im2, tmp0, tmp1, tmp2;
282
283     if (invmass != NULL)
284     {
285         for (bi = 0; bi < ncons; bi++)
286         {
287             b        = ind[bi];
288             i        = bla[2*b];
289             j        = bla[2*b+1];
290             mvb      = prefac*fac[b];
291             im1      = invmass[i];
292             im2      = invmass[j];
293             tmp0     = r[b][0]*mvb;
294             tmp1     = r[b][1]*mvb;
295             tmp2     = r[b][2]*mvb;
296             x[i][0] -= tmp0*im1;
297             x[i][1] -= tmp1*im1;
298             x[i][2] -= tmp2*im1;
299             x[j][0] += tmp0*im2;
300             x[j][1] += tmp1*im2;
301             x[j][2] += tmp2*im2;
302         } /* 16 ncons flops */
303     }
304     else
305     {
306         for (bi = 0; bi < ncons; bi++)
307         {
308             b        = ind[bi];
309             i        = bla[2*b];
310             j        = bla[2*b+1];
311             mvb      = prefac*fac[b];
312             tmp0     = r[b][0]*mvb;
313             tmp1     = r[b][1]*mvb;
314             tmp2     = r[b][2]*mvb;
315             x[i][0] -= tmp0;
316             x[i][1] -= tmp1;
317             x[i][2] -= tmp2;
318             x[j][0] += tmp0;
319             x[j][1] += tmp1;
320             x[j][2] += tmp2;
321         } /* 16 ncons flops */
322     }
323 }
324
325 static void lincs_update_atoms(struct gmx_lincsdata *li, int th,
326                                real prefac,
327                                const real *fac, rvec *r,
328                                const real *invmass,
329                                rvec *x)
330 {
331     if (li->nth == 1)
332     {
333         /* Single thread, we simply update for all constraints */
334         lincs_update_atoms_noind(li->nc, li->bla, prefac, fac, r, invmass, x);
335     }
336     else
337     {
338         /* Update the atom vector components for our thread local
339          * constraints that only access our local atom range.
340          * This can be done without a barrier.
341          */
342         lincs_update_atoms_ind(li->th[th].nind, li->th[th].ind,
343                                li->bla, prefac, fac, r, invmass, x);
344
345         if (li->th[li->nth].nind > 0)
346         {
347             /* Update the constraints that operate on atoms
348              * in multiple thread atom blocks on the master thread.
349              */
350 #pragma omp barrier
351 #pragma omp master
352             {
353                 lincs_update_atoms_ind(li->th[li->nth].nind,
354                                        li->th[li->nth].ind,
355                                        li->bla, prefac, fac, r, invmass, x);
356             }
357         }
358     }
359 }
360
361 /* LINCS projection, works on derivatives of the coordinates */
362 static void do_lincsp(rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
363                       struct gmx_lincsdata *lincsd, int th,
364                       real *invmass,
365                       int econq, real *dvdlambda,
366                       gmx_bool bCalcVir, tensor rmdf)
367 {
368     int      b0, b1, b, i, j, k, n;
369     real     tmp0, tmp1, tmp2, im1, im2, mvb, rlen, len, wfac, lam;
370     rvec     dx;
371     int     *bla, *blnr, *blbnb;
372     rvec    *r;
373     real    *blc, *blmf, *blcc, *rhs1, *rhs2, *sol;
374
375     b0 = lincsd->th[th].b0;
376     b1 = lincsd->th[th].b1;
377
378     bla    = lincsd->bla;
379     r      = lincsd->tmpv;
380     blnr   = lincsd->blnr;
381     blbnb  = lincsd->blbnb;
382     if (econq != econqForce)
383     {
384         /* Use mass-weighted parameters */
385         blc  = lincsd->blc;
386         blmf = lincsd->blmf;
387     }
388     else
389     {
390         /* Use non mass-weighted parameters */
391         blc  = lincsd->blc1;
392         blmf = lincsd->blmf1;
393     }
394     blcc   = lincsd->tmpncc;
395     rhs1   = lincsd->tmp1;
396     rhs2   = lincsd->tmp2;
397     sol    = lincsd->tmp3;
398
399     /* Compute normalized i-j vectors */
400     if (pbc)
401     {
402         for (b = b0; b < b1; b++)
403         {
404             pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
405             unitv(dx, r[b]);
406         }
407     }
408     else
409     {
410         for (b = b0; b < b1; b++)
411         {
412             rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
413             unitv(dx, r[b]);
414         } /* 16 ncons flops */
415     }
416
417 #pragma omp barrier
418     for (b = b0; b < b1; b++)
419     {
420         tmp0 = r[b][0];
421         tmp1 = r[b][1];
422         tmp2 = r[b][2];
423         i    = bla[2*b];
424         j    = bla[2*b+1];
425         for (n = blnr[b]; n < blnr[b+1]; n++)
426         {
427             k       = blbnb[n];
428             blcc[n] = blmf[n]*(tmp0*r[k][0] + tmp1*r[k][1] + tmp2*r[k][2]);
429         } /* 6 nr flops */
430         mvb = blc[b]*(tmp0*(f[i][0] - f[j][0]) +
431                       tmp1*(f[i][1] - f[j][1]) +
432                       tmp2*(f[i][2] - f[j][2]));
433         rhs1[b] = mvb;
434         sol[b]  = mvb;
435         /* 7 flops */
436     }
437     /* Together: 23*ncons + 6*nrtot flops */
438
439     lincs_matrix_expand(lincsd, b0, b1, blcc, rhs1, rhs2, sol);
440     /* nrec*(ncons+2*nrtot) flops */
441
442     if (econq == econqDeriv_FlexCon)
443     {
444         /* We only want to constraint the flexible constraints,
445          * so we mask out the normal ones by setting sol to 0.
446          */
447         for (b = b0; b < b1; b++)
448         {
449             if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
450             {
451                 sol[b] = 0;
452             }
453         }
454     }
455
456     /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
457     for (b = b0; b < b1; b++)
458     {
459         sol[b] *= blc[b];
460     }
461
462     /* When constraining forces, we should not use mass weighting,
463      * so we pass invmass=NULL, which results in the use of 1 for all atoms.
464      */
465     lincs_update_atoms(lincsd, th, 1.0, sol, r,
466                        (econq != econqForce) ? invmass : NULL, fp);
467
468     if (dvdlambda != NULL)
469     {
470 #pragma omp barrier
471         for (b = b0; b < b1; b++)
472         {
473             *dvdlambda -= sol[b]*lincsd->ddist[b];
474         }
475         /* 10 ncons flops */
476     }
477
478     if (bCalcVir)
479     {
480         /* Constraint virial,
481          * determines sum r_bond x delta f,
482          * where delta f is the constraint correction
483          * of the quantity that is being constrained.
484          */
485         for (b = b0; b < b1; b++)
486         {
487             mvb = lincsd->bllen[b]*sol[b];
488             for (i = 0; i < DIM; i++)
489             {
490                 tmp1 = mvb*r[b][i];
491                 for (j = 0; j < DIM; j++)
492                 {
493                     rmdf[i][j] += tmp1*r[b][j];
494                 }
495             }
496         } /* 23 ncons flops */
497     }
498 }
499
500 static void do_lincs(rvec *x, rvec *xp, matrix box, t_pbc *pbc,
501                      struct gmx_lincsdata *lincsd, int th,
502                      real *invmass,
503                      t_commrec *cr,
504                      gmx_bool bCalcLambda,
505                      real wangle, int *warn,
506                      real invdt, rvec *v,
507                      gmx_bool bCalcVir, tensor vir_r_m_dr)
508 {
509     int      b0, b1, b, i, j, k, n, iter;
510     real     tmp0, tmp1, tmp2, im1, im2, mvb, rlen, len, len2, dlen2, wfac;
511     rvec     dx;
512     int     *bla, *blnr, *blbnb;
513     rvec    *r;
514     real    *blc, *blmf, *bllen, *blcc, *rhs1, *rhs2, *sol, *blc_sol, *mlambda;
515     int     *nlocat;
516
517     b0 = lincsd->th[th].b0;
518     b1 = lincsd->th[th].b1;
519
520     bla     = lincsd->bla;
521     r       = lincsd->tmpv;
522     blnr    = lincsd->blnr;
523     blbnb   = lincsd->blbnb;
524     blc     = lincsd->blc;
525     blmf    = lincsd->blmf;
526     bllen   = lincsd->bllen;
527     blcc    = lincsd->tmpncc;
528     rhs1    = lincsd->tmp1;
529     rhs2    = lincsd->tmp2;
530     sol     = lincsd->tmp3;
531     blc_sol = lincsd->tmp4;
532     mlambda = lincsd->mlambda;
533
534     if (DOMAINDECOMP(cr) && cr->dd->constraints)
535     {
536         nlocat = dd_constraints_nlocalatoms(cr->dd);
537     }
538     else
539     {
540         nlocat = NULL;
541     }
542
543     if (pbc)
544     {
545         /* Compute normalized i-j vectors */
546         for (b = b0; b < b1; b++)
547         {
548             pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
549             unitv(dx, r[b]);
550         }
551 #pragma omp barrier
552         for (b = b0; b < b1; b++)
553         {
554             for (n = blnr[b]; n < blnr[b+1]; n++)
555             {
556                 blcc[n] = blmf[n]*iprod(r[b], r[blbnb[n]]);
557             }
558             pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
559             mvb     = blc[b]*(iprod(r[b], dx) - bllen[b]);
560             rhs1[b] = mvb;
561             sol[b]  = mvb;
562         }
563     }
564     else
565     {
566         /* Compute normalized i-j vectors */
567         for (b = b0; b < b1; b++)
568         {
569             i       = bla[2*b];
570             j       = bla[2*b+1];
571             tmp0    = x[i][0] - x[j][0];
572             tmp1    = x[i][1] - x[j][1];
573             tmp2    = x[i][2] - x[j][2];
574             rlen    = gmx_invsqrt(tmp0*tmp0+tmp1*tmp1+tmp2*tmp2);
575             r[b][0] = rlen*tmp0;
576             r[b][1] = rlen*tmp1;
577             r[b][2] = rlen*tmp2;
578         } /* 16 ncons flops */
579
580 #pragma omp barrier
581         for (b = b0; b < b1; b++)
582         {
583             tmp0 = r[b][0];
584             tmp1 = r[b][1];
585             tmp2 = r[b][2];
586             len  = bllen[b];
587             i    = bla[2*b];
588             j    = bla[2*b+1];
589             for (n = blnr[b]; n < blnr[b+1]; n++)
590             {
591                 k       = blbnb[n];
592                 blcc[n] = blmf[n]*(tmp0*r[k][0] + tmp1*r[k][1] + tmp2*r[k][2]);
593             } /* 6 nr flops */
594             mvb = blc[b]*(tmp0*(xp[i][0] - xp[j][0]) +
595                           tmp1*(xp[i][1] - xp[j][1]) +
596                           tmp2*(xp[i][2] - xp[j][2]) - len);
597             rhs1[b] = mvb;
598             sol[b]  = mvb;
599             /* 10 flops */
600         }
601         /* Together: 26*ncons + 6*nrtot flops */
602     }
603
604     lincs_matrix_expand(lincsd, b0, b1, blcc, rhs1, rhs2, sol);
605     /* nrec*(ncons+2*nrtot) flops */
606
607     for (b = b0; b < b1; b++)
608     {
609         mlambda[b] = blc[b]*sol[b];
610     }
611
612     /* Update the coordinates */
613     lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
614
615     /*
616      ********  Correction for centripetal effects  ********
617      */
618
619     wfac = cos(DEG2RAD*wangle);
620     wfac = wfac*wfac;
621
622     for (iter = 0; iter < lincsd->nIter; iter++)
623     {
624         if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
625         {
626 #pragma omp barrier
627 #pragma omp master
628             {
629                 /* Communicate the corrected non-local coordinates */
630                 if (DOMAINDECOMP(cr))
631                 {
632                     dd_move_x_constraints(cr->dd, box, xp, NULL, FALSE);
633                 }
634             }
635         }
636
637 #pragma omp barrier
638         for (b = b0; b < b1; b++)
639         {
640             len = bllen[b];
641             if (pbc)
642             {
643                 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
644             }
645             else
646             {
647                 rvec_sub(xp[bla[2*b]], xp[bla[2*b+1]], dx);
648             }
649             len2  = len*len;
650             dlen2 = 2*len2 - norm2(dx);
651             if (dlen2 < wfac*len2 && (nlocat == NULL || nlocat[b]))
652             {
653                 *warn = b;
654             }
655             if (dlen2 > 0)
656             {
657                 mvb = blc[b]*(len - dlen2*gmx_invsqrt(dlen2));
658             }
659             else
660             {
661                 mvb = blc[b]*len;
662             }
663             rhs1[b] = mvb;
664             sol[b]  = mvb;
665         } /* 20*ncons flops */
666
667         lincs_matrix_expand(lincsd, b0, b1, blcc, rhs1, rhs2, sol);
668         /* nrec*(ncons+2*nrtot) flops */
669
670         for (b = b0; b < b1; b++)
671         {
672             mvb         = blc[b]*sol[b];
673             blc_sol[b]  = mvb;
674             mlambda[b] += mvb;
675         }
676
677         /* Update the coordinates */
678         lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
679     }
680     /* nit*ncons*(37+9*nrec) flops */
681
682     if (v != NULL)
683     {
684         /* Update the velocities */
685         lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
686         /* 16 ncons flops */
687     }
688
689     if (nlocat != NULL && bCalcLambda)
690     {
691         /* In lincs_update_atoms thread might cross-read mlambda */
692 #pragma omp barrier
693
694         /* Only account for local atoms */
695         for (b = b0; b < b1; b++)
696         {
697             mlambda[b] *= 0.5*nlocat[b];
698         }
699     }
700
701     if (bCalcVir)
702     {
703         /* Constraint virial */
704         for (b = b0; b < b1; b++)
705         {
706             tmp0 = -bllen[b]*mlambda[b];
707             for (i = 0; i < DIM; i++)
708             {
709                 tmp1 = tmp0*r[b][i];
710                 for (j = 0; j < DIM; j++)
711                 {
712                     vir_r_m_dr[i][j] -= tmp1*r[b][j];
713                 }
714             }
715         } /* 22 ncons flops */
716     }
717
718     /* Total:
719      * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
720      * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
721      *
722      * (26+nrec)*ncons + (6+2*nrec)*nrtot
723      * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
724      * if nit=1
725      * (63+nrec)*ncons + (6+4*nrec)*nrtot
726      */
727 }
728
729 void set_lincs_matrix(struct gmx_lincsdata *li, real *invmass, real lambda)
730 {
731     int        i, a1, a2, n, k, sign, center;
732     int        end, nk, kk;
733     const real invsqrt2 = 0.7071067811865475244;
734
735     for (i = 0; (i < li->nc); i++)
736     {
737         a1          = li->bla[2*i];
738         a2          = li->bla[2*i+1];
739         li->blc[i]  = gmx_invsqrt(invmass[a1] + invmass[a2]);
740         li->blc1[i] = invsqrt2;
741     }
742
743     /* Construct the coupling coefficient matrix blmf */
744     li->ntriangle    = 0;
745     li->ncc_triangle = 0;
746     for (i = 0; (i < li->nc); i++)
747     {
748         a1 = li->bla[2*i];
749         a2 = li->bla[2*i+1];
750         for (n = li->blnr[i]; (n < li->blnr[i+1]); n++)
751         {
752             k = li->blbnb[n];
753             if (a1 == li->bla[2*k] || a2 == li->bla[2*k+1])
754             {
755                 sign = -1;
756             }
757             else
758             {
759                 sign = 1;
760             }
761             if (a1 == li->bla[2*k] || a1 == li->bla[2*k+1])
762             {
763                 center = a1;
764                 end    = a2;
765             }
766             else
767             {
768                 center = a2;
769                 end    = a1;
770             }
771             li->blmf[n]  = sign*invmass[center]*li->blc[i]*li->blc[k];
772             li->blmf1[n] = sign*0.5;
773             if (li->ncg_triangle > 0)
774             {
775                 /* Look for constraint triangles */
776                 for (nk = li->blnr[k]; (nk < li->blnr[k+1]); nk++)
777                 {
778                     kk = li->blbnb[nk];
779                     if (kk != i && kk != k &&
780                         (li->bla[2*kk] == end || li->bla[2*kk+1] == end))
781                     {
782                         if (li->ntriangle == 0 ||
783                             li->triangle[li->ntriangle-1] < i)
784                         {
785                             /* Add this constraint to the triangle list */
786                             li->triangle[li->ntriangle] = i;
787                             li->tri_bits[li->ntriangle] = 0;
788                             li->ntriangle++;
789                             if (li->blnr[i+1] - li->blnr[i] > sizeof(li->tri_bits[0])*8 - 1)
790                             {
791                                 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %d allowed for constraints participating in triangles",
792                                           li->blnr[i+1] - li->blnr[i],
793                                           sizeof(li->tri_bits[0])*8-1);
794                             }
795                         }
796                         li->tri_bits[li->ntriangle-1] |= (1<<(n-li->blnr[i]));
797                         li->ncc_triangle++;
798                     }
799                 }
800             }
801         }
802     }
803
804     if (debug)
805     {
806         fprintf(debug, "Of the %d constraints %d participate in triangles\n",
807                 li->nc, li->ntriangle);
808         fprintf(debug, "There are %d couplings of which %d in triangles\n",
809                 li->ncc, li->ncc_triangle);
810     }
811
812     /* Set matlam,
813      * so we know with which lambda value the masses have been set.
814      */
815     li->matlam = lambda;
816 }
817
818 static int count_triangle_constraints(t_ilist *ilist, t_blocka *at2con)
819 {
820     int      ncon1, ncon_tot;
821     int      c0, a00, a01, n1, c1, a10, a11, ac1, n2, c2, a20, a21;
822     int      ncon_triangle;
823     gmx_bool bTriangle;
824     t_iatom *ia1, *ia2, *iap;
825
826     ncon1    = ilist[F_CONSTR].nr/3;
827     ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
828
829     ia1 = ilist[F_CONSTR].iatoms;
830     ia2 = ilist[F_CONSTRNC].iatoms;
831
832     ncon_triangle = 0;
833     for (c0 = 0; c0 < ncon_tot; c0++)
834     {
835         bTriangle = FALSE;
836         iap       = constr_iatomptr(ncon1, ia1, ia2, c0);
837         a00       = iap[1];
838         a01       = iap[2];
839         for (n1 = at2con->index[a01]; n1 < at2con->index[a01+1]; n1++)
840         {
841             c1 = at2con->a[n1];
842             if (c1 != c0)
843             {
844                 iap = constr_iatomptr(ncon1, ia1, ia2, c1);
845                 a10 = iap[1];
846                 a11 = iap[2];
847                 if (a10 == a01)
848                 {
849                     ac1 = a11;
850                 }
851                 else
852                 {
853                     ac1 = a10;
854                 }
855                 for (n2 = at2con->index[ac1]; n2 < at2con->index[ac1+1]; n2++)
856                 {
857                     c2 = at2con->a[n2];
858                     if (c2 != c0 && c2 != c1)
859                     {
860                         iap = constr_iatomptr(ncon1, ia1, ia2, c2);
861                         a20 = iap[1];
862                         a21 = iap[2];
863                         if (a20 == a00 || a21 == a00)
864                         {
865                             bTriangle = TRUE;
866                         }
867                     }
868                 }
869             }
870         }
871         if (bTriangle)
872         {
873             ncon_triangle++;
874         }
875     }
876
877     return ncon_triangle;
878 }
879
880 static gmx_bool more_than_two_sequential_constraints(const t_ilist  *ilist,
881                                                      const t_blocka *at2con)
882 {
883     t_iatom  *ia1, *ia2, *iap;
884     int       ncon1, ncon_tot, c;
885     int       a1, a2;
886     gmx_bool  bMoreThanTwoSequentialConstraints;
887
888     ncon1    = ilist[F_CONSTR].nr/3;
889     ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
890
891     ia1 = ilist[F_CONSTR].iatoms;
892     ia2 = ilist[F_CONSTRNC].iatoms;
893
894     bMoreThanTwoSequentialConstraints = FALSE;
895     for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
896     {
897         iap = constr_iatomptr(ncon1, ia1, ia2, c);
898         a1  = iap[1];
899         a2  = iap[2];
900         /* Check if this constraint has constraints connected at both atoms */
901         if (at2con->index[a1+1] - at2con->index[a1] > 1 &&
902             at2con->index[a2+1] - at2con->index[a2] > 1)
903         {
904             bMoreThanTwoSequentialConstraints = TRUE;
905         }
906     }
907
908     return bMoreThanTwoSequentialConstraints;
909 }
910
911 static int int_comp(const void *a, const void *b)
912 {
913     return (*(int *)a) - (*(int *)b);
914 }
915
916 gmx_lincsdata_t init_lincs(FILE *fplog, gmx_mtop_t *mtop,
917                            int nflexcon_global, t_blocka *at2con,
918                            gmx_bool bPLINCS, int nIter, int nProjOrder)
919 {
920     struct gmx_lincsdata *li;
921     int                   mb;
922     gmx_moltype_t        *molt;
923
924     if (fplog)
925     {
926         fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
927                 bPLINCS ? " Parallel" : "");
928     }
929
930     snew(li, 1);
931
932     li->ncg      =
933         gmx_mtop_ftype_count(mtop, F_CONSTR) +
934         gmx_mtop_ftype_count(mtop, F_CONSTRNC);
935     li->ncg_flex = nflexcon_global;
936
937     li->nIter  = nIter;
938     li->nOrder = nProjOrder;
939
940     li->ncg_triangle = 0;
941     li->bCommIter    = FALSE;
942     for (mb = 0; mb < mtop->nmolblock; mb++)
943     {
944         molt              = &mtop->moltype[mtop->molblock[mb].type];
945         li->ncg_triangle +=
946             mtop->molblock[mb].nmol*
947             count_triangle_constraints(molt->ilist,
948                                        &at2con[mtop->molblock[mb].type]);
949         if (bPLINCS && li->bCommIter == FALSE)
950         {
951             /* Check if we need to communicate not only before LINCS,
952              * but also before each iteration.
953              * The check for only two sequential constraints is only
954              * useful for the common case of H-bond only constraints.
955              * With more effort we could also make it useful for small
956              * molecules with nr. sequential constraints <= nOrder-1.
957              */
958             li->bCommIter = (li->nOrder < 1 || more_than_two_sequential_constraints(molt->ilist, &at2con[mtop->molblock[mb].type]));
959         }
960     }
961     if (debug && bPLINCS)
962     {
963         fprintf(debug, "PLINCS communication before each iteration: %d\n",
964                 li->bCommIter);
965     }
966
967     /* LINCS can run on any number of threads.
968      * Currently the number is fixed for the whole simulation,
969      * but it could be set in set_lincs().
970      */
971     li->nth = gmx_omp_nthreads_get(emntLINCS);
972     if (li->nth == 1)
973     {
974         snew(li->th, 1);
975     }
976     else
977     {
978         /* Allocate an extra elements for "thread-overlap" constraints */
979         snew(li->th, li->nth+1);
980     }
981     if (debug)
982     {
983         fprintf(debug, "LINCS: using %d threads\n", li->nth);
984     }
985
986     if (bPLINCS || li->ncg_triangle > 0)
987     {
988         please_cite(fplog, "Hess2008a");
989     }
990     else
991     {
992         please_cite(fplog, "Hess97a");
993     }
994
995     if (fplog)
996     {
997         fprintf(fplog, "The number of constraints is %d\n", li->ncg);
998         if (bPLINCS)
999         {
1000             fprintf(fplog, "There are inter charge-group constraints,\n"
1001                     "will communicate selected coordinates each lincs iteration\n");
1002         }
1003         if (li->ncg_triangle > 0)
1004         {
1005             fprintf(fplog,
1006                     "%d constraints are involved in constraint triangles,\n"
1007                     "will apply an additional matrix expansion of order %d for couplings\n"
1008                     "between constraints inside triangles\n",
1009                     li->ncg_triangle, li->nOrder);
1010         }
1011     }
1012
1013     return li;
1014 }
1015
1016 /* Sets up the work division over the threads */
1017 static void lincs_thread_setup(struct gmx_lincsdata *li, int natoms)
1018 {
1019     lincs_thread_t *li_m;
1020     int             th;
1021     unsigned       *atf;
1022     int             a;
1023
1024     if (natoms > li->atf_nalloc)
1025     {
1026         li->atf_nalloc = over_alloc_large(natoms);
1027         srenew(li->atf, li->atf_nalloc);
1028     }
1029
1030     atf = li->atf;
1031     /* Clear the atom flags */
1032     for (a = 0; a < natoms; a++)
1033     {
1034         atf[a] = 0;
1035     }
1036
1037     for (th = 0; th < li->nth; th++)
1038     {
1039         lincs_thread_t *li_th;
1040         int             b;
1041
1042         li_th = &li->th[th];
1043
1044         /* The constraints are divided equally over the threads */
1045         li_th->b0 = (li->nc* th   )/li->nth;
1046         li_th->b1 = (li->nc*(th+1))/li->nth;
1047
1048         if (th < sizeof(*atf)*8)
1049         {
1050             /* For each atom set a flag for constraints from each */
1051             for (b = li_th->b0; b < li_th->b1; b++)
1052             {
1053                 atf[li->bla[b*2]  ] |= (1U<<th);
1054                 atf[li->bla[b*2+1]] |= (1U<<th);
1055             }
1056         }
1057     }
1058
1059 #pragma omp parallel for num_threads(li->nth) schedule(static)
1060     for (th = 0; th < li->nth; th++)
1061     {
1062         lincs_thread_t *li_th;
1063         unsigned        mask;
1064         int             b;
1065
1066         li_th = &li->th[th];
1067
1068         if (li_th->b1 - li_th->b0 > li_th->ind_nalloc)
1069         {
1070             li_th->ind_nalloc = over_alloc_large(li_th->b1-li_th->b0);
1071             srenew(li_th->ind, li_th->ind_nalloc);
1072             srenew(li_th->ind_r, li_th->ind_nalloc);
1073         }
1074
1075         if (th < sizeof(*atf)*8)
1076         {
1077             mask = (1U<<th) - 1U;
1078
1079             li_th->nind   = 0;
1080             li_th->nind_r = 0;
1081             for (b = li_th->b0; b < li_th->b1; b++)
1082             {
1083                 /* We let the constraint with the lowest thread index
1084                  * operate on atoms with constraints from multiple threads.
1085                  */
1086                 if (((atf[li->bla[b*2]]   & mask) == 0) &&
1087                     ((atf[li->bla[b*2+1]] & mask) == 0))
1088                 {
1089                     /* Add the constraint to the local atom update index */
1090                     li_th->ind[li_th->nind++] = b;
1091                 }
1092                 else
1093                 {
1094                     /* Add the constraint to the rest block */
1095                     li_th->ind_r[li_th->nind_r++] = b;
1096                 }
1097             }
1098         }
1099         else
1100         {
1101             /* We are out of bits, assign all constraints to rest */
1102             for (b = li_th->b0; b < li_th->b1; b++)
1103             {
1104                 li_th->ind_r[li_th->nind_r++] = b;
1105             }
1106         }
1107     }
1108
1109     /* We need to copy all constraints which have not be assigned
1110      * to a thread to a separate list which will be handled by one thread.
1111      */
1112     li_m = &li->th[li->nth];
1113
1114     li_m->nind = 0;
1115     for (th = 0; th < li->nth; th++)
1116     {
1117         lincs_thread_t *li_th;
1118         int             b;
1119
1120         li_th   = &li->th[th];
1121
1122         if (li_m->nind + li_th->nind_r > li_m->ind_nalloc)
1123         {
1124             li_m->ind_nalloc = over_alloc_large(li_m->nind+li_th->nind_r);
1125             srenew(li_m->ind, li_m->ind_nalloc);
1126         }
1127
1128         for (b = 0; b < li_th->nind_r; b++)
1129         {
1130             li_m->ind[li_m->nind++] = li_th->ind_r[b];
1131         }
1132
1133         if (debug)
1134         {
1135             fprintf(debug, "LINCS thread %d: %d constraints\n",
1136                     th, li_th->nind);
1137         }
1138     }
1139
1140     if (debug)
1141     {
1142         fprintf(debug, "LINCS thread r: %d constraints\n",
1143                 li_m->nind);
1144     }
1145 }
1146
1147
1148 void set_lincs(t_idef *idef, t_mdatoms *md,
1149                gmx_bool bDynamics, t_commrec *cr,
1150                struct gmx_lincsdata *li)
1151 {
1152     int          start, natoms, nflexcon;
1153     t_blocka     at2con;
1154     t_iatom     *iatom;
1155     int          i, k, ncc_alloc, ni, con, nconnect, concon;
1156     int          type, a1, a2;
1157     real         lenA = 0, lenB;
1158     gmx_bool     bLocal;
1159
1160     li->nc  = 0;
1161     li->ncc = 0;
1162     /* Zero the thread index ranges.
1163      * Otherwise without local constraints we could return with old ranges.
1164      */
1165     for (i = 0; i < li->nth; i++)
1166     {
1167         li->th[i].b0   = 0;
1168         li->th[i].b1   = 0;
1169         li->th[i].nind = 0;
1170     }
1171     if (li->nth > 1)
1172     {
1173         li->th[li->nth].nind = 0;
1174     }
1175
1176     /* This is the local topology, so there are only F_CONSTR constraints */
1177     if (idef->il[F_CONSTR].nr == 0)
1178     {
1179         /* There are no constraints,
1180          * we do not need to fill any data structures.
1181          */
1182         return;
1183     }
1184
1185     if (debug)
1186     {
1187         fprintf(debug, "Building the LINCS connectivity\n");
1188     }
1189
1190     if (DOMAINDECOMP(cr))
1191     {
1192         if (cr->dd->constraints)
1193         {
1194             dd_get_constraint_range(cr->dd, &start, &natoms);
1195         }
1196         else
1197         {
1198             natoms = cr->dd->nat_home;
1199         }
1200         start = 0;
1201     }
1202     else
1203     {
1204         start  = 0;
1205         natoms = md->homenr;
1206     }
1207     at2con = make_at2con(start, natoms, idef->il, idef->iparams, bDynamics,
1208                          &nflexcon);
1209
1210
1211     if (idef->il[F_CONSTR].nr/3 > li->nc_alloc || li->nc_alloc == 0)
1212     {
1213         li->nc_alloc = over_alloc_dd(idef->il[F_CONSTR].nr/3);
1214         srenew(li->bllen0, li->nc_alloc);
1215         srenew(li->ddist, li->nc_alloc);
1216         srenew(li->bla, 2*li->nc_alloc);
1217         srenew(li->blc, li->nc_alloc);
1218         srenew(li->blc1, li->nc_alloc);
1219         srenew(li->blnr, li->nc_alloc+1);
1220         srenew(li->bllen, li->nc_alloc);
1221         srenew(li->tmpv, li->nc_alloc);
1222         srenew(li->tmp1, li->nc_alloc);
1223         srenew(li->tmp2, li->nc_alloc);
1224         srenew(li->tmp3, li->nc_alloc);
1225         srenew(li->tmp4, li->nc_alloc);
1226         srenew(li->mlambda, li->nc_alloc);
1227         if (li->ncg_triangle > 0)
1228         {
1229             /* This is allocating too much, but it is difficult to improve */
1230             srenew(li->triangle, li->nc_alloc);
1231             srenew(li->tri_bits, li->nc_alloc);
1232         }
1233     }
1234
1235     iatom = idef->il[F_CONSTR].iatoms;
1236
1237     ncc_alloc   = li->ncc_alloc;
1238     li->blnr[0] = 0;
1239
1240     ni = idef->il[F_CONSTR].nr/3;
1241
1242     con           = 0;
1243     nconnect      = 0;
1244     li->blnr[con] = nconnect;
1245     for (i = 0; i < ni; i++)
1246     {
1247         bLocal = TRUE;
1248         type   = iatom[3*i];
1249         a1     = iatom[3*i+1];
1250         a2     = iatom[3*i+2];
1251         lenA   = idef->iparams[type].constr.dA;
1252         lenB   = idef->iparams[type].constr.dB;
1253         /* Skip the flexible constraints when not doing dynamics */
1254         if (bDynamics || lenA != 0 || lenB != 0)
1255         {
1256             li->bllen0[con]  = lenA;
1257             li->ddist[con]   = lenB - lenA;
1258             /* Set the length to the topology A length */
1259             li->bllen[con]   = li->bllen0[con];
1260             li->bla[2*con]   = a1;
1261             li->bla[2*con+1] = a2;
1262             /* Construct the constraint connection matrix blbnb */
1263             for (k = at2con.index[a1-start]; k < at2con.index[a1-start+1]; k++)
1264             {
1265                 concon = at2con.a[k];
1266                 if (concon != i)
1267                 {
1268                     if (nconnect >= ncc_alloc)
1269                     {
1270                         ncc_alloc = over_alloc_small(nconnect+1);
1271                         srenew(li->blbnb, ncc_alloc);
1272                     }
1273                     li->blbnb[nconnect++] = concon;
1274                 }
1275             }
1276             for (k = at2con.index[a2-start]; k < at2con.index[a2-start+1]; k++)
1277             {
1278                 concon = at2con.a[k];
1279                 if (concon != i)
1280                 {
1281                     if (nconnect+1 > ncc_alloc)
1282                     {
1283                         ncc_alloc = over_alloc_small(nconnect+1);
1284                         srenew(li->blbnb, ncc_alloc);
1285                     }
1286                     li->blbnb[nconnect++] = concon;
1287                 }
1288             }
1289             li->blnr[con+1] = nconnect;
1290
1291             if (cr->dd == NULL)
1292             {
1293                 /* Order the blbnb matrix to optimize memory access */
1294                 qsort(&(li->blbnb[li->blnr[con]]), li->blnr[con+1]-li->blnr[con],
1295                       sizeof(li->blbnb[0]), int_comp);
1296             }
1297             /* Increase the constraint count */
1298             con++;
1299         }
1300     }
1301
1302     done_blocka(&at2con);
1303
1304     /* This is the real number of constraints,
1305      * without dynamics the flexible constraints are not present.
1306      */
1307     li->nc = con;
1308
1309     li->ncc = li->blnr[con];
1310     if (cr->dd == NULL)
1311     {
1312         /* Since the matrix is static, we can free some memory */
1313         ncc_alloc = li->ncc;
1314         srenew(li->blbnb, ncc_alloc);
1315     }
1316
1317     if (ncc_alloc > li->ncc_alloc)
1318     {
1319         li->ncc_alloc = ncc_alloc;
1320         srenew(li->blmf, li->ncc_alloc);
1321         srenew(li->blmf1, li->ncc_alloc);
1322         srenew(li->tmpncc, li->ncc_alloc);
1323     }
1324
1325     if (debug)
1326     {
1327         fprintf(debug, "Number of constraints is %d, couplings %d\n",
1328                 li->nc, li->ncc);
1329     }
1330
1331     if (li->nth == 1)
1332     {
1333         li->th[0].b0 = 0;
1334         li->th[0].b1 = li->nc;
1335     }
1336     else
1337     {
1338         lincs_thread_setup(li, md->nr);
1339     }
1340
1341     set_lincs_matrix(li, md->invmass, md->lambda);
1342 }
1343
1344 static void lincs_warning(FILE *fplog,
1345                           gmx_domdec_t *dd, rvec *x, rvec *xprime, t_pbc *pbc,
1346                           int ncons, int *bla, real *bllen, real wangle,
1347                           int maxwarn, int *warncount)
1348 {
1349     int  b, i, j;
1350     rvec v0, v1;
1351     real wfac, d0, d1, cosine;
1352     char buf[STRLEN];
1353
1354     wfac = cos(DEG2RAD*wangle);
1355
1356     sprintf(buf, "bonds that rotated more than %g degrees:\n"
1357             " atom 1 atom 2  angle  previous, current, constraint length\n",
1358             wangle);
1359     fprintf(stderr, "%s", buf);
1360     if (fplog)
1361     {
1362         fprintf(fplog, "%s", buf);
1363     }
1364
1365     for (b = 0; b < ncons; b++)
1366     {
1367         i = bla[2*b];
1368         j = bla[2*b+1];
1369         if (pbc)
1370         {
1371             pbc_dx_aiuc(pbc, x[i], x[j], v0);
1372             pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
1373         }
1374         else
1375         {
1376             rvec_sub(x[i], x[j], v0);
1377             rvec_sub(xprime[i], xprime[j], v1);
1378         }
1379         d0     = norm(v0);
1380         d1     = norm(v1);
1381         cosine = iprod(v0, v1)/(d0*d1);
1382         if (cosine < wfac)
1383         {
1384             sprintf(buf, " %6d %6d  %5.1f  %8.4f %8.4f    %8.4f\n",
1385                     ddglatnr(dd, i), ddglatnr(dd, j),
1386                     RAD2DEG*acos(cosine), d0, d1, bllen[b]);
1387             fprintf(stderr, "%s", buf);
1388             if (fplog)
1389             {
1390                 fprintf(fplog, "%s", buf);
1391             }
1392             if (!gmx_isfinite(d1))
1393             {
1394                 gmx_fatal(FARGS, "Bond length not finite.");
1395             }
1396
1397             (*warncount)++;
1398         }
1399     }
1400     if (*warncount > maxwarn)
1401     {
1402         too_many_constraint_warnings(econtLINCS, *warncount);
1403     }
1404 }
1405
1406 static void cconerr(gmx_domdec_t *dd,
1407                     int ncons, int *bla, real *bllen, rvec *x, t_pbc *pbc,
1408                     real *ncons_loc, real *ssd, real *max, int *imax)
1409 {
1410     real       len, d, ma, ssd2, r2;
1411     int       *nlocat, count, b, im;
1412     rvec       dx;
1413
1414     if (dd && dd->constraints)
1415     {
1416         nlocat = dd_constraints_nlocalatoms(dd);
1417     }
1418     else
1419     {
1420         nlocat = 0;
1421     }
1422
1423     ma    = 0;
1424     ssd2  = 0;
1425     im    = 0;
1426     count = 0;
1427     for (b = 0; b < ncons; b++)
1428     {
1429         if (pbc)
1430         {
1431             pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
1432         }
1433         else
1434         {
1435             rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
1436         }
1437         r2  = norm2(dx);
1438         len = r2*gmx_invsqrt(r2);
1439         d   = fabs(len/bllen[b]-1);
1440         if (d > ma && (nlocat == NULL || nlocat[b]))
1441         {
1442             ma = d;
1443             im = b;
1444         }
1445         if (nlocat == NULL)
1446         {
1447             ssd2 += d*d;
1448             count++;
1449         }
1450         else
1451         {
1452             ssd2  += nlocat[b]*d*d;
1453             count += nlocat[b];
1454         }
1455     }
1456
1457     *ncons_loc = (nlocat ? 0.5 : 1)*count;
1458     *ssd       = (nlocat ? 0.5 : 1)*ssd2;
1459     *max       = ma;
1460     *imax      = im;
1461 }
1462
1463 gmx_bool constrain_lincs(FILE *fplog, gmx_bool bLog, gmx_bool bEner,
1464                          t_inputrec *ir,
1465                          gmx_int64_t step,
1466                          struct gmx_lincsdata *lincsd, t_mdatoms *md,
1467                          t_commrec *cr,
1468                          rvec *x, rvec *xprime, rvec *min_proj,
1469                          matrix box, t_pbc *pbc,
1470                          real lambda, real *dvdlambda,
1471                          real invdt, rvec *v,
1472                          gmx_bool bCalcVir, tensor vir_r_m_dr,
1473                          int econq,
1474                          t_nrnb *nrnb,
1475                          int maxwarn, int *warncount)
1476 {
1477     char      buf[STRLEN], buf2[22], buf3[STRLEN];
1478     int       i, warn, p_imax, error;
1479     real      ncons_loc, p_ssd, p_max = 0;
1480     rvec      dx;
1481     gmx_bool  bOK;
1482
1483     bOK = TRUE;
1484
1485     if (lincsd->nc == 0 && cr->dd == NULL)
1486     {
1487         if (bLog || bEner)
1488         {
1489             lincsd->rmsd_data[0] = 0;
1490             if (ir->eI == eiSD2 && v == NULL)
1491             {
1492                 i = 2;
1493             }
1494             else
1495             {
1496                 i = 1;
1497             }
1498             lincsd->rmsd_data[i] = 0;
1499         }
1500
1501         return bOK;
1502     }
1503
1504     if (econq == econqCoord)
1505     {
1506         if (ir->efep != efepNO)
1507         {
1508             if (md->nMassPerturbed && lincsd->matlam != md->lambda)
1509             {
1510                 set_lincs_matrix(lincsd, md->invmass, md->lambda);
1511             }
1512
1513             for (i = 0; i < lincsd->nc; i++)
1514             {
1515                 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
1516             }
1517         }
1518
1519         if (lincsd->ncg_flex)
1520         {
1521             /* Set the flexible constraint lengths to the old lengths */
1522             if (pbc != NULL)
1523             {
1524                 for (i = 0; i < lincsd->nc; i++)
1525                 {
1526                     if (lincsd->bllen[i] == 0)
1527                     {
1528                         pbc_dx_aiuc(pbc, x[lincsd->bla[2*i]], x[lincsd->bla[2*i+1]], dx);
1529                         lincsd->bllen[i] = norm(dx);
1530                     }
1531                 }
1532             }
1533             else
1534             {
1535                 for (i = 0; i < lincsd->nc; i++)
1536                 {
1537                     if (lincsd->bllen[i] == 0)
1538                     {
1539                         lincsd->bllen[i] =
1540                             sqrt(distance2(x[lincsd->bla[2*i]],
1541                                            x[lincsd->bla[2*i+1]]));
1542                     }
1543                 }
1544             }
1545         }
1546
1547         if (bLog && fplog)
1548         {
1549             cconerr(cr->dd, lincsd->nc, lincsd->bla, lincsd->bllen, xprime, pbc,
1550                     &ncons_loc, &p_ssd, &p_max, &p_imax);
1551         }
1552
1553         /* This warn var can be updated by multiple threads
1554          * at the same time. But as we only need to detect
1555          * if a warning occured or not, this is not an issue.
1556          */
1557         warn = -1;
1558
1559         /* The OpenMP parallel region of constrain_lincs for coords */
1560 #pragma omp parallel num_threads(lincsd->nth)
1561         {
1562             int th = gmx_omp_get_thread_num();
1563
1564             clear_mat(lincsd->th[th].vir_r_m_dr);
1565
1566             do_lincs(x, xprime, box, pbc, lincsd, th,
1567                      md->invmass, cr,
1568                      bCalcVir || (ir->efep != efepNO),
1569                      ir->LincsWarnAngle, &warn,
1570                      invdt, v, bCalcVir,
1571                      th == 0 ? vir_r_m_dr : lincsd->th[th].vir_r_m_dr);
1572         }
1573
1574         if (ir->efep != efepNO)
1575         {
1576             real dt_2, dvdl = 0;
1577
1578             dt_2 = 1.0/(ir->delta_t*ir->delta_t);
1579             for (i = 0; (i < lincsd->nc); i++)
1580             {
1581                 dvdl -= lincsd->mlambda[i]*dt_2*lincsd->ddist[i];
1582             }
1583             *dvdlambda += dvdl;
1584         }
1585
1586         if (bLog && fplog && lincsd->nc > 0)
1587         {
1588             fprintf(fplog, "   Rel. Constraint Deviation:  RMS         MAX     between atoms\n");
1589             fprintf(fplog, "       Before LINCS          %.6f    %.6f %6d %6d\n",
1590                     sqrt(p_ssd/ncons_loc), p_max,
1591                     ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
1592                     ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
1593         }
1594         if (bLog || bEner)
1595         {
1596             cconerr(cr->dd, lincsd->nc, lincsd->bla, lincsd->bllen, xprime, pbc,
1597                     &ncons_loc, &p_ssd, &p_max, &p_imax);
1598             /* Check if we are doing the second part of SD */
1599             if (ir->eI == eiSD2 && v == NULL)
1600             {
1601                 i = 2;
1602             }
1603             else
1604             {
1605                 i = 1;
1606             }
1607             lincsd->rmsd_data[0] = ncons_loc;
1608             lincsd->rmsd_data[i] = p_ssd;
1609         }
1610         else
1611         {
1612             lincsd->rmsd_data[0] = 0;
1613             lincsd->rmsd_data[1] = 0;
1614             lincsd->rmsd_data[2] = 0;
1615         }
1616         if (bLog && fplog && lincsd->nc > 0)
1617         {
1618             fprintf(fplog,
1619                     "        After LINCS          %.6f    %.6f %6d %6d\n\n",
1620                     sqrt(p_ssd/ncons_loc), p_max,
1621                     ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
1622                     ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
1623         }
1624
1625         if (warn >= 0)
1626         {
1627             if (maxwarn >= 0)
1628             {
1629                 cconerr(cr->dd, lincsd->nc, lincsd->bla, lincsd->bllen, xprime, pbc,
1630                         &ncons_loc, &p_ssd, &p_max, &p_imax);
1631                 if (MULTISIM(cr))
1632                 {
1633                     sprintf(buf3, " in simulation %d", cr->ms->sim);
1634                 }
1635                 else
1636                 {
1637                     buf3[0] = 0;
1638                 }
1639                 sprintf(buf, "\nStep %s, time %g (ps)  LINCS WARNING%s\n"
1640                         "relative constraint deviation after LINCS:\n"
1641                         "rms %.6f, max %.6f (between atoms %d and %d)\n",
1642                         gmx_step_str(step, buf2), ir->init_t+step*ir->delta_t,
1643                         buf3,
1644                         sqrt(p_ssd/ncons_loc), p_max,
1645                         ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
1646                         ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
1647                 if (fplog)
1648                 {
1649                     fprintf(fplog, "%s", buf);
1650                 }
1651                 fprintf(stderr, "%s", buf);
1652                 lincs_warning(fplog, cr->dd, x, xprime, pbc,
1653                               lincsd->nc, lincsd->bla, lincsd->bllen,
1654                               ir->LincsWarnAngle, maxwarn, warncount);
1655             }
1656             bOK = (p_max < 0.5);
1657         }
1658
1659         if (lincsd->ncg_flex)
1660         {
1661             for (i = 0; (i < lincsd->nc); i++)
1662             {
1663                 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
1664                 {
1665                     lincsd->bllen[i] = 0;
1666                 }
1667             }
1668         }
1669     }
1670     else
1671     {
1672         /* The OpenMP parallel region of constrain_lincs for derivatives */
1673 #pragma omp parallel num_threads(lincsd->nth)
1674         {
1675             int th = gmx_omp_get_thread_num();
1676
1677             do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
1678                       md->invmass, econq, ir->efep != efepNO ? dvdlambda : NULL,
1679                       bCalcVir, th == 0 ? vir_r_m_dr : lincsd->th[th].vir_r_m_dr);
1680         }
1681     }
1682
1683     if (bCalcVir && lincsd->nth > 1)
1684     {
1685         for (i = 1; i < lincsd->nth; i++)
1686         {
1687             m_add(vir_r_m_dr, lincsd->th[i].vir_r_m_dr, vir_r_m_dr);
1688         }
1689     }
1690
1691     /* count assuming nit=1 */
1692     inc_nrnb(nrnb, eNR_LINCS, lincsd->nc);
1693     inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
1694     if (lincsd->ntriangle > 0)
1695     {
1696         inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
1697     }
1698     if (v)
1699     {
1700         inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc*2);
1701     }
1702     if (bCalcVir)
1703     {
1704         inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc);
1705     }
1706
1707     return bOK;
1708 }