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