d8b1d96b62d13a6c122e5dee8c808483b4c6ae2e
[alexxy/gromacs.git] / src / gromacs / mdlib / lincs.cpp
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,2015,2016,2017,2018, 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 /*! \internal \file
38  * \brief Defines LINCS code.
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \author Mark Abraham <mark.j.abraham@gmail.com>
42  * \ingroup module_mdlib
43  */
44 #include "gmxpre.h"
45
46 #include "lincs.h"
47
48 #include "config.h"
49
50 #include <cassert>
51 #include <cmath>
52 #include <cstdlib>
53
54 #include <algorithm>
55 #include <vector>
56
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/gmxlib/nrnb.h"
60 #include "gromacs/math/functions.h"
61 #include "gromacs/math/units.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/constr.h"
64 #include "gromacs/mdlib/gmx_omp_nthreads.h"
65 #include "gromacs/mdlib/mdrun.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/mdatom.h"
70 #include "gromacs/pbcutil/pbc.h"
71 #include "gromacs/pbcutil/pbc-simd.h"
72 #include "gromacs/simd/simd.h"
73 #include "gromacs/simd/simd_math.h"
74 #include "gromacs/simd/vector_operations.h"
75 #include "gromacs/topology/block.h"
76 #include "gromacs/topology/mtop_util.h"
77 #include "gromacs/utility/arrayref.h"
78 #include "gromacs/utility/basedefinitions.h"
79 #include "gromacs/utility/bitmask.h"
80 #include "gromacs/utility/cstringutil.h"
81 #include "gromacs/utility/exceptions.h"
82 #include "gromacs/utility/fatalerror.h"
83 #include "gromacs/utility/gmxomp.h"
84 #include "gromacs/utility/pleasecite.h"
85 #include "gromacs/utility/smalloc.h"
86
87 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
88
89 namespace gmx
90 {
91
92 //! Unit of work within LINCS.
93 struct Task
94 {
95     //! First constraint for this task.
96     int    b0 = 0;
97     //! b1-1 is the last constraint for this task.
98     int    b1 = 0;
99     //! The number of constraints in triangles.
100     int    ntriangle = 0;
101     //! The list of triangle constraints.
102     int   *triangle = nullptr;
103     //! The bits tell if the matrix element should be used.
104     int   *tri_bits = nullptr;
105     //! Allocation size of triangle and tri_bits.
106     int    tri_alloc = 0;
107     //! Number of indices.
108     int    nind = 0;
109     //! Constraint index for updating atom data.
110     int   *ind = nullptr;
111     //! Number of indices.
112     int    nind_r = 0;
113     //! Constraint index for updating atom data.
114     int   *ind_r = nullptr;
115     //! Allocation size of ind and ind_r.
116     int    ind_nalloc = 0;
117     //! Temporary variable for virial calculation.
118     tensor vir_r_m_dr = {{0}};
119     //! Temporary variable for lambda derivative.
120     real   dhdlambda;
121 };
122
123 /*! \brief Data for LINCS algorithm.
124  */
125 class Lincs
126 {
127     public:
128         //! The global number of constraints.
129         int             ncg = 0;
130         //! The global number of flexible constraints.
131         int             ncg_flex = 0;
132         //! The global number of constraints in triangles.
133         int             ncg_triangle = 0;
134         //! The number of iterations.
135         int             nIter = 0;
136         //! The order of the matrix expansion.
137         int             nOrder = 0;
138         //! The maximum number of constraints connected to a single atom.
139         int             max_connect = 0;
140
141         //! The number of real constraints.
142         int             nc_real = 0;
143         //! The number of constraints including padding for SIMD.
144         int             nc = 0;
145         //! The number we allocated memory for.
146         int             nc_alloc = 0;
147         //! The number of constraint connections.
148         int             ncc = 0;
149         //! The number we allocated memory for.
150         int             ncc_alloc = 0;
151         //! The FE lambda value used for filling blc and blmf.
152         real            matlam = 0;
153         //! mapping from topology to LINCS constraints.
154         int            *con_index = nullptr;
155         //! The reference distance in topology A.
156         real           *bllen0 = nullptr;
157         //! The reference distance in top B - the r.d. in top A.
158         real           *ddist = nullptr;
159         //! The atom pairs involved in the constraints.
160         int            *bla = nullptr;
161         //! 1/sqrt(invmass1  invmass2).
162         real           *blc = nullptr;
163         //! As blc, but with all masses 1.
164         real           *blc1 = nullptr;
165         //! Index into blbnb and blmf.
166         int            *blnr = nullptr;
167         //! List of constraint connections.
168         int            *blbnb = nullptr;
169         //! The local number of constraints in triangles.
170         int             ntriangle = 0;
171         //! The number of constraint connections in triangles.
172         int             ncc_triangle = 0;
173         //! Communicate before each LINCS interation.
174         bool            bCommIter = false;
175         //! Matrix of mass factors for constraint connections.
176         real           *blmf = nullptr;
177         //! As blmf, but with all masses 1.
178         real           *blmf1 = nullptr;
179         //! The reference bond length.
180         real           *bllen = nullptr;
181         //! The local atom count per constraint, can be NULL.
182         int            *nlocat = nullptr;
183
184         /*! \brief The number of tasks used for LINCS work.
185          *
186          * \todo This is mostly used to loop over \c task, which would
187          * be nicer to do with range-based for loops, but the thread
188          * index is used for constructing bit masks and organizing the
189          * virial output buffer, so other things need to change,
190          * first. */
191         int               ntask = 0;
192         /*! \brief LINCS thread division */
193         std::vector<Task> task;
194         //! Atom flags for thread parallelization.
195         gmx_bitmask_t    *atf = nullptr;
196         //! Allocation size of atf
197         int               atf_nalloc = 0;
198         //! Are the LINCS tasks interdependent?
199         bool              bTaskDep = false;
200         //! Are there triangle constraints that cross task borders?
201         bool              bTaskDepTri = false;
202         //! Arrays for temporary storage in the LINCS algorithm.
203         /*! @{ */
204         rvec           *tmpv   = nullptr;
205         real           *tmpncc = nullptr;
206         real           *tmp1   = nullptr;
207         real           *tmp2   = nullptr;
208         real           *tmp3   = nullptr;
209         real           *tmp4   = nullptr;
210         /*! @} */
211         //! The Lagrange multipliers times -1.
212         real               *mlambda = nullptr;
213         //! Storage for the constraint RMS relative deviation output.
214         std::array<real, 2> rmsdData = {{0}};
215 };
216
217 /*! \brief Define simd_width for memory allocation used for SIMD code */
218 #if GMX_SIMD_HAVE_REAL
219 static const int simd_width = GMX_SIMD_REAL_WIDTH;
220 #else
221 static const int simd_width = 1;
222 #endif
223
224 /*! \brief Align to 128 bytes, consistent with the current implementation of
225    AlignedAllocator, which currently forces 128 byte alignment. */
226 static const int align_bytes = 128;
227
228 ArrayRef<real> lincs_rmsdData(Lincs *lincsd)
229 {
230     return lincsd->rmsdData;
231 }
232
233 real lincs_rmsd(const Lincs *lincsd)
234 {
235     if (lincsd->rmsdData[0] > 0)
236     {
237         return std::sqrt(lincsd->rmsdData[1]/lincsd->rmsdData[0]);
238     }
239     else
240     {
241         return 0;
242     }
243 }
244
245 /*! \brief Do a set of nrec LINCS matrix multiplications.
246  *
247  * This function will return with up to date thread-local
248  * constraint data, without an OpenMP barrier.
249  */
250 static void lincs_matrix_expand(const Lincs *lincsd,
251                                 const Task *li_task,
252                                 const real *blcc,
253                                 real *rhs1, real *rhs2, real *sol)
254 {
255     int        b0, b1, nrec, rec;
256     const int *blnr  = lincsd->blnr;
257     const int *blbnb = lincsd->blbnb;
258
259     b0   = li_task->b0;
260     b1   = li_task->b1;
261     nrec = lincsd->nOrder;
262
263     for (rec = 0; rec < nrec; rec++)
264     {
265         int b;
266
267         if (lincsd->bTaskDep)
268         {
269 #pragma omp barrier
270         }
271         for (b = b0; b < b1; b++)
272         {
273             real mvb;
274             int  n;
275
276             mvb = 0;
277             for (n = blnr[b]; n < blnr[b+1]; n++)
278             {
279                 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
280             }
281             rhs2[b] = mvb;
282             sol[b]  = sol[b] + mvb;
283         }
284
285         real *swap;
286
287         swap = rhs1;
288         rhs1 = rhs2;
289         rhs2 = swap;
290     }   /* nrec*(ncons+2*nrtot) flops */
291
292     if (lincsd->ntriangle > 0)
293     {
294         /* Perform an extra nrec recursions for only the constraints
295          * involved in rigid triangles.
296          * In this way their accuracy should come close to those of the other
297          * constraints, since traingles of constraints can produce eigenvalues
298          * around 0.7, while the effective eigenvalue for bond constraints
299          * is around 0.4 (and 0.7*0.7=0.5).
300          */
301
302         if (lincsd->bTaskDep)
303         {
304             /* We need a barrier here, since other threads might still be
305              * reading the contents of rhs1 and/o rhs2.
306              * We could avoid this barrier by introducing two extra rhs
307              * arrays for the triangle constraints only.
308              */
309 #pragma omp barrier
310         }
311
312         /* Constraints involved in a triangle are ensured to be in the same
313          * LINCS task. This means no barriers are required during the extra
314          * iterations for the triangle constraints.
315          */
316         const int *triangle = li_task->triangle;
317         const int *tri_bits = li_task->tri_bits;
318
319         for (rec = 0; rec < nrec; rec++)
320         {
321             int tb;
322
323             for (tb = 0; tb < li_task->ntriangle; tb++)
324             {
325                 int  b, bits, nr0, nr1, n;
326                 real mvb;
327
328                 b    = triangle[tb];
329                 bits = tri_bits[tb];
330                 mvb  = 0;
331                 nr0  = blnr[b];
332                 nr1  = blnr[b+1];
333                 for (n = nr0; n < nr1; n++)
334                 {
335                     if (bits & (1 << (n - nr0)))
336                     {
337                         mvb = mvb + blcc[n]*rhs1[blbnb[n]];
338                     }
339                 }
340                 rhs2[b] = mvb;
341                 sol[b]  = sol[b] + mvb;
342             }
343
344             real *swap;
345
346             swap = rhs1;
347             rhs1 = rhs2;
348             rhs2 = swap;
349         }   /* nrec*(ntriangle + ncc_triangle*2) flops */
350
351         if (lincsd->bTaskDepTri)
352         {
353             /* The constraints triangles are decoupled from each other,
354              * but constraints in one triangle cross thread task borders.
355              * We could probably avoid this with more advanced setup code.
356              */
357 #pragma omp barrier
358         }
359     }
360 }
361
362 //! Update atomic coordinates when an index is not required.
363 static void lincs_update_atoms_noind(int ncons, const int *bla,
364                                      real prefac,
365                                      const real *fac, rvec *r,
366                                      const real *invmass,
367                                      rvec *x)
368 {
369     int  b, i, j;
370     real mvb, im1, im2, tmp0, tmp1, tmp2;
371
372     if (invmass != nullptr)
373     {
374         for (b = 0; b < ncons; b++)
375         {
376             i        = bla[2*b];
377             j        = bla[2*b+1];
378             mvb      = prefac*fac[b];
379             im1      = invmass[i];
380             im2      = invmass[j];
381             tmp0     = r[b][0]*mvb;
382             tmp1     = r[b][1]*mvb;
383             tmp2     = r[b][2]*mvb;
384             x[i][0] -= tmp0*im1;
385             x[i][1] -= tmp1*im1;
386             x[i][2] -= tmp2*im1;
387             x[j][0] += tmp0*im2;
388             x[j][1] += tmp1*im2;
389             x[j][2] += tmp2*im2;
390         }   /* 16 ncons flops */
391     }
392     else
393     {
394         for (b = 0; b < ncons; b++)
395         {
396             i        = bla[2*b];
397             j        = bla[2*b+1];
398             mvb      = prefac*fac[b];
399             tmp0     = r[b][0]*mvb;
400             tmp1     = r[b][1]*mvb;
401             tmp2     = r[b][2]*mvb;
402             x[i][0] -= tmp0;
403             x[i][1] -= tmp1;
404             x[i][2] -= tmp2;
405             x[j][0] += tmp0;
406             x[j][1] += tmp1;
407             x[j][2] += tmp2;
408         }
409     }
410 }
411
412 //! Update atomic coordinates when an index is required.
413 static void lincs_update_atoms_ind(int ncons, const int *ind, const int *bla,
414                                    real prefac,
415                                    const real *fac, rvec *r,
416                                    const real *invmass,
417                                    rvec *x)
418 {
419     int  bi, b, i, j;
420     real mvb, im1, im2, tmp0, tmp1, tmp2;
421
422     if (invmass != nullptr)
423     {
424         for (bi = 0; bi < ncons; bi++)
425         {
426             b        = ind[bi];
427             i        = bla[2*b];
428             j        = bla[2*b+1];
429             mvb      = prefac*fac[b];
430             im1      = invmass[i];
431             im2      = invmass[j];
432             tmp0     = r[b][0]*mvb;
433             tmp1     = r[b][1]*mvb;
434             tmp2     = r[b][2]*mvb;
435             x[i][0] -= tmp0*im1;
436             x[i][1] -= tmp1*im1;
437             x[i][2] -= tmp2*im1;
438             x[j][0] += tmp0*im2;
439             x[j][1] += tmp1*im2;
440             x[j][2] += tmp2*im2;
441         }   /* 16 ncons flops */
442     }
443     else
444     {
445         for (bi = 0; bi < ncons; bi++)
446         {
447             b        = ind[bi];
448             i        = bla[2*b];
449             j        = bla[2*b+1];
450             mvb      = prefac*fac[b];
451             tmp0     = r[b][0]*mvb;
452             tmp1     = r[b][1]*mvb;
453             tmp2     = r[b][2]*mvb;
454             x[i][0] -= tmp0;
455             x[i][1] -= tmp1;
456             x[i][2] -= tmp2;
457             x[j][0] += tmp0;
458             x[j][1] += tmp1;
459             x[j][2] += tmp2;
460         }   /* 16 ncons flops */
461     }
462 }
463
464 //! Update coordinates for atoms.
465 static void lincs_update_atoms(Lincs *li, int th,
466                                real prefac,
467                                const real *fac, rvec *r,
468                                const real *invmass,
469                                rvec *x)
470 {
471     if (li->ntask == 1)
472     {
473         /* Single thread, we simply update for all constraints */
474         lincs_update_atoms_noind(li->nc_real,
475                                  li->bla, prefac, fac, r, invmass, x);
476     }
477     else
478     {
479         /* Update the atom vector components for our thread local
480          * constraints that only access our local atom range.
481          * This can be done without a barrier.
482          */
483         lincs_update_atoms_ind(li->task[th].nind, li->task[th].ind,
484                                li->bla, prefac, fac, r, invmass, x);
485
486         if (li->task[li->ntask].nind > 0)
487         {
488             /* Update the constraints that operate on atoms
489              * in multiple thread atom blocks on the master thread.
490              */
491 #pragma omp barrier
492 #pragma omp master
493             {
494                 lincs_update_atoms_ind(li->task[li->ntask].nind,
495                                        li->task[li->ntask].ind,
496                                        li->bla, prefac, fac, r, invmass, x);
497             }
498         }
499     }
500 }
501
502 #if GMX_SIMD_HAVE_REAL
503 /*! \brief Calculate the constraint distance vectors r to project on from x.
504  *
505  * Determine the right-hand side of the matrix equation using quantity f.
506  * This function only differs from calc_dr_x_xp_simd below in that
507  * no constraint length is subtracted and no PBC is used for f. */
508 static void gmx_simdcall
509 calc_dr_x_f_simd(int                       b0,
510                  int                       b1,
511                  const int *               bla,
512                  const rvec * gmx_restrict x,
513                  const rvec * gmx_restrict f,
514                  const real * gmx_restrict blc,
515                  const real *              pbc_simd,
516                  rvec * gmx_restrict       r,
517                  real * gmx_restrict       rhs,
518                  real * gmx_restrict       sol)
519 {
520     assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
521
522     alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
523
524     for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
525     {
526         offset2[i] = i;
527     }
528
529     for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
530     {
531         SimdReal x0_S, y0_S, z0_S;
532         SimdReal x1_S, y1_S, z1_S;
533         SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
534         SimdReal fx_S, fy_S, fz_S, ip_S, rhs_S;
535         alignas(GMX_SIMD_ALIGNMENT) std::int32_t      offset0[GMX_SIMD_REAL_WIDTH];
536         alignas(GMX_SIMD_ALIGNMENT) std::int32_t      offset1[GMX_SIMD_REAL_WIDTH];
537
538         for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
539         {
540             offset0[i] = bla[bs*2 + i*2];
541             offset1[i] = bla[bs*2 + i*2 + 1];
542         }
543
544         gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
545         gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
546         rx_S = x0_S - x1_S;
547         ry_S = y0_S - y1_S;
548         rz_S = z0_S - z1_S;
549
550         pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
551
552         n2_S  = norm2(rx_S, ry_S, rz_S);
553         il_S  = invsqrt(n2_S);
554
555         rx_S  = rx_S * il_S;
556         ry_S  = ry_S * il_S;
557         rz_S  = rz_S * il_S;
558
559         transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
560
561         gatherLoadUTranspose<3>(reinterpret_cast<const real *>(f), offset0, &x0_S, &y0_S, &z0_S);
562         gatherLoadUTranspose<3>(reinterpret_cast<const real *>(f), offset1, &x1_S, &y1_S, &z1_S);
563         fx_S = x0_S - x1_S;
564         fy_S = y0_S - y1_S;
565         fz_S = z0_S - z1_S;
566
567         ip_S  = iprod(rx_S, ry_S, rz_S, fx_S, fy_S, fz_S);
568
569         rhs_S = load<SimdReal>(blc + bs) * ip_S;
570
571         store(rhs + bs, rhs_S);
572         store(sol + bs, rhs_S);
573     }
574 }
575 #endif // GMX_SIMD_HAVE_REAL
576
577 /*! \brief LINCS projection, works on derivatives of the coordinates. */
578 static void do_lincsp(const rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
579                       Lincs *lincsd, int th,
580                       real *invmass,
581                       ConstraintVariable econq, bool bCalcDHDL,
582                       bool bCalcVir, tensor rmdf)
583 {
584     int      b0, b1, b;
585     int     *bla, *blnr, *blbnb;
586     rvec    *r;
587     real    *blc, *blmf, *blcc, *rhs1, *rhs2, *sol;
588
589     b0 = lincsd->task[th].b0;
590     b1 = lincsd->task[th].b1;
591
592     bla    = lincsd->bla;
593     r      = lincsd->tmpv;
594     blnr   = lincsd->blnr;
595     blbnb  = lincsd->blbnb;
596     if (econq != ConstraintVariable::Force)
597     {
598         /* Use mass-weighted parameters */
599         blc  = lincsd->blc;
600         blmf = lincsd->blmf;
601     }
602     else
603     {
604         /* Use non mass-weighted parameters */
605         blc  = lincsd->blc1;
606         blmf = lincsd->blmf1;
607     }
608     blcc   = lincsd->tmpncc;
609     rhs1   = lincsd->tmp1;
610     rhs2   = lincsd->tmp2;
611     sol    = lincsd->tmp3;
612
613 #if GMX_SIMD_HAVE_REAL
614     /* This SIMD code does the same as the plain-C code after the #else.
615      * The only difference is that we always call pbc code, as with SIMD
616      * the overhead of pbc computation (when not needed) is small.
617      */
618     alignas(GMX_SIMD_ALIGNMENT) real    pbc_simd[9*GMX_SIMD_REAL_WIDTH];
619
620     /* Convert the pbc struct for SIMD */
621     set_pbc_simd(pbc, pbc_simd);
622
623     /* Compute normalized x i-j vectors, store in r.
624      * Compute the inner product of r and xp i-j and store in rhs1.
625      */
626     calc_dr_x_f_simd(b0, b1, bla, x, f, blc,
627                      pbc_simd,
628                      r, rhs1, sol);
629
630 #else   // GMX_SIMD_HAVE_REAL
631
632     /* Compute normalized i-j vectors */
633     if (pbc)
634     {
635         for (b = b0; b < b1; b++)
636         {
637             rvec dx;
638
639             pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
640             unitv(dx, r[b]);
641         }
642     }
643     else
644     {
645         for (b = b0; b < b1; b++)
646         {
647             rvec dx;
648
649             rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
650             unitv(dx, r[b]);
651         }   /* 16 ncons flops */
652     }
653
654     for (b = b0; b < b1; b++)
655     {
656         int  i, j;
657         real mvb;
658
659         i       = bla[2*b];
660         j       = bla[2*b+1];
661         mvb     = blc[b]*(r[b][0]*(f[i][0] - f[j][0]) +
662                           r[b][1]*(f[i][1] - f[j][1]) +
663                           r[b][2]*(f[i][2] - f[j][2]));
664         rhs1[b] = mvb;
665         sol[b]  = mvb;
666         /* 7 flops */
667     }
668
669 #endif  // GMX_SIMD_HAVE_REAL
670
671     if (lincsd->bTaskDep)
672     {
673         /* We need a barrier, since the matrix construction below
674          * can access entries in r of other threads.
675          */
676 #pragma omp barrier
677     }
678
679     /* Construct the (sparse) LINCS matrix */
680     for (b = b0; b < b1; b++)
681     {
682         int n;
683
684         for (n = blnr[b]; n < blnr[b+1]; n++)
685         {
686             blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
687         }   /* 6 nr flops */
688     }
689     /* Together: 23*ncons + 6*nrtot flops */
690
691     lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
692     /* nrec*(ncons+2*nrtot) flops */
693
694     if (econq == ConstraintVariable::Deriv_FlexCon)
695     {
696         /* We only want to constraint the flexible constraints,
697          * so we mask out the normal ones by setting sol to 0.
698          */
699         for (b = b0; b < b1; b++)
700         {
701             if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
702             {
703                 sol[b] = 0;
704             }
705         }
706     }
707
708     /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
709     for (b = b0; b < b1; b++)
710     {
711         sol[b] *= blc[b];
712     }
713
714     /* When constraining forces, we should not use mass weighting,
715      * so we pass invmass=NULL, which results in the use of 1 for all atoms.
716      */
717     lincs_update_atoms(lincsd, th, 1.0, sol, r,
718                        (econq != ConstraintVariable::Force) ? invmass : nullptr, fp);
719
720     if (bCalcDHDL)
721     {
722         real dhdlambda;
723
724         dhdlambda = 0;
725         for (b = b0; b < b1; b++)
726         {
727             dhdlambda -= sol[b]*lincsd->ddist[b];
728         }
729
730         lincsd->task[th].dhdlambda = dhdlambda;
731     }
732
733     if (bCalcVir)
734     {
735         /* Constraint virial,
736          * determines sum r_bond x delta f,
737          * where delta f is the constraint correction
738          * of the quantity that is being constrained.
739          */
740         for (b = b0; b < b1; b++)
741         {
742             real mvb, tmp1;
743             int  i, j;
744
745             mvb = lincsd->bllen[b]*sol[b];
746             for (i = 0; i < DIM; i++)
747             {
748                 tmp1 = mvb*r[b][i];
749                 for (j = 0; j < DIM; j++)
750                 {
751                     rmdf[i][j] += tmp1*r[b][j];
752                 }
753             }
754         }   /* 23 ncons flops */
755     }
756 }
757
758 #if GMX_SIMD_HAVE_REAL
759 /*! \brief Calculate the constraint distance vectors r to project on from x.
760  *
761  * Determine the right-hand side of the matrix equation using coordinates xp. */
762 static void gmx_simdcall
763 calc_dr_x_xp_simd(int                       b0,
764                   int                       b1,
765                   const int *               bla,
766                   const rvec * gmx_restrict x,
767                   const rvec * gmx_restrict xp,
768                   const real * gmx_restrict bllen,
769                   const real * gmx_restrict blc,
770                   const real *              pbc_simd,
771                   rvec * gmx_restrict       r,
772                   real * gmx_restrict       rhs,
773                   real * gmx_restrict       sol)
774 {
775     assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
776     alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
777
778     for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
779     {
780         offset2[i] = i;
781     }
782
783     for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
784     {
785         SimdReal x0_S, y0_S, z0_S;
786         SimdReal x1_S, y1_S, z1_S;
787         SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
788         SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
789         alignas(GMX_SIMD_ALIGNMENT) std::int32_t      offset0[GMX_SIMD_REAL_WIDTH];
790         alignas(GMX_SIMD_ALIGNMENT) std::int32_t      offset1[GMX_SIMD_REAL_WIDTH];
791
792         for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
793         {
794             offset0[i] = bla[bs*2 + i*2];
795             offset1[i] = bla[bs*2 + i*2 + 1];
796         }
797
798         gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
799         gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
800         rx_S = x0_S - x1_S;
801         ry_S = y0_S - y1_S;
802         rz_S = z0_S - z1_S;
803
804         pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
805
806         n2_S  = norm2(rx_S, ry_S, rz_S);
807         il_S  = invsqrt(n2_S);
808
809         rx_S  = rx_S * il_S;
810         ry_S  = ry_S * il_S;
811         rz_S  = rz_S * il_S;
812
813         transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
814
815         gatherLoadUTranspose<3>(reinterpret_cast<const real *>(xp), offset0, &x0_S, &y0_S, &z0_S);
816         gatherLoadUTranspose<3>(reinterpret_cast<const real *>(xp), offset1, &x1_S, &y1_S, &z1_S);
817         rxp_S = x0_S - x1_S;
818         ryp_S = y0_S - y1_S;
819         rzp_S = z0_S - z1_S;
820
821         pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
822
823         ip_S  = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
824
825         rhs_S = load<SimdReal>(blc + bs) * (ip_S - load<SimdReal>(bllen + bs));
826
827         store(rhs + bs, rhs_S);
828         store(sol + bs, rhs_S);
829     }
830 }
831 #endif // GMX_SIMD_HAVE_REAL
832
833 /*! \brief Determine the distances and right-hand side for the next iteration. */
834 gmx_unused static void calc_dist_iter(
835         int                       b0,
836         int                       b1,
837         const int *               bla,
838         const rvec * gmx_restrict xp,
839         const real * gmx_restrict bllen,
840         const real * gmx_restrict blc,
841         const t_pbc *             pbc,
842         real                      wfac,
843         real * gmx_restrict       rhs,
844         real * gmx_restrict       sol,
845         bool *                    bWarn)
846 {
847     int b;
848
849     for (b = b0; b < b1; b++)
850     {
851         real len, len2, dlen2, mvb;
852         rvec dx;
853
854         len = bllen[b];
855         if (pbc)
856         {
857             pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
858         }
859         else
860         {
861             rvec_sub(xp[bla[2*b]], xp[bla[2*b+1]], dx);
862         }
863         len2  = len*len;
864         dlen2 = 2*len2 - ::norm2(dx);
865         if (dlen2 < wfac*len2)
866         {
867             /* not race free - see detailed comment in caller */
868             *bWarn = TRUE;
869         }
870         if (dlen2 > 0)
871         {
872             mvb = blc[b]*(len - dlen2*gmx::invsqrt(dlen2));
873         }
874         else
875         {
876             mvb = blc[b]*len;
877         }
878         rhs[b]  = mvb;
879         sol[b]  = mvb;
880     }   /* 20*ncons flops */
881 }
882
883 #if GMX_SIMD_HAVE_REAL
884 /*! \brief As calc_dist_iter(), but using SIMD intrinsics. */
885 static void gmx_simdcall
886 calc_dist_iter_simd(int                       b0,
887                     int                       b1,
888                     const int *               bla,
889                     const rvec * gmx_restrict x,
890                     const real * gmx_restrict bllen,
891                     const real * gmx_restrict blc,
892                     const real *              pbc_simd,
893                     real                      wfac,
894                     real * gmx_restrict       rhs,
895                     real * gmx_restrict       sol,
896                     bool *                    bWarn)
897 {
898     SimdReal        min_S(GMX_REAL_MIN);
899     SimdReal        two_S(2.0);
900     SimdReal        wfac_S(wfac);
901     SimdBool        warn_S;
902
903     int             bs;
904
905     assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
906
907     /* Initialize all to FALSE */
908     warn_S = (two_S < setZero());
909
910     for (bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
911     {
912         SimdReal x0_S, y0_S, z0_S;
913         SimdReal x1_S, y1_S, z1_S;
914         SimdReal rx_S, ry_S, rz_S, n2_S;
915         SimdReal len_S, len2_S, dlen2_S, lc_S, blc_S;
916         alignas(GMX_SIMD_ALIGNMENT) std::int32_t      offset0[GMX_SIMD_REAL_WIDTH];
917         alignas(GMX_SIMD_ALIGNMENT) std::int32_t      offset1[GMX_SIMD_REAL_WIDTH];
918
919         for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
920         {
921             offset0[i] = bla[bs*2 + i*2];
922             offset1[i] = bla[bs*2 + i*2 + 1];
923         }
924
925         gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
926         gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
927         rx_S = x0_S - x1_S;
928         ry_S = y0_S - y1_S;
929         rz_S = z0_S - z1_S;
930
931         pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
932
933         n2_S    = norm2(rx_S, ry_S, rz_S);
934
935         len_S   = load<SimdReal>(bllen + bs);
936         len2_S  = len_S * len_S;
937
938         dlen2_S = fms(two_S, len2_S, n2_S);
939
940         warn_S  = warn_S || (dlen2_S < (wfac_S * len2_S));
941
942         /* Avoid 1/0 by taking the max with REAL_MIN.
943          * Note: when dlen2 is close to zero (90 degree constraint rotation),
944          * the accuracy of the algorithm is no longer relevant.
945          */
946         dlen2_S = max(dlen2_S, min_S);
947
948         lc_S    = fnma(dlen2_S, invsqrt(dlen2_S), len_S);
949
950         blc_S   = load<SimdReal>(blc + bs);
951
952         lc_S    = blc_S * lc_S;
953
954         store(rhs + bs, lc_S);
955         store(sol + bs, lc_S);
956     }
957
958     if (anyTrue(warn_S))
959     {
960         *bWarn = TRUE;
961     }
962 }
963 #endif // GMX_SIMD_HAVE_REAL
964
965 //! Implements LINCS constraining.
966 static void do_lincs(const rvec *x, rvec *xp, matrix box, t_pbc *pbc,
967                      Lincs *lincsd, int th,
968                      const real *invmass,
969                      const t_commrec *cr,
970                      bool bCalcDHDL,
971                      real wangle, bool *bWarn,
972                      real invdt, rvec * gmx_restrict v,
973                      bool bCalcVir, tensor vir_r_m_dr)
974 {
975     int      b0, b1, b, i, j, n, iter;
976     int     *bla, *blnr, *blbnb;
977     rvec    *r;
978     real    *blc, *blmf, *bllen, *blcc, *rhs1, *rhs2, *sol, *blc_sol, *mlambda;
979     int     *nlocat;
980
981     b0 = lincsd->task[th].b0;
982     b1 = lincsd->task[th].b1;
983
984     bla     = lincsd->bla;
985     r       = lincsd->tmpv;
986     blnr    = lincsd->blnr;
987     blbnb   = lincsd->blbnb;
988     blc     = lincsd->blc;
989     blmf    = lincsd->blmf;
990     bllen   = lincsd->bllen;
991     blcc    = lincsd->tmpncc;
992     rhs1    = lincsd->tmp1;
993     rhs2    = lincsd->tmp2;
994     sol     = lincsd->tmp3;
995     blc_sol = lincsd->tmp4;
996     mlambda = lincsd->mlambda;
997     nlocat  = lincsd->nlocat;
998
999 #if GMX_SIMD_HAVE_REAL
1000
1001     /* This SIMD code does the same as the plain-C code after the #else.
1002      * The only difference is that we always call pbc code, as with SIMD
1003      * the overhead of pbc computation (when not needed) is small.
1004      */
1005     alignas(GMX_SIMD_ALIGNMENT) real    pbc_simd[9*GMX_SIMD_REAL_WIDTH];
1006
1007     /* Convert the pbc struct for SIMD */
1008     set_pbc_simd(pbc, pbc_simd);
1009
1010     /* Compute normalized x i-j vectors, store in r.
1011      * Compute the inner product of r and xp i-j and store in rhs1.
1012      */
1013     calc_dr_x_xp_simd(b0, b1, bla, x, xp, bllen, blc,
1014                       pbc_simd,
1015                       r, rhs1, sol);
1016
1017 #else   // GMX_SIMD_HAVE_REAL
1018
1019     if (pbc)
1020     {
1021         /* Compute normalized i-j vectors */
1022         for (b = b0; b < b1; b++)
1023         {
1024             rvec dx;
1025             real mvb;
1026
1027             pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
1028             unitv(dx, r[b]);
1029
1030             pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
1031             mvb     = blc[b]*(::iprod(r[b], dx) - bllen[b]);
1032             rhs1[b] = mvb;
1033             sol[b]  = mvb;
1034         }
1035     }
1036     else
1037     {
1038         /* Compute normalized i-j vectors */
1039         for (b = b0; b < b1; b++)
1040         {
1041             real tmp0, tmp1, tmp2, rlen, mvb;
1042
1043             i       = bla[2*b];
1044             j       = bla[2*b+1];
1045             tmp0    = x[i][0] - x[j][0];
1046             tmp1    = x[i][1] - x[j][1];
1047             tmp2    = x[i][2] - x[j][2];
1048             rlen    = gmx::invsqrt(tmp0*tmp0 + tmp1*tmp1 + tmp2*tmp2);
1049             r[b][0] = rlen*tmp0;
1050             r[b][1] = rlen*tmp1;
1051             r[b][2] = rlen*tmp2;
1052             /* 16 ncons flops */
1053
1054             i       = bla[2*b];
1055             j       = bla[2*b+1];
1056             mvb     = blc[b]*(r[b][0]*(xp[i][0] - xp[j][0]) +
1057                               r[b][1]*(xp[i][1] - xp[j][1]) +
1058                               r[b][2]*(xp[i][2] - xp[j][2]) - bllen[b]);
1059             rhs1[b] = mvb;
1060             sol[b]  = mvb;
1061             /* 10 flops */
1062         }
1063         /* Together: 26*ncons + 6*nrtot flops */
1064     }
1065
1066 #endif  // GMX_SIMD_HAVE_REAL
1067
1068     if (lincsd->bTaskDep)
1069     {
1070         /* We need a barrier, since the matrix construction below
1071          * can access entries in r of other threads.
1072          */
1073 #pragma omp barrier
1074     }
1075
1076     /* Construct the (sparse) LINCS matrix */
1077     for (b = b0; b < b1; b++)
1078     {
1079         for (n = blnr[b]; n < blnr[b+1]; n++)
1080         {
1081             blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
1082         }
1083     }
1084     /* Together: 26*ncons + 6*nrtot flops */
1085
1086     lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1087     /* nrec*(ncons+2*nrtot) flops */
1088
1089 #if GMX_SIMD_HAVE_REAL
1090     for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1091     {
1092         SimdReal t1 = load<SimdReal>(blc + b);
1093         SimdReal t2 = load<SimdReal>(sol + b);
1094         store(mlambda + b, t1 * t2);
1095     }
1096 #else
1097     for (b = b0; b < b1; b++)
1098     {
1099         mlambda[b] = blc[b]*sol[b];
1100     }
1101 #endif  // GMX_SIMD_HAVE_REAL
1102
1103     /* Update the coordinates */
1104     lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1105
1106     /*
1107      ********  Correction for centripetal effects  ********
1108      */
1109
1110     real wfac;
1111
1112     wfac = std::cos(DEG2RAD*wangle);
1113     wfac = wfac*wfac;
1114
1115     for (iter = 0; iter < lincsd->nIter; iter++)
1116     {
1117         if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
1118         {
1119 #pragma omp barrier
1120 #pragma omp master
1121             {
1122                 /* Communicate the corrected non-local coordinates */
1123                 if (DOMAINDECOMP(cr))
1124                 {
1125                     dd_move_x_constraints(cr->dd, box, xp, nullptr, FALSE);
1126                 }
1127             }
1128 #pragma omp barrier
1129         }
1130         else if (lincsd->bTaskDep)
1131         {
1132 #pragma omp barrier
1133         }
1134
1135 #if GMX_SIMD_HAVE_REAL
1136         calc_dist_iter_simd(b0, b1, bla, xp, bllen, blc, pbc_simd, wfac,
1137                             rhs1, sol, bWarn);
1138 #else
1139         calc_dist_iter(b0, b1, bla, xp, bllen, blc, pbc, wfac,
1140                        rhs1, sol, bWarn);
1141         /* 20*ncons flops */
1142 #endif      // GMX_SIMD_HAVE_REAL
1143
1144         lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1145         /* nrec*(ncons+2*nrtot) flops */
1146
1147 #if GMX_SIMD_HAVE_REAL
1148         for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1149         {
1150             SimdReal t1  = load<SimdReal>(blc + b);
1151             SimdReal t2  = load<SimdReal>(sol + b);
1152             SimdReal mvb = t1 * t2;
1153             store(blc_sol + b, mvb);
1154             store(mlambda + b, load<SimdReal>(mlambda + b) + mvb);
1155         }
1156 #else
1157         for (b = b0; b < b1; b++)
1158         {
1159             real mvb;
1160
1161             mvb         = blc[b]*sol[b];
1162             blc_sol[b]  = mvb;
1163             mlambda[b] += mvb;
1164         }
1165 #endif      // GMX_SIMD_HAVE_REAL
1166
1167         /* Update the coordinates */
1168         lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1169     }
1170     /* nit*ncons*(37+9*nrec) flops */
1171
1172     if (v != nullptr)
1173     {
1174         /* Update the velocities */
1175         lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1176         /* 16 ncons flops */
1177     }
1178
1179     if (nlocat != nullptr && (bCalcDHDL || bCalcVir))
1180     {
1181         if (lincsd->bTaskDep)
1182         {
1183             /* In lincs_update_atoms threads might cross-read mlambda */
1184 #pragma omp barrier
1185         }
1186
1187         /* Only account for local atoms */
1188         for (b = b0; b < b1; b++)
1189         {
1190             mlambda[b] *= 0.5*nlocat[b];
1191         }
1192     }
1193
1194     if (bCalcDHDL)
1195     {
1196         real dhdl;
1197
1198         dhdl = 0;
1199         for (b = b0; b < b1; b++)
1200         {
1201             /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1202              * later after the contributions are reduced over the threads.
1203              */
1204             dhdl -= lincsd->mlambda[b]*lincsd->ddist[b];
1205         }
1206         lincsd->task[th].dhdlambda = dhdl;
1207     }
1208
1209     if (bCalcVir)
1210     {
1211         /* Constraint virial */
1212         for (b = b0; b < b1; b++)
1213         {
1214             real tmp0, tmp1;
1215
1216             tmp0 = -bllen[b]*mlambda[b];
1217             for (i = 0; i < DIM; i++)
1218             {
1219                 tmp1 = tmp0*r[b][i];
1220                 for (j = 0; j < DIM; j++)
1221                 {
1222                     vir_r_m_dr[i][j] -= tmp1*r[b][j];
1223                 }
1224             }
1225         }   /* 22 ncons flops */
1226     }
1227
1228     /* Total:
1229      * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1230      * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1231      *
1232      * (26+nrec)*ncons + (6+2*nrec)*nrtot
1233      * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1234      * if nit=1
1235      * (63+nrec)*ncons + (6+4*nrec)*nrtot
1236      */
1237 }
1238
1239 /*! \brief Sets the elements in the LINCS matrix for task task. */
1240 static void set_lincs_matrix_task(Lincs                *li,
1241                                   Task                 *li_task,
1242                                   const real           *invmass,
1243                                   int                  *ncc_triangle,
1244                                   int                  *nCrossTaskTriangles)
1245 {
1246     int        i;
1247
1248     /* Construct the coupling coefficient matrix blmf */
1249     li_task->ntriangle   = 0;
1250     *ncc_triangle        = 0;
1251     *nCrossTaskTriangles = 0;
1252     for (i = li_task->b0; i < li_task->b1; i++)
1253     {
1254         int a1, a2, n;
1255
1256         a1 = li->bla[2*i];
1257         a2 = li->bla[2*i+1];
1258         for (n = li->blnr[i]; (n < li->blnr[i+1]); n++)
1259         {
1260             int k, sign, center, end;
1261
1262             k = li->blbnb[n];
1263
1264             /* If we are using multiple, independent tasks for LINCS,
1265              * the calls to check_assign_connected should have
1266              * put all connected constraints in our task.
1267              */
1268             assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1269
1270             if (a1 == li->bla[2*k] || a2 == li->bla[2*k+1])
1271             {
1272                 sign = -1;
1273             }
1274             else
1275             {
1276                 sign = 1;
1277             }
1278             if (a1 == li->bla[2*k] || a1 == li->bla[2*k+1])
1279             {
1280                 center = a1;
1281                 end    = a2;
1282             }
1283             else
1284             {
1285                 center = a2;
1286                 end    = a1;
1287             }
1288             li->blmf[n]  = sign*invmass[center]*li->blc[i]*li->blc[k];
1289             li->blmf1[n] = sign*0.5;
1290             if (li->ncg_triangle > 0)
1291             {
1292                 int nk, kk;
1293
1294                 /* Look for constraint triangles */
1295                 for (nk = li->blnr[k]; (nk < li->blnr[k+1]); nk++)
1296                 {
1297                     kk = li->blbnb[nk];
1298                     if (kk != i && kk != k &&
1299                         (li->bla[2*kk] == end || li->bla[2*kk+1] == end))
1300                     {
1301                         /* Check if the constraints in this triangle actually
1302                          * belong to a different task. We still assign them
1303                          * here, since it's convenient for the triangle
1304                          * iterations, but we then need an extra barrier.
1305                          */
1306                         if (k  < li_task->b0 || k  >= li_task->b1 ||
1307                             kk < li_task->b0 || kk >= li_task->b1)
1308                         {
1309                             (*nCrossTaskTriangles)++;
1310                         }
1311
1312                         if (li_task->ntriangle == 0 ||
1313                             li_task->triangle[li_task->ntriangle - 1] < i)
1314                         {
1315                             /* Add this constraint to the triangle list */
1316                             li_task->triangle[li_task->ntriangle] = i;
1317                             li_task->tri_bits[li_task->ntriangle] = 0;
1318                             li_task->ntriangle++;
1319                             if (li->blnr[i+1] - li->blnr[i] > static_cast<int>(sizeof(li_task->tri_bits[0])*8 - 1))
1320                             {
1321                                 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %lu allowed for constraints participating in triangles",
1322                                           li->blnr[i+1] - li->blnr[i],
1323                                           sizeof(li_task->tri_bits[0])*8-1);
1324                             }
1325                         }
1326                         li_task->tri_bits[li_task->ntriangle-1] |= (1 << (n - li->blnr[i]));
1327                         (*ncc_triangle)++;
1328                     }
1329                 }
1330             }
1331         }
1332     }
1333 }
1334
1335 /*! \brief Sets the elements in the LINCS matrix. */
1336 static void set_lincs_matrix(Lincs *li, real *invmass, real lambda)
1337 {
1338     int        i;
1339     const real invsqrt2 = 0.7071067811865475244;
1340
1341     for (i = 0; (i < li->nc); i++)
1342     {
1343         int a1, a2;
1344
1345         a1          = li->bla[2*i];
1346         a2          = li->bla[2*i+1];
1347         li->blc[i]  = gmx::invsqrt(invmass[a1] + invmass[a2]);
1348         li->blc1[i] = invsqrt2;
1349     }
1350
1351     /* Construct the coupling coefficient matrix blmf */
1352     int th, ntriangle = 0, ncc_triangle = 0, nCrossTaskTriangles = 0;
1353 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle, nCrossTaskTriangles) num_threads(li->ntask) schedule(static)
1354     for (th = 0; th < li->ntask; th++)
1355     {
1356         try
1357         {
1358             set_lincs_matrix_task(li, &li->task[th], invmass,
1359                                   &ncc_triangle, &nCrossTaskTriangles);
1360             ntriangle = li->task[th].ntriangle;
1361         }
1362         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1363     }
1364     li->ntriangle    = ntriangle;
1365     li->ncc_triangle = ncc_triangle;
1366     li->bTaskDepTri  = (nCrossTaskTriangles > 0);
1367
1368     if (debug)
1369     {
1370         fprintf(debug, "The %d constraints participate in %d triangles\n",
1371                 li->nc, li->ntriangle);
1372         fprintf(debug, "There are %d constraint couplings, of which %d in triangles\n",
1373                 li->ncc, li->ncc_triangle);
1374         if (li->ntriangle > 0 && li->ntask > 1)
1375         {
1376             fprintf(debug, "%d constraint triangles contain constraints assigned to different tasks\n",
1377                     nCrossTaskTriangles);
1378         }
1379     }
1380
1381     /* Set matlam,
1382      * so we know with which lambda value the masses have been set.
1383      */
1384     li->matlam = lambda;
1385 }
1386
1387 //! Finds all triangles of atoms that share constraints to a central atom.
1388 static int count_triangle_constraints(const t_ilist  *ilist,
1389                                       const t_blocka &at2con)
1390 {
1391     int      ncon1, ncon_tot;
1392     int      c0, a00, a01, n1, c1, a10, a11, ac1, n2, c2, a20, a21;
1393     int      ncon_triangle;
1394     bool     bTriangle;
1395     t_iatom *ia1, *ia2, *iap;
1396
1397     ncon1    = ilist[F_CONSTR].nr/3;
1398     ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
1399
1400     ia1 = ilist[F_CONSTR].iatoms;
1401     ia2 = ilist[F_CONSTRNC].iatoms;
1402
1403     ncon_triangle = 0;
1404     for (c0 = 0; c0 < ncon_tot; c0++)
1405     {
1406         bTriangle = FALSE;
1407         iap       = constr_iatomptr(ncon1, ia1, ia2, c0);
1408         a00       = iap[1];
1409         a01       = iap[2];
1410         for (n1 = at2con.index[a01]; n1 < at2con.index[a01+1]; n1++)
1411         {
1412             c1 = at2con.a[n1];
1413             if (c1 != c0)
1414             {
1415                 iap = constr_iatomptr(ncon1, ia1, ia2, c1);
1416                 a10 = iap[1];
1417                 a11 = iap[2];
1418                 if (a10 == a01)
1419                 {
1420                     ac1 = a11;
1421                 }
1422                 else
1423                 {
1424                     ac1 = a10;
1425                 }
1426                 for (n2 = at2con.index[ac1]; n2 < at2con.index[ac1+1]; n2++)
1427                 {
1428                     c2 = at2con.a[n2];
1429                     if (c2 != c0 && c2 != c1)
1430                     {
1431                         iap = constr_iatomptr(ncon1, ia1, ia2, c2);
1432                         a20 = iap[1];
1433                         a21 = iap[2];
1434                         if (a20 == a00 || a21 == a00)
1435                         {
1436                             bTriangle = TRUE;
1437                         }
1438                     }
1439                 }
1440             }
1441         }
1442         if (bTriangle)
1443         {
1444             ncon_triangle++;
1445         }
1446     }
1447
1448     return ncon_triangle;
1449 }
1450
1451 //! Finds sequences of sequential constraints.
1452 static bool more_than_two_sequential_constraints(const t_ilist  *ilist,
1453                                                  const t_blocka &at2con)
1454 {
1455     t_iatom  *ia1, *ia2, *iap;
1456     int       ncon1, ncon_tot, c;
1457     int       a1, a2;
1458     bool      bMoreThanTwoSequentialConstraints;
1459
1460     ncon1    = ilist[F_CONSTR].nr/3;
1461     ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
1462
1463     ia1 = ilist[F_CONSTR].iatoms;
1464     ia2 = ilist[F_CONSTRNC].iatoms;
1465
1466     bMoreThanTwoSequentialConstraints = FALSE;
1467     for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
1468     {
1469         iap = constr_iatomptr(ncon1, ia1, ia2, c);
1470         a1  = iap[1];
1471         a2  = iap[2];
1472         /* Check if this constraint has constraints connected at both atoms */
1473         if (at2con.index[a1+1] - at2con.index[a1] > 1 &&
1474             at2con.index[a2+1] - at2con.index[a2] > 1)
1475         {
1476             bMoreThanTwoSequentialConstraints = TRUE;
1477         }
1478     }
1479
1480     return bMoreThanTwoSequentialConstraints;
1481 }
1482
1483 Lincs *init_lincs(FILE *fplog, const gmx_mtop_t &mtop,
1484                   int nflexcon_global, ArrayRef<const t_blocka> at2con,
1485                   bool bPLINCS, int nIter, int nProjOrder)
1486 {
1487     // TODO this should become a unique_ptr
1488     Lincs                *li;
1489     bool                  bMoreThanTwoSeq;
1490
1491     if (fplog)
1492     {
1493         fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
1494                 bPLINCS ? " Parallel" : "");
1495     }
1496
1497     li = new Lincs;
1498
1499     li->ncg      =
1500         gmx_mtop_ftype_count(mtop, F_CONSTR) +
1501         gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1502     li->ncg_flex = nflexcon_global;
1503
1504     li->nIter  = nIter;
1505     li->nOrder = nProjOrder;
1506
1507     li->max_connect = 0;
1508     for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1509     {
1510         for (int a = 0; a < mtop.moltype[mt].atoms.nr; a++)
1511         {
1512             li->max_connect = std::max(li->max_connect,
1513                                        at2con[mt].index[a + 1] - at2con[mt].index[a]);
1514         }
1515     }
1516
1517     li->ncg_triangle = 0;
1518     bMoreThanTwoSeq  = FALSE;
1519     for (const gmx_molblock_t &molb : mtop.molblock)
1520     {
1521         const gmx_moltype_t &molt = mtop.moltype[molb.type];
1522
1523         li->ncg_triangle +=
1524             molb.nmol*
1525             count_triangle_constraints(molt.ilist, at2con[molb.type]);
1526
1527         if (!bMoreThanTwoSeq &&
1528             more_than_two_sequential_constraints(molt.ilist, at2con[molb.type]))
1529         {
1530             bMoreThanTwoSeq = TRUE;
1531         }
1532     }
1533
1534     /* Check if we need to communicate not only before LINCS,
1535      * but also before each iteration.
1536      * The check for only two sequential constraints is only
1537      * useful for the common case of H-bond only constraints.
1538      * With more effort we could also make it useful for small
1539      * molecules with nr. sequential constraints <= nOrder-1.
1540      */
1541     li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1542
1543     if (debug && bPLINCS)
1544     {
1545         fprintf(debug, "PLINCS communication before each iteration: %d\n",
1546                 static_cast<int>(li->bCommIter));
1547     }
1548
1549     /* LINCS can run on any number of threads.
1550      * Currently the number is fixed for the whole simulation,
1551      * but it could be set in set_lincs().
1552      * The current constraint to task assignment code can create independent
1553      * tasks only when not more than two constraints are connected sequentially.
1554      */
1555     li->ntask    = gmx_omp_nthreads_get(emntLINCS);
1556     li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1557     if (debug)
1558     {
1559         fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n",
1560                 li->ntask, li->bTaskDep ? "" : "in");
1561     }
1562     if (li->ntask == 1)
1563     {
1564         li->task.resize(1);
1565     }
1566     else
1567     {
1568         /* Allocate an extra elements for "task-overlap" constraints */
1569         li->task.resize(li->ntask + 1);
1570     }
1571
1572     if (bPLINCS || li->ncg_triangle > 0)
1573     {
1574         please_cite(fplog, "Hess2008a");
1575     }
1576     else
1577     {
1578         please_cite(fplog, "Hess97a");
1579     }
1580
1581     if (fplog)
1582     {
1583         fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1584         if (bPLINCS)
1585         {
1586             fprintf(fplog, "There are inter charge-group constraints,\n"
1587                     "will communicate selected coordinates each lincs iteration\n");
1588         }
1589         if (li->ncg_triangle > 0)
1590         {
1591             fprintf(fplog,
1592                     "%d constraints are involved in constraint triangles,\n"
1593                     "will apply an additional matrix expansion of order %d for couplings\n"
1594                     "between constraints inside triangles\n",
1595                     li->ncg_triangle, li->nOrder);
1596         }
1597     }
1598
1599     return li;
1600 }
1601
1602 void done_lincs(Lincs *li)
1603 {
1604     delete li;
1605 }
1606
1607 /*! \brief Sets up the work division over the threads. */
1608 static void lincs_thread_setup(Lincs *li, int natoms)
1609 {
1610     Task           *li_m;
1611     int             th;
1612     gmx_bitmask_t  *atf;
1613     int             a;
1614
1615     if (natoms > li->atf_nalloc)
1616     {
1617         li->atf_nalloc = over_alloc_large(natoms);
1618         srenew(li->atf, li->atf_nalloc);
1619     }
1620
1621     atf = li->atf;
1622     /* Clear the atom flags */
1623     for (a = 0; a < natoms; a++)
1624     {
1625         bitmask_clear(&atf[a]);
1626     }
1627
1628     if (li->ntask > BITMASK_SIZE)
1629     {
1630         gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1631     }
1632
1633     for (th = 0; th < li->ntask; th++)
1634     {
1635         Task         *li_task;
1636         int           b;
1637
1638         li_task = &li->task[th];
1639
1640         /* For each atom set a flag for constraints from each */
1641         for (b = li_task->b0; b < li_task->b1; b++)
1642         {
1643             bitmask_set_bit(&atf[li->bla[b*2    ]], th);
1644             bitmask_set_bit(&atf[li->bla[b*2 + 1]], th);
1645         }
1646     }
1647
1648 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1649     for (th = 0; th < li->ntask; th++)
1650     {
1651         try
1652         {
1653             Task          *li_task;
1654             gmx_bitmask_t  mask;
1655             int            b;
1656
1657             li_task = &li->task[th];
1658
1659             if (li_task->b1 - li_task->b0 > li_task->ind_nalloc)
1660             {
1661                 li_task->ind_nalloc = over_alloc_large(li_task->b1-li_task->b0);
1662                 srenew(li_task->ind, li_task->ind_nalloc);
1663                 srenew(li_task->ind_r, li_task->ind_nalloc);
1664             }
1665
1666             bitmask_init_low_bits(&mask, th);
1667
1668             li_task->nind   = 0;
1669             li_task->nind_r = 0;
1670             for (b = li_task->b0; b < li_task->b1; b++)
1671             {
1672                 /* We let the constraint with the lowest thread index
1673                  * operate on atoms with constraints from multiple threads.
1674                  */
1675                 if (bitmask_is_disjoint(atf[li->bla[b*2]], mask) &&
1676                     bitmask_is_disjoint(atf[li->bla[b*2+1]], mask))
1677                 {
1678                     /* Add the constraint to the local atom update index */
1679                     li_task->ind[li_task->nind++] = b;
1680                 }
1681                 else
1682                 {
1683                     /* Add the constraint to the rest block */
1684                     li_task->ind_r[li_task->nind_r++] = b;
1685                 }
1686             }
1687         }
1688         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1689     }
1690
1691     /* We need to copy all constraints which have not be assigned
1692      * to a thread to a separate list which will be handled by one thread.
1693      */
1694     li_m = &li->task[li->ntask];
1695
1696     li_m->nind = 0;
1697     for (th = 0; th < li->ntask; th++)
1698     {
1699         Task         *li_task;
1700         int           b;
1701
1702         li_task   = &li->task[th];
1703
1704         if (li_m->nind + li_task->nind_r > li_m->ind_nalloc)
1705         {
1706             li_m->ind_nalloc = over_alloc_large(li_m->nind+li_task->nind_r);
1707             srenew(li_m->ind, li_m->ind_nalloc);
1708         }
1709
1710         for (b = 0; b < li_task->nind_r; b++)
1711         {
1712             li_m->ind[li_m->nind++] = li_task->ind_r[b];
1713         }
1714
1715         if (debug)
1716         {
1717             fprintf(debug, "LINCS thread %d: %d constraints\n",
1718                     th, li_task->nind);
1719         }
1720     }
1721
1722     if (debug)
1723     {
1724         fprintf(debug, "LINCS thread r: %d constraints\n",
1725                 li_m->nind);
1726     }
1727 }
1728
1729 /*! \brief There is no realloc with alignment, so here we make one for reals.
1730  * Note that this function does not preserve the contents of the memory.
1731  */
1732 static void resize_real_aligned(real **ptr, int nelem)
1733 {
1734     sfree_aligned(*ptr);
1735     snew_aligned(*ptr, nelem, align_bytes);
1736 }
1737
1738 //! Assign a constraint.
1739 static void assign_constraint(Lincs *li,
1740                               int constraint_index,
1741                               int a1, int a2,
1742                               real lenA, real lenB,
1743                               const t_blocka *at2con)
1744 {
1745     int con;
1746
1747     con = li->nc;
1748
1749     /* Make an mapping of local topology constraint index to LINCS index */
1750     li->con_index[constraint_index] = con;
1751
1752     li->bllen0[con]  = lenA;
1753     li->ddist[con]   = lenB - lenA;
1754     /* Set the length to the topology A length */
1755     li->bllen[con]   = lenA;
1756     li->bla[2*con]   = a1;
1757     li->bla[2*con+1] = a2;
1758
1759     /* Make space in the constraint connection matrix for constraints
1760      * connected to both end of the current constraint.
1761      */
1762     li->ncc +=
1763         at2con->index[a1 + 1] - at2con->index[a1] - 1 +
1764         at2con->index[a2 + 1] - at2con->index[a2] - 1;
1765
1766     li->blnr[con + 1] = li->ncc;
1767
1768     /* Increase the constraint count */
1769     li->nc++;
1770 }
1771
1772 /*! \brief Check if constraint with topology index constraint_index is connected
1773  * to other constraints, and if so add those connected constraints to our task. */
1774 static void check_assign_connected(Lincs *li,
1775                                    const t_iatom *iatom,
1776                                    const t_idef &idef,
1777                                    bool bDynamics,
1778                                    int a1, int a2,
1779                                    const t_blocka *at2con)
1780 {
1781     /* Currently this function only supports constraint groups
1782      * in which all constraints share at least one atom
1783      * (e.g. H-bond constraints).
1784      * Check both ends of the current constraint for
1785      * connected constraints. We need to assign those
1786      * to the same task.
1787      */
1788     int end;
1789
1790     for (end = 0; end < 2; end++)
1791     {
1792         int a, k;
1793
1794         a = (end == 0 ? a1 : a2);
1795
1796         for (k = at2con->index[a]; k < at2con->index[a + 1]; k++)
1797         {
1798             int cc;
1799
1800             cc = at2con->a[k];
1801             /* Check if constraint cc has not yet been assigned */
1802             if (li->con_index[cc] == -1)
1803             {
1804                 int  type;
1805                 real lenA, lenB;
1806
1807                 type = iatom[cc*3];
1808                 lenA = idef.iparams[type].constr.dA;
1809                 lenB = idef.iparams[type].constr.dB;
1810
1811                 if (bDynamics || lenA != 0 || lenB != 0)
1812                 {
1813                     assign_constraint(li, cc, iatom[3*cc + 1], iatom[3*cc + 2], lenA, lenB, at2con);
1814                 }
1815             }
1816         }
1817     }
1818 }
1819
1820 /*! \brief Check if constraint with topology index constraint_index is involved
1821  * in a constraint triangle, and if so add the other two constraints
1822  * in the triangle to our task. */
1823 static void check_assign_triangle(Lincs *li,
1824                                   const t_iatom *iatom,
1825                                   const t_idef &idef,
1826                                   bool bDynamics,
1827                                   int constraint_index,
1828                                   int a1, int a2,
1829                                   const t_blocka *at2con)
1830 {
1831     int nca, cc[32], ca[32], k;
1832     int c_triangle[2] = { -1, -1 };
1833
1834     nca = 0;
1835     for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1836     {
1837         int c;
1838
1839         c = at2con->a[k];
1840         if (c != constraint_index)
1841         {
1842             int aa1, aa2;
1843
1844             aa1 = iatom[c*3 + 1];
1845             aa2 = iatom[c*3 + 2];
1846             if (aa1 != a1)
1847             {
1848                 cc[nca] = c;
1849                 ca[nca] = aa1;
1850                 nca++;
1851             }
1852             if (aa2 != a1)
1853             {
1854                 cc[nca] = c;
1855                 ca[nca] = aa2;
1856                 nca++;
1857             }
1858         }
1859     }
1860
1861     for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1862     {
1863         int c;
1864
1865         c = at2con->a[k];
1866         if (c != constraint_index)
1867         {
1868             int aa1, aa2, i;
1869
1870             aa1 = iatom[c*3 + 1];
1871             aa2 = iatom[c*3 + 2];
1872             if (aa1 != a2)
1873             {
1874                 for (i = 0; i < nca; i++)
1875                 {
1876                     if (aa1 == ca[i])
1877                     {
1878                         c_triangle[0] = cc[i];
1879                         c_triangle[1] = c;
1880                     }
1881                 }
1882             }
1883             if (aa2 != a2)
1884             {
1885                 for (i = 0; i < nca; i++)
1886                 {
1887                     if (aa2 == ca[i])
1888                     {
1889                         c_triangle[0] = cc[i];
1890                         c_triangle[1] = c;
1891                     }
1892                 }
1893             }
1894         }
1895     }
1896
1897     if (c_triangle[0] >= 0)
1898     {
1899         int end;
1900
1901         for (end = 0; end < 2; end++)
1902         {
1903             /* Check if constraint c_triangle[end] has not yet been assigned */
1904             if (li->con_index[c_triangle[end]] == -1)
1905             {
1906                 int  i, type;
1907                 real lenA, lenB;
1908
1909                 i    = c_triangle[end]*3;
1910                 type = iatom[i];
1911                 lenA = idef.iparams[type].constr.dA;
1912                 lenB = idef.iparams[type].constr.dB;
1913
1914                 if (bDynamics || lenA != 0 || lenB != 0)
1915                 {
1916                     assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1917                 }
1918             }
1919         }
1920     }
1921 }
1922
1923 //! Sets matrix indices.
1924 static void set_matrix_indices(Lincs                *li,
1925                                const Task           *li_task,
1926                                const t_blocka       *at2con,
1927                                bool                  bSortMatrix)
1928 {
1929     int b;
1930
1931     for (b = li_task->b0; b < li_task->b1; b++)
1932     {
1933         int a1, a2, i, k;
1934
1935         a1 = li->bla[b*2];
1936         a2 = li->bla[b*2 + 1];
1937
1938         i = li->blnr[b];
1939         for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1940         {
1941             int concon;
1942
1943             concon = li->con_index[at2con->a[k]];
1944             if (concon != b)
1945             {
1946                 li->blbnb[i++] = concon;
1947             }
1948         }
1949         for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1950         {
1951             int concon;
1952
1953             concon = li->con_index[at2con->a[k]];
1954             if (concon != b)
1955             {
1956                 li->blbnb[i++] = concon;
1957             }
1958         }
1959
1960         if (bSortMatrix)
1961         {
1962             /* Order the blbnb matrix to optimize memory access */
1963             std::sort(&(li->blbnb[li->blnr[b]]), &(li->blbnb[li->blnr[b+1]]));
1964         }
1965     }
1966 }
1967
1968 void set_lincs(const t_idef         &idef,
1969                const t_mdatoms      &md,
1970                bool                  bDynamics,
1971                const t_commrec      *cr,
1972                Lincs                *li)
1973 {
1974     int          natoms;
1975     t_blocka     at2con;
1976     t_iatom     *iatom;
1977     int          i, ncc_alloc_old, ncon_tot;
1978
1979     li->nc_real = 0;
1980     li->nc      = 0;
1981     li->ncc     = 0;
1982     /* Zero the thread index ranges.
1983      * Otherwise without local constraints we could return with old ranges.
1984      */
1985     for (i = 0; i < li->ntask; i++)
1986     {
1987         li->task[i].b0   = 0;
1988         li->task[i].b1   = 0;
1989         li->task[i].nind = 0;
1990     }
1991     if (li->ntask > 1)
1992     {
1993         li->task[li->ntask].nind = 0;
1994     }
1995
1996     /* This is the local topology, so there are only F_CONSTR constraints */
1997     if (idef.il[F_CONSTR].nr == 0)
1998     {
1999         /* There are no constraints,
2000          * we do not need to fill any data structures.
2001          */
2002         return;
2003     }
2004
2005     if (debug)
2006     {
2007         fprintf(debug, "Building the LINCS connectivity\n");
2008     }
2009
2010     if (DOMAINDECOMP(cr))
2011     {
2012         if (cr->dd->constraints)
2013         {
2014             int start;
2015
2016             dd_get_constraint_range(cr->dd, &start, &natoms);
2017         }
2018         else
2019         {
2020             natoms = dd_numHomeAtoms(*cr->dd);
2021         }
2022     }
2023     else
2024     {
2025         natoms = md.homenr;
2026     }
2027
2028     at2con = make_at2con(natoms, idef.il, idef.iparams,
2029                          flexibleConstraintTreatment(bDynamics));
2030
2031     ncon_tot = idef.il[F_CONSTR].nr/3;
2032
2033     /* Ensure we have enough padding for aligned loads for each thread */
2034     if (ncon_tot + li->ntask*simd_width > li->nc_alloc || li->nc_alloc == 0)
2035     {
2036         li->nc_alloc = over_alloc_dd(ncon_tot + li->ntask*simd_width);
2037         srenew(li->con_index, li->nc_alloc);
2038         resize_real_aligned(&li->bllen0, li->nc_alloc);
2039         resize_real_aligned(&li->ddist, li->nc_alloc);
2040         srenew(li->bla, 2*li->nc_alloc);
2041         resize_real_aligned(&li->blc, li->nc_alloc);
2042         resize_real_aligned(&li->blc1, li->nc_alloc);
2043         srenew(li->blnr, li->nc_alloc + 1);
2044         resize_real_aligned(&li->bllen, li->nc_alloc);
2045         srenew(li->tmpv, li->nc_alloc);
2046         if (DOMAINDECOMP(cr))
2047         {
2048             srenew(li->nlocat, li->nc_alloc);
2049         }
2050         resize_real_aligned(&li->tmp1, li->nc_alloc);
2051         resize_real_aligned(&li->tmp2, li->nc_alloc);
2052         resize_real_aligned(&li->tmp3, li->nc_alloc);
2053         resize_real_aligned(&li->tmp4, li->nc_alloc);
2054         resize_real_aligned(&li->mlambda, li->nc_alloc);
2055     }
2056
2057     iatom = idef.il[F_CONSTR].iatoms;
2058
2059     ncc_alloc_old = li->ncc_alloc;
2060     li->blnr[0]   = li->ncc;
2061
2062     /* Assign the constraints for li->ntask LINCS tasks.
2063      * We target a uniform distribution of constraints over the tasks.
2064      * Note that when flexible constraints are present, but are removed here
2065      * (e.g. because we are doing EM) we get imbalance, but since that doesn't
2066      * happen during normal MD, that's ok.
2067      */
2068     int ncon_assign, ncon_target, con, th;
2069
2070     /* Determine the number of constraints we need to assign here */
2071     ncon_assign      = ncon_tot;
2072     if (!bDynamics)
2073     {
2074         /* With energy minimization, flexible constraints are ignored
2075          * (and thus minimized, as they should be).
2076          */
2077         ncon_assign -= countFlexibleConstraints(idef.il, idef.iparams);
2078     }
2079
2080     /* Set the target constraint count per task to exactly uniform,
2081      * this might be overridden below.
2082      */
2083     ncon_target = (ncon_assign + li->ntask - 1)/li->ntask;
2084
2085     /* Mark all constraints as unassigned by setting their index to -1 */
2086     for (con = 0; con < ncon_tot; con++)
2087     {
2088         li->con_index[con] = -1;
2089     }
2090
2091     con = 0;
2092     for (th = 0; th < li->ntask; th++)
2093     {
2094         Task *li_task;
2095
2096         li_task = &li->task[th];
2097
2098 #if GMX_SIMD_HAVE_REAL
2099         /* With indepedent tasks we likely have H-bond constraints or constraint
2100          * pairs. The connected constraints will be pulled into the task, so the
2101          * constraints per task will often exceed ncon_target.
2102          * Triangle constraints can also increase the count, but there are
2103          * relatively few of those, so we usually expect to get ncon_target.
2104          */
2105         if (li->bTaskDep)
2106         {
2107             /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
2108              * since otherwise a lot of operations can be wasted.
2109              * There are several ways to round here, we choose the one
2110              * that alternates block sizes, which helps with Intel HT.
2111              */
2112             ncon_target = ((ncon_assign*(th + 1))/li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1) & ~(GMX_SIMD_REAL_WIDTH - 1);
2113         }
2114 #endif      // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
2115
2116         /* Continue filling the arrays where we left off with the previous task,
2117          * including padding for SIMD.
2118          */
2119         li_task->b0 = li->nc;
2120
2121         while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2122         {
2123             if (li->con_index[con] == -1)
2124             {
2125                 int  type, a1, a2;
2126                 real lenA, lenB;
2127
2128                 type   = iatom[3*con];
2129                 a1     = iatom[3*con + 1];
2130                 a2     = iatom[3*con + 2];
2131                 lenA   = idef.iparams[type].constr.dA;
2132                 lenB   = idef.iparams[type].constr.dB;
2133                 /* Skip the flexible constraints when not doing dynamics */
2134                 if (bDynamics || lenA != 0 || lenB != 0)
2135                 {
2136                     assign_constraint(li, con, a1, a2, lenA, lenB, &at2con);
2137
2138                     if (li->ntask > 1 && !li->bTaskDep)
2139                     {
2140                         /* We can generate independent tasks. Check if we
2141                          * need to assign connected constraints to our task.
2142                          */
2143                         check_assign_connected(li, iatom, idef, bDynamics,
2144                                                a1, a2, &at2con);
2145                     }
2146                     if (li->ntask > 1 && li->ncg_triangle > 0)
2147                     {
2148                         /* Ensure constraints in one triangle are assigned
2149                          * to the same task.
2150                          */
2151                         check_assign_triangle(li, iatom, idef, bDynamics,
2152                                               con, a1, a2, &at2con);
2153                     }
2154                 }
2155             }
2156
2157             con++;
2158         }
2159
2160         li_task->b1 = li->nc;
2161
2162         if (simd_width > 1)
2163         {
2164             /* Copy the last atom pair indices and lengths for constraints
2165              * up to a multiple of simd_width, such that we can do all
2166              * SIMD operations without having to worry about end effects.
2167              */
2168             int i, last;
2169
2170             li->nc = ((li_task->b1 + simd_width - 1)/simd_width)*simd_width;
2171             last   = li_task->b1 - 1;
2172             for (i = li_task->b1; i < li->nc; i++)
2173             {
2174                 li->bla[i*2    ] = li->bla[last*2    ];
2175                 li->bla[i*2 + 1] = li->bla[last*2 + 1];
2176                 li->bllen0[i]    = li->bllen0[last];
2177                 li->ddist[i]     = li->ddist[last];
2178                 li->bllen[i]     = li->bllen[last];
2179                 li->blnr[i + 1]  = li->blnr[last + 1];
2180             }
2181         }
2182
2183         /* Keep track of how many constraints we assigned */
2184         li->nc_real += li_task->b1 - li_task->b0;
2185
2186         if (debug)
2187         {
2188             fprintf(debug, "LINCS task %d constraints %d - %d\n",
2189                     th, li_task->b0, li_task->b1);
2190         }
2191     }
2192
2193     assert(li->nc_real == ncon_assign);
2194
2195     bool bSortMatrix;
2196
2197     /* Without DD we order the blbnb matrix to optimize memory access.
2198      * With DD the overhead of sorting is more than the gain during access.
2199      */
2200     bSortMatrix = !DOMAINDECOMP(cr);
2201
2202     if (li->ncc > li->ncc_alloc)
2203     {
2204         li->ncc_alloc = over_alloc_small(li->ncc);
2205         srenew(li->blbnb, li->ncc_alloc);
2206     }
2207
2208 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2209     for (th = 0; th < li->ntask; th++)
2210     {
2211         try
2212         {
2213             Task *li_task;
2214
2215             li_task = &li->task[th];
2216
2217             if (li->ncg_triangle > 0 &&
2218                 li_task->b1 - li_task->b0 > li_task->tri_alloc)
2219             {
2220                 /* This is allocating too much, but it is difficult to improve */
2221                 li_task->tri_alloc = over_alloc_dd(li_task->b1 - li_task->b0);
2222                 srenew(li_task->triangle, li_task->tri_alloc);
2223                 srenew(li_task->tri_bits, li_task->tri_alloc);
2224             }
2225
2226             set_matrix_indices(li, li_task, &at2con, bSortMatrix);
2227         }
2228         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2229     }
2230
2231     done_blocka(&at2con);
2232
2233     if (cr->dd == nullptr)
2234     {
2235         /* Since the matrix is static, we should free some memory */
2236         li->ncc_alloc = li->ncc;
2237         srenew(li->blbnb, li->ncc_alloc);
2238     }
2239
2240     if (li->ncc_alloc > ncc_alloc_old)
2241     {
2242         srenew(li->blmf, li->ncc_alloc);
2243         srenew(li->blmf1, li->ncc_alloc);
2244         srenew(li->tmpncc, li->ncc_alloc);
2245     }
2246
2247     gmx::ArrayRef<const int> nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2248     if (!nlocat_dd.empty())
2249     {
2250         /* Convert nlocat from local topology to LINCS constraint indexing */
2251         for (con = 0; con < ncon_tot; con++)
2252         {
2253             li->nlocat[li->con_index[con]] = nlocat_dd[con];
2254         }
2255     }
2256     else
2257     {
2258         li->nlocat = nullptr;
2259     }
2260
2261     if (debug)
2262     {
2263         fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n",
2264                 li->nc_real, li->nc, li->ncc);
2265     }
2266
2267     if (li->ntask > 1)
2268     {
2269         lincs_thread_setup(li, md.nr);
2270     }
2271
2272     set_lincs_matrix(li, md.invmass, md.lambda);
2273 }
2274
2275 //! Issues a warning when LINCS constraints cannot be satisfied.
2276 static void lincs_warning(gmx_domdec_t *dd, const rvec *x, rvec *xprime, t_pbc *pbc,
2277                           int ncons, const int *bla, real *bllen, real wangle,
2278                           int maxwarn, int *warncount)
2279 {
2280     int  b, i, j;
2281     rvec v0, v1;
2282     real wfac, d0, d1, cosine;
2283
2284     wfac = std::cos(DEG2RAD*wangle);
2285
2286     fprintf(stderr,
2287             "bonds that rotated more than %g degrees:\n"
2288             " atom 1 atom 2  angle  previous, current, constraint length\n",
2289             wangle);
2290
2291     for (b = 0; b < ncons; b++)
2292     {
2293         i = bla[2*b];
2294         j = bla[2*b+1];
2295         if (pbc)
2296         {
2297             pbc_dx_aiuc(pbc, x[i], x[j], v0);
2298             pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2299         }
2300         else
2301         {
2302             rvec_sub(x[i], x[j], v0);
2303             rvec_sub(xprime[i], xprime[j], v1);
2304         }
2305         d0     = norm(v0);
2306         d1     = norm(v1);
2307         cosine = ::iprod(v0, v1)/(d0*d1);
2308         if (cosine < wfac)
2309         {
2310             fprintf(stderr,
2311                     " %6d %6d  %5.1f  %8.4f %8.4f    %8.4f\n",
2312                     ddglatnr(dd, i), ddglatnr(dd, j),
2313                     RAD2DEG*std::acos(cosine), d0, d1, bllen[b]);
2314             if (!std::isfinite(d1))
2315             {
2316                 gmx_fatal(FARGS, "Bond length not finite.");
2317             }
2318
2319             (*warncount)++;
2320         }
2321     }
2322     if (*warncount > maxwarn)
2323     {
2324         too_many_constraint_warnings(econtLINCS, *warncount);
2325     }
2326 }
2327
2328 //! Determine how well the constraints have been satisfied.
2329 static void cconerr(const Lincs *lincsd,
2330                     rvec *x, t_pbc *pbc,
2331                     real *ncons_loc, real *ssd, real *max, int *imax)
2332 {
2333     const int  *bla, *nlocat;
2334     const real *bllen;
2335     real        ma, ssd2;
2336     int         count, im, task;
2337
2338     bla    = lincsd->bla;
2339     bllen  = lincsd->bllen;
2340     nlocat = lincsd->nlocat;
2341
2342     ma    = 0;
2343     ssd2  = 0;
2344     im    = 0;
2345     count = 0;
2346     for (task = 0; task < lincsd->ntask; task++)
2347     {
2348         int b;
2349
2350         for (b = lincsd->task[task].b0; b < lincsd->task[task].b1; b++)
2351         {
2352             real len, d, r2;
2353             rvec dx;
2354
2355             if (pbc)
2356             {
2357                 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
2358             }
2359             else
2360             {
2361                 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
2362             }
2363             r2  = ::norm2(dx);
2364             len = r2*gmx::invsqrt(r2);
2365             d   = std::abs(len/bllen[b]-1);
2366             if (d > ma && (nlocat == nullptr || nlocat[b]))
2367             {
2368                 ma = d;
2369                 im = b;
2370             }
2371             if (nlocat == nullptr)
2372             {
2373                 ssd2 += d*d;
2374                 count++;
2375             }
2376             else
2377             {
2378                 ssd2  += nlocat[b]*d*d;
2379                 count += nlocat[b];
2380             }
2381         }
2382     }
2383
2384     *ncons_loc = (nlocat ? 0.5 : 1)*count;
2385     *ssd       = (nlocat ? 0.5 : 1)*ssd2;
2386     *max       = ma;
2387     *imax      = im;
2388 }
2389
2390 bool constrain_lincs(bool computeRmsd,
2391                      const t_inputrec &ir,
2392                      int64_t step,
2393                      Lincs *lincsd, const t_mdatoms &md,
2394                      const t_commrec *cr,
2395                      const gmx_multisim_t &ms,
2396                      const rvec *x, rvec *xprime, rvec *min_proj,
2397                      matrix box, t_pbc *pbc,
2398                      real lambda, real *dvdlambda,
2399                      real invdt, rvec *v,
2400                      bool bCalcVir, tensor vir_r_m_dr,
2401                      ConstraintVariable econq,
2402                      t_nrnb *nrnb,
2403                      int maxwarn, int *warncount)
2404 {
2405     gmx_bool  bCalcDHDL;
2406     char      buf2[22], buf3[STRLEN];
2407     int       i, p_imax;
2408     real      ncons_loc, p_ssd, p_max = 0;
2409     rvec      dx;
2410     bool      bOK, bWarn;
2411
2412     bOK = TRUE;
2413
2414     /* This boolean should be set by a flag passed to this routine.
2415      * We can also easily check if any constraint length is changed,
2416      * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2417      */
2418     bCalcDHDL = (ir.efep != efepNO && dvdlambda != nullptr);
2419
2420     if (lincsd->nc == 0 && cr->dd == nullptr)
2421     {
2422         if (computeRmsd)
2423         {
2424             lincsd->rmsdData = {{0}};
2425         }
2426
2427         return bOK;
2428     }
2429
2430     if (econq == ConstraintVariable::Positions)
2431     {
2432         /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2433          * also with efep!=fepNO.
2434          */
2435         if (ir.efep != efepNO)
2436         {
2437             if (md.nMassPerturbed && lincsd->matlam != md.lambda)
2438             {
2439                 set_lincs_matrix(lincsd, md.invmass, md.lambda);
2440             }
2441
2442             for (i = 0; i < lincsd->nc; i++)
2443             {
2444                 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
2445             }
2446         }
2447
2448         if (lincsd->ncg_flex)
2449         {
2450             /* Set the flexible constraint lengths to the old lengths */
2451             if (pbc != nullptr)
2452             {
2453                 for (i = 0; i < lincsd->nc; i++)
2454                 {
2455                     if (lincsd->bllen[i] == 0)
2456                     {
2457                         pbc_dx_aiuc(pbc, x[lincsd->bla[2*i]], x[lincsd->bla[2*i+1]], dx);
2458                         lincsd->bllen[i] = norm(dx);
2459                     }
2460                 }
2461             }
2462             else
2463             {
2464                 for (i = 0; i < lincsd->nc; i++)
2465                 {
2466                     if (lincsd->bllen[i] == 0)
2467                     {
2468                         lincsd->bllen[i] =
2469                             std::sqrt(distance2(x[lincsd->bla[2*i]],
2470                                                 x[lincsd->bla[2*i+1]]));
2471                     }
2472                 }
2473             }
2474         }
2475
2476         if (debug)
2477         {
2478             cconerr(lincsd, xprime, pbc,
2479                     &ncons_loc, &p_ssd, &p_max, &p_imax);
2480         }
2481
2482         /* This bWarn var can be updated by multiple threads
2483          * at the same time. But as we only need to detect
2484          * if a warning occurred or not, this is not an issue.
2485          */
2486         bWarn = FALSE;
2487
2488         /* The OpenMP parallel region of constrain_lincs for coords */
2489 #pragma omp parallel num_threads(lincsd->ntask)
2490         {
2491             try
2492             {
2493                 int th = gmx_omp_get_thread_num();
2494
2495                 clear_mat(lincsd->task[th].vir_r_m_dr);
2496
2497                 do_lincs(x, xprime, box, pbc, lincsd, th,
2498                          md.invmass, cr,
2499                          bCalcDHDL,
2500                          ir.LincsWarnAngle, &bWarn,
2501                          invdt, v, bCalcVir,
2502                          th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2503             }
2504             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2505         }
2506
2507         if (debug && lincsd->nc > 0)
2508         {
2509             fprintf(debug, "   Rel. Constraint Deviation:  RMS         MAX     between atoms\n");
2510             fprintf(debug, "       Before LINCS          %.6f    %.6f %6d %6d\n",
2511                     std::sqrt(p_ssd/ncons_loc), p_max,
2512                     ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2513                     ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2514         }
2515         if (computeRmsd || debug)
2516         {
2517             cconerr(lincsd, xprime, pbc,
2518                     &ncons_loc, &p_ssd, &p_max, &p_imax);
2519             lincsd->rmsdData[0] = ncons_loc;
2520             lincsd->rmsdData[1] = p_ssd;
2521         }
2522         else
2523         {
2524             lincsd->rmsdData = {{0}};
2525         }
2526         if (debug && lincsd->nc > 0)
2527         {
2528             fprintf(debug,
2529                     "        After LINCS          %.6f    %.6f %6d %6d\n\n",
2530                     std::sqrt(p_ssd/ncons_loc), p_max,
2531                     ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2532                     ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2533         }
2534
2535         if (bWarn)
2536         {
2537             if (maxwarn < INT_MAX)
2538             {
2539                 cconerr(lincsd, xprime, pbc,
2540                         &ncons_loc, &p_ssd, &p_max, &p_imax);
2541                 if (isMultiSim(&ms))
2542                 {
2543                     sprintf(buf3, " in simulation %d", ms.sim);
2544                 }
2545                 else
2546                 {
2547                     buf3[0] = 0;
2548                 }
2549                 fprintf(stderr,
2550                         "\nStep %s, time %g (ps)  LINCS WARNING%s\n"
2551                         "relative constraint deviation after LINCS:\n"
2552                         "rms %.6f, max %.6f (between atoms %d and %d)\n",
2553                         gmx_step_str(step, buf2), ir.init_t+step*ir.delta_t,
2554                         buf3,
2555                         std::sqrt(p_ssd/ncons_loc), p_max,
2556                         ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2557                         ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2558
2559                 lincs_warning(cr->dd, x, xprime, pbc,
2560                               lincsd->nc, lincsd->bla, lincsd->bllen,
2561                               ir.LincsWarnAngle, maxwarn, warncount);
2562             }
2563             bOK = (p_max < 0.5);
2564         }
2565
2566         if (lincsd->ncg_flex)
2567         {
2568             for (i = 0; (i < lincsd->nc); i++)
2569             {
2570                 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2571                 {
2572                     lincsd->bllen[i] = 0;
2573                 }
2574             }
2575         }
2576     }
2577     else
2578     {
2579         /* The OpenMP parallel region of constrain_lincs for derivatives */
2580 #pragma omp parallel num_threads(lincsd->ntask)
2581         {
2582             try
2583             {
2584                 int th = gmx_omp_get_thread_num();
2585
2586                 do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
2587                           md.invmass, econq, bCalcDHDL,
2588                           bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2589             }
2590             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2591         }
2592     }
2593
2594     if (bCalcDHDL)
2595     {
2596         /* Reduce the dH/dlambda contributions over the threads */
2597         real dhdlambda;
2598         int  th;
2599
2600         dhdlambda = 0;
2601         for (th = 0; th < lincsd->ntask; th++)
2602         {
2603             dhdlambda += lincsd->task[th].dhdlambda;
2604         }
2605         if (econq == ConstraintVariable::Positions)
2606         {
2607             /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2608             /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2609             dhdlambda /= ir.delta_t*ir.delta_t;
2610         }
2611         *dvdlambda += dhdlambda;
2612     }
2613
2614     if (bCalcVir && lincsd->ntask > 1)
2615     {
2616         for (i = 1; i < lincsd->ntask; i++)
2617         {
2618             m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2619         }
2620     }
2621
2622     /* count assuming nit=1 */
2623     inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2624     inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
2625     if (lincsd->ntriangle > 0)
2626     {
2627         inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
2628     }
2629     if (v)
2630     {
2631         inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real*2);
2632     }
2633     if (bCalcVir)
2634     {
2635         inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);
2636     }
2637
2638     return bOK;
2639 }
2640
2641 }  // namespace gmx