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