2 * This file is part of the GROMACS molecular simulation package.
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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Defines LINCS code.
41 * \author Berk Hess <hess@kth.se>
42 * \author Mark Abraham <mark.j.abraham@gmail.com>
43 * \ingroup module_mdlib
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/paddedvector.h"
63 #include "gromacs/math/units.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/mdlib/constr.h"
66 #include "gromacs/mdlib/gmx_omp_nthreads.h"
67 #include "gromacs/mdrunutility/multisim.h"
68 #include "gromacs/mdtypes/commrec.h"
69 #include "gromacs/mdtypes/inputrec.h"
70 #include "gromacs/mdtypes/md_enums.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/mtop_util.h"
77 #include "gromacs/utility/alignedallocator.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/listoflists.h"
86 #include "gromacs/utility/pleasecite.h"
91 //! Indices of the two atoms involved in a single constraint
94 //! \brief Constructor, does not initialize to catch bugs and faster construction
103 //! Unit of work within LINCS.
106 //! First constraint for this task.
108 //! b1-1 is the last constraint for this task.
110 //! The number of constraints in triangles.
112 //! The list of triangle constraints.
113 std::vector<int> triangle;
114 //! The bits tell if the matrix element should be used.
115 std::vector<int> tri_bits;
116 //! Constraint index for updating atom data.
117 std::vector<int> ind;
118 //! Constraint index for updating atom data.
119 std::vector<int> ind_r;
120 //! Temporary variable for virial calculation.
121 tensor vir_r_m_dr = { { 0 } };
122 //! Temporary variable for lambda derivative.
126 /*! \brief Data for LINCS algorithm.
131 //! The global number of constraints.
133 //! The global number of flexible constraints.
135 //! The global number of constraints in triangles.
136 int ncg_triangle = 0;
137 //! The number of iterations.
139 //! The order of the matrix expansion.
141 //! The maximum number of constraints connected to a single atom.
144 //! The number of real constraints.
146 //! The number of constraints including padding for SIMD.
148 //! The number of constraint connections.
150 //! The FE lambda value used for filling blc and blmf.
152 //! mapping from topology to LINCS constraints.
153 std::vector<int> con_index;
154 //! The reference distance in topology A.
155 std::vector<real, AlignedAllocator<real>> bllen0;
156 //! The reference distance in top B - the r.d. in top A.
157 std::vector<real, AlignedAllocator<real>> ddist;
158 //! The atom pairs involved in the constraints.
159 std::vector<AtomPair> atoms;
160 //! 1/sqrt(invmass1 invmass2).
161 std::vector<real, AlignedAllocator<real>> blc;
162 //! As blc, but with all masses 1.
163 std::vector<real, AlignedAllocator<real>> blc1;
164 //! Index into blbnb and blmf.
165 std::vector<int> blnr;
166 //! List of constraint connections.
167 std::vector<int> blbnb;
168 //! The local number of constraints in triangles.
170 //! The number of constraint connections in triangles.
171 int ncc_triangle = 0;
172 //! Communicate before each LINCS interation.
173 bool bCommIter = false;
174 //! Matrix of mass factors for constraint connections.
175 std::vector<real> blmf;
176 //! As blmf, but with all masses 1.
177 std::vector<real> blmf1;
178 //! The reference bond length.
179 std::vector<real, AlignedAllocator<real>> bllen;
180 //! The local atom count per constraint, can be NULL.
181 std::vector<int> nlocat;
183 /*! \brief The number of tasks used for LINCS work.
185 * \todo This is mostly used to loop over \c task, which would
186 * be nicer to do with range-based for loops, but the thread
187 * index is used for constructing bit masks and organizing the
188 * virial output buffer, so other things need to change,
191 /*! \brief LINCS thread division */
192 std::vector<Task> task;
193 //! Atom flags for thread parallelization.
194 std::vector<gmx_bitmask_t> atf;
195 //! Are the LINCS tasks interdependent?
196 bool bTaskDep = false;
197 //! Are there triangle constraints that cross task borders?
198 bool bTaskDepTri = false;
199 //! Arrays for temporary storage in the LINCS algorithm.
201 PaddedVector<gmx::RVec> tmpv;
202 std::vector<real> tmpncc;
203 std::vector<real, AlignedAllocator<real>> tmp1;
204 std::vector<real, AlignedAllocator<real>> tmp2;
205 std::vector<real, AlignedAllocator<real>> tmp3;
206 std::vector<real, AlignedAllocator<real>> tmp4;
208 //! The Lagrange multipliers times -1.
209 std::vector<real, AlignedAllocator<real>> mlambda;
210 //! Storage for the constraint RMS relative deviation output.
211 std::array<real, 2> rmsdData = { { 0 } };
214 /*! \brief Define simd_width for memory allocation used for SIMD code */
215 #if GMX_SIMD_HAVE_REAL
216 static const int simd_width = GMX_SIMD_REAL_WIDTH;
218 static const int simd_width = 1;
221 ArrayRef<real> lincs_rmsdData(Lincs* lincsd)
223 return lincsd->rmsdData;
226 real lincs_rmsd(const Lincs* lincsd)
228 if (lincsd->rmsdData[0] > 0)
230 return std::sqrt(lincsd->rmsdData[1] / lincsd->rmsdData[0]);
238 /*! \brief Do a set of nrec LINCS matrix multiplications.
240 * This function will return with up to date thread-local
241 * constraint data, without an OpenMP barrier.
243 static void lincs_matrix_expand(const Lincs& lincsd,
245 gmx::ArrayRef<const real> blcc,
246 gmx::ArrayRef<real> rhs1,
247 gmx::ArrayRef<real> rhs2,
248 gmx::ArrayRef<real> sol)
250 gmx::ArrayRef<const int> blnr = lincsd.blnr;
251 gmx::ArrayRef<const int> blbnb = lincsd.blbnb;
253 const int b0 = li_task.b0;
254 const int b1 = li_task.b1;
255 const int nrec = lincsd.nOrder;
257 for (int rec = 0; rec < nrec; rec++)
263 for (int b = b0; b < b1; b++)
269 for (n = blnr[b]; n < blnr[b + 1]; n++)
271 mvb = mvb + blcc[n] * rhs1[blbnb[n]];
274 sol[b] = sol[b] + mvb;
277 std::swap(rhs1, rhs2);
278 } /* nrec*(ncons+2*nrtot) flops */
280 if (lincsd.ntriangle > 0)
282 /* Perform an extra nrec recursions for only the constraints
283 * involved in rigid triangles.
284 * In this way their accuracy should come close to those of the other
285 * constraints, since traingles of constraints can produce eigenvalues
286 * around 0.7, while the effective eigenvalue for bond constraints
287 * is around 0.4 (and 0.7*0.7=0.5).
292 /* We need a barrier here, since other threads might still be
293 * reading the contents of rhs1 and/o rhs2.
294 * We could avoid this barrier by introducing two extra rhs
295 * arrays for the triangle constraints only.
300 /* Constraints involved in a triangle are ensured to be in the same
301 * LINCS task. This means no barriers are required during the extra
302 * iterations for the triangle constraints.
304 gmx::ArrayRef<const int> triangle = li_task.triangle;
305 gmx::ArrayRef<const int> tri_bits = li_task.tri_bits;
307 for (int rec = 0; rec < nrec; rec++)
309 for (int tb = 0; tb < li_task.ntriangle; tb++)
311 int b, bits, nr0, nr1, n;
319 for (n = nr0; n < nr1; n++)
321 if (bits & (1 << (n - nr0)))
323 mvb = mvb + blcc[n] * rhs1[blbnb[n]];
327 sol[b] = sol[b] + mvb;
330 std::swap(rhs1, rhs2);
331 } /* nrec*(ntriangle + ncc_triangle*2) flops */
333 if (lincsd.bTaskDepTri)
335 /* The constraints triangles are decoupled from each other,
336 * but constraints in one triangle cross thread task borders.
337 * We could probably avoid this with more advanced setup code.
344 //! Update atomic coordinates when an index is not required.
345 static void lincs_update_atoms_noind(int ncons,
346 gmx::ArrayRef<const AtomPair> atoms,
348 gmx::ArrayRef<const real> fac,
349 gmx::ArrayRef<const gmx::RVec> r,
350 gmx::ArrayRef<const real> invmass,
353 if (!invmass.empty())
355 for (int b = 0; b < ncons; b++)
357 int i = atoms[b].index1;
358 int j = atoms[b].index2;
359 real mvb = preFactor * fac[b];
360 real im1 = invmass[i];
361 real im2 = invmass[j];
362 real tmp0 = r[b][0] * mvb;
363 real tmp1 = r[b][1] * mvb;
364 real tmp2 = r[b][2] * mvb;
365 x[i][0] -= tmp0 * im1;
366 x[i][1] -= tmp1 * im1;
367 x[i][2] -= tmp2 * im1;
368 x[j][0] += tmp0 * im2;
369 x[j][1] += tmp1 * im2;
370 x[j][2] += tmp2 * im2;
371 } /* 16 ncons flops */
375 for (int b = 0; b < ncons; b++)
377 int i = atoms[b].index1;
378 int j = atoms[b].index2;
379 real mvb = preFactor * fac[b];
380 real tmp0 = r[b][0] * mvb;
381 real tmp1 = r[b][1] * mvb;
382 real tmp2 = r[b][2] * mvb;
393 //! Update atomic coordinates when an index is required.
394 static void lincs_update_atoms_ind(gmx::ArrayRef<const int> ind,
395 gmx::ArrayRef<const AtomPair> atoms,
397 gmx::ArrayRef<const real> fac,
398 gmx::ArrayRef<const gmx::RVec> r,
399 gmx::ArrayRef<const real> invmass,
402 if (!invmass.empty())
406 int i = atoms[b].index1;
407 int j = atoms[b].index2;
408 real mvb = preFactor * fac[b];
409 real im1 = invmass[i];
410 real im2 = invmass[j];
411 real tmp0 = r[b][0] * mvb;
412 real tmp1 = r[b][1] * mvb;
413 real tmp2 = r[b][2] * mvb;
414 x[i][0] -= tmp0 * im1;
415 x[i][1] -= tmp1 * im1;
416 x[i][2] -= tmp2 * im1;
417 x[j][0] += tmp0 * im2;
418 x[j][1] += tmp1 * im2;
419 x[j][2] += tmp2 * im2;
420 } /* 16 ncons flops */
426 int i = atoms[b].index1;
427 int j = atoms[b].index2;
428 real mvb = preFactor * fac[b];
429 real tmp0 = r[b][0] * mvb;
430 real tmp1 = r[b][1] * mvb;
431 real tmp2 = r[b][2] * mvb;
438 } /* 16 ncons flops */
442 //! Update coordinates for atoms.
443 static void lincs_update_atoms(Lincs* li,
446 gmx::ArrayRef<const real> fac,
447 gmx::ArrayRef<const gmx::RVec> r,
448 gmx::ArrayRef<const real> invmass,
453 /* Single thread, we simply update for all constraints */
454 lincs_update_atoms_noind(li->nc_real, li->atoms, preFactor, fac, r, invmass, x);
458 /* Update the atom vector components for our thread local
459 * constraints that only access our local atom range.
460 * This can be done without a barrier.
462 lincs_update_atoms_ind(li->task[th].ind, li->atoms, preFactor, fac, r, invmass, x);
464 if (!li->task[li->ntask].ind.empty())
466 /* Update the constraints that operate on atoms
467 * in multiple thread atom blocks on the master thread.
472 lincs_update_atoms_ind(li->task[li->ntask].ind, li->atoms, preFactor, fac, r, invmass, x);
478 #if GMX_SIMD_HAVE_REAL
479 //! Helper function so that we can run TSAN with SIMD support (where implemented).
481 static inline void gmx_simdcall gatherLoadUTransposeTSANSafe(const real* base,
482 const std::int32_t* offset,
487 # if (CMAKE_BUILD_TYPE == CMAKE_BUILD_TYPE_TSAN) && GMX_SIMD_X86_AVX2_256
488 // This function is only implemented in this case
489 gatherLoadUTransposeSafe<align>(base, offset, v0, v1, v2);
491 gatherLoadUTranspose<align>(base, offset, v0, v1, v2);
495 /*! \brief Calculate the constraint distance vectors r to project on from x.
497 * Determine the right-hand side of the matrix equation using quantity f.
498 * This function only differs from calc_dr_x_xp_simd below in that
499 * no constraint length is subtracted and no PBC is used for f. */
500 static void gmx_simdcall calc_dr_x_f_simd(int b0,
502 gmx::ArrayRef<const AtomPair> atoms,
503 const rvec* gmx_restrict x,
504 const rvec* gmx_restrict f,
505 const real* gmx_restrict blc,
506 const real* pbc_simd,
507 rvec* gmx_restrict r,
508 real* gmx_restrict rhs,
509 real* gmx_restrict sol)
511 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
513 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
515 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
520 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
522 SimdReal x0_S, y0_S, z0_S;
523 SimdReal x1_S, y1_S, z1_S;
524 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
525 SimdReal fx_S, fy_S, fz_S, ip_S, rhs_S;
526 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
527 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
529 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
531 offset0[i] = atoms[bs + i].index1;
532 offset1[i] = atoms[bs + i].index2;
535 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset0, &x0_S, &y0_S, &z0_S);
536 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset1, &x1_S, &y1_S, &z1_S);
541 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
543 n2_S = norm2(rx_S, ry_S, rz_S);
544 il_S = invsqrt(n2_S);
550 transposeScatterStoreU<3>(reinterpret_cast<real*>(r + bs), offset2, rx_S, ry_S, rz_S);
552 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(f), offset0, &x0_S, &y0_S, &z0_S);
553 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(f), offset1, &x1_S, &y1_S, &z1_S);
558 ip_S = iprod(rx_S, ry_S, rz_S, fx_S, fy_S, fz_S);
560 rhs_S = load<SimdReal>(blc + bs) * ip_S;
562 store(rhs + bs, rhs_S);
563 store(sol + bs, rhs_S);
566 #endif // GMX_SIMD_HAVE_REAL
568 /*! \brief LINCS projection, works on derivatives of the coordinates. */
569 static void do_lincsp(ArrayRefWithPadding<const RVec> xPadded,
570 ArrayRefWithPadding<RVec> fPadded,
575 ArrayRef<const real> invmass,
576 ConstraintVariable econq,
581 const int b0 = lincsd->task[th].b0;
582 const int b1 = lincsd->task[th].b1;
584 gmx::ArrayRef<const AtomPair> atoms = lincsd->atoms;
585 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
586 gmx::ArrayRef<const int> blnr = lincsd->blnr;
587 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
589 gmx::ArrayRef<const real> blc;
590 gmx::ArrayRef<const real> blmf;
591 if (econq != ConstraintVariable::Force)
593 /* Use mass-weighted parameters */
599 /* Use non mass-weighted parameters */
601 blmf = lincsd->blmf1;
603 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
604 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
605 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
606 gmx::ArrayRef<real> sol = lincsd->tmp3;
608 const rvec* x = as_rvec_array(xPadded.paddedArrayRef().data());
609 rvec* f = as_rvec_array(fPadded.paddedArrayRef().data());
611 #if GMX_SIMD_HAVE_REAL
612 /* This SIMD code does the same as the plain-C code after the #else.
613 * The only difference is that we always call pbc code, as with SIMD
614 * the overhead of pbc computation (when not needed) is small.
616 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9 * GMX_SIMD_REAL_WIDTH];
618 /* Convert the pbc struct for SIMD */
619 set_pbc_simd(pbc, pbc_simd);
621 /* Compute normalized x i-j vectors, store in r.
622 * Compute the inner product of r and xp i-j and store in rhs1.
625 b0, b1, atoms, x, f, blc.data(), pbc_simd, as_rvec_array(r.data()), rhs1.data(), sol.data());
627 #else // GMX_SIMD_HAVE_REAL
629 /* Compute normalized i-j vectors */
632 for (int b = b0; b < b1; b++)
636 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
642 for (int b = b0; b < b1; b++)
646 rvec_sub(x[atoms[b].index1], x[atoms[b].index2], dx);
648 } /* 16 ncons flops */
651 for (int b = b0; b < b1; b++)
653 int i = atoms[b].index1;
654 int j = atoms[b].index2;
656 * (r[b][0] * (f[i][0] - f[j][0]) + r[b][1] * (f[i][1] - f[j][1])
657 + r[b][2] * (f[i][2] - f[j][2]));
663 #endif // GMX_SIMD_HAVE_REAL
665 if (lincsd->bTaskDep)
667 /* We need a barrier, since the matrix construction below
668 * can access entries in r of other threads.
673 /* Construct the (sparse) LINCS matrix */
674 for (int b = b0; b < b1; b++)
676 for (int n = blnr[b]; n < blnr[b + 1]; n++)
678 blcc[n] = blmf[n] * ::iprod(r[b], r[blbnb[n]]);
681 /* Together: 23*ncons + 6*nrtot flops */
683 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
684 /* nrec*(ncons+2*nrtot) flops */
686 if (econq == ConstraintVariable::Deriv_FlexCon)
688 /* We only want to constraint the flexible constraints,
689 * so we mask out the normal ones by setting sol to 0.
691 for (int b = b0; b < b1; b++)
693 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
700 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
701 for (int b = b0; b < b1; b++)
706 /* When constraining forces, we should not use mass weighting,
707 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
709 lincs_update_atoms(lincsd,
714 (econq != ConstraintVariable::Force) ? invmass : gmx::ArrayRef<real>(),
715 as_rvec_array(fp.data()));
720 for (int b = b0; b < b1; b++)
722 dhdlambda -= sol[b] * lincsd->ddist[b];
725 lincsd->task[th].dhdlambda = dhdlambda;
730 /* Constraint virial,
731 * determines sum r_bond x delta f,
732 * where delta f is the constraint correction
733 * of the quantity that is being constrained.
735 for (int b = b0; b < b1; b++)
737 const real mvb = lincsd->bllen[b] * sol[b];
738 for (int i = 0; i < DIM; i++)
740 const real tmp1 = mvb * r[b][i];
741 for (int j = 0; j < DIM; j++)
743 rmdf[i][j] += tmp1 * r[b][j];
746 } /* 23 ncons flops */
750 #if GMX_SIMD_HAVE_REAL
752 /*! \brief Calculate the constraint distance vectors r to project on from x.
754 * Determine the right-hand side of the matrix equation using coordinates xp. */
755 static void gmx_simdcall calc_dr_x_xp_simd(int b0,
757 gmx::ArrayRef<const AtomPair> atoms,
758 const rvec* gmx_restrict x,
759 const rvec* gmx_restrict xp,
760 const real* gmx_restrict bllen,
761 const real* gmx_restrict blc,
762 const real* pbc_simd,
763 rvec* gmx_restrict r,
764 real* gmx_restrict rhs,
765 real* gmx_restrict sol)
767 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
768 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
770 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
775 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
777 SimdReal x0_S, y0_S, z0_S;
778 SimdReal x1_S, y1_S, z1_S;
779 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
780 SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
781 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
782 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
784 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
786 offset0[i] = atoms[bs + i].index1;
787 offset1[i] = atoms[bs + i].index2;
790 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset0, &x0_S, &y0_S, &z0_S);
791 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset1, &x1_S, &y1_S, &z1_S);
796 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
798 n2_S = norm2(rx_S, ry_S, rz_S);
799 il_S = invsqrt(n2_S);
805 transposeScatterStoreU<3>(reinterpret_cast<real*>(r + bs), offset2, rx_S, ry_S, rz_S);
807 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(xp), offset0, &x0_S, &y0_S, &z0_S);
808 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(xp), offset1, &x1_S, &y1_S, &z1_S);
814 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
816 ip_S = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
818 rhs_S = load<SimdReal>(blc + bs) * (ip_S - load<SimdReal>(bllen + bs));
820 store(rhs + bs, rhs_S);
821 store(sol + bs, rhs_S);
824 #endif // GMX_SIMD_HAVE_REAL
826 /*! \brief Determine the distances and right-hand side for the next iteration. */
827 gmx_unused static void calc_dist_iter(int b0,
829 gmx::ArrayRef<const AtomPair> atoms,
830 const rvec* gmx_restrict xp,
831 const real* gmx_restrict bllen,
832 const real* gmx_restrict blc,
835 real* gmx_restrict rhs,
836 real* gmx_restrict sol,
839 for (int b = b0; b < b1; b++)
845 pbc_dx_aiuc(pbc, xp[atoms[b].index1], xp[atoms[b].index2], dx);
849 rvec_sub(xp[atoms[b].index1], xp[atoms[b].index2], dx);
851 real len2 = len * len;
852 real dlen2 = 2 * len2 - ::norm2(dx);
853 if (dlen2 < wfac * len2)
855 /* not race free - see detailed comment in caller */
861 mvb = blc[b] * (len - dlen2 * gmx::invsqrt(dlen2));
869 } /* 20*ncons flops */
872 #if GMX_SIMD_HAVE_REAL
873 /*! \brief As calc_dist_iter(), but using SIMD intrinsics. */
874 static void gmx_simdcall calc_dist_iter_simd(int b0,
876 gmx::ArrayRef<const AtomPair> atoms,
877 const rvec* gmx_restrict x,
878 const real* gmx_restrict bllen,
879 const real* gmx_restrict blc,
880 const real* pbc_simd,
882 real* gmx_restrict rhs,
883 real* gmx_restrict sol,
886 SimdReal min_S(GMX_REAL_MIN);
888 SimdReal wfac_S(wfac);
891 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
893 /* Initialize all to FALSE */
894 warn_S = (two_S < setZero());
896 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
898 SimdReal x0_S, y0_S, z0_S;
899 SimdReal x1_S, y1_S, z1_S;
900 SimdReal rx_S, ry_S, rz_S, n2_S;
901 SimdReal len_S, len2_S, dlen2_S, lc_S, blc_S;
902 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
903 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
905 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
907 offset0[i] = atoms[bs + i].index1;
908 offset1[i] = atoms[bs + i].index2;
911 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset0, &x0_S, &y0_S, &z0_S);
912 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset1, &x1_S, &y1_S, &z1_S);
918 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
920 n2_S = norm2(rx_S, ry_S, rz_S);
922 len_S = load<SimdReal>(bllen + bs);
923 len2_S = len_S * len_S;
925 dlen2_S = fms(two_S, len2_S, n2_S);
927 warn_S = warn_S || (dlen2_S < (wfac_S * len2_S));
929 /* Avoid 1/0 by taking the max with REAL_MIN.
930 * Note: when dlen2 is close to zero (90 degree constraint rotation),
931 * the accuracy of the algorithm is no longer relevant.
933 dlen2_S = max(dlen2_S, min_S);
935 lc_S = fnma(dlen2_S, invsqrt(dlen2_S), len_S);
937 blc_S = load<SimdReal>(blc + bs);
941 store(rhs + bs, lc_S);
942 store(sol + bs, lc_S);
950 #endif // GMX_SIMD_HAVE_REAL
952 //! Implements LINCS constraining.
953 static void do_lincs(ArrayRefWithPadding<const RVec> xPadded,
954 ArrayRefWithPadding<RVec> xpPadded,
959 ArrayRef<const real> invmass,
969 const rvec* x = as_rvec_array(xPadded.paddedArrayRef().data());
970 rvec* xp = as_rvec_array(xpPadded.paddedArrayRef().data());
971 rvec* gmx_restrict v = as_rvec_array(vRef.data());
973 const int b0 = lincsd->task[th].b0;
974 const int b1 = lincsd->task[th].b1;
976 gmx::ArrayRef<const AtomPair> atoms = lincsd->atoms;
977 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
978 gmx::ArrayRef<const int> blnr = lincsd->blnr;
979 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
980 gmx::ArrayRef<const real> blc = lincsd->blc;
981 gmx::ArrayRef<const real> blmf = lincsd->blmf;
982 gmx::ArrayRef<const real> bllen = lincsd->bllen;
983 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
984 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
985 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
986 gmx::ArrayRef<real> sol = lincsd->tmp3;
987 gmx::ArrayRef<real> blc_sol = lincsd->tmp4;
988 gmx::ArrayRef<real> mlambda = lincsd->mlambda;
989 gmx::ArrayRef<const int> nlocat = lincsd->nlocat;
991 #if GMX_SIMD_HAVE_REAL
993 /* This SIMD code does the same as the plain-C code after the #else.
994 * The only difference is that we always call pbc code, as with SIMD
995 * the overhead of pbc computation (when not needed) is small.
997 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9 * GMX_SIMD_REAL_WIDTH];
999 /* Convert the pbc struct for SIMD */
1000 set_pbc_simd(pbc, pbc_simd);
1002 /* Compute normalized x i-j vectors, store in r.
1003 * Compute the inner product of r and xp i-j and store in rhs1.
1006 b0, b1, atoms, x, xp, bllen.data(), blc.data(), pbc_simd, as_rvec_array(r.data()), rhs1.data(), sol.data());
1008 #else // GMX_SIMD_HAVE_REAL
1012 /* Compute normalized i-j vectors */
1013 for (int b = b0; b < b1; b++)
1016 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
1019 pbc_dx_aiuc(pbc, xp[atoms[b].index1], xp[atoms[b].index2], dx);
1020 real mvb = blc[b] * (::iprod(r[b], dx) - bllen[b]);
1027 /* Compute normalized i-j vectors */
1028 for (int b = b0; b < b1; b++)
1030 int i = atoms[b].index1;
1031 int j = atoms[b].index2;
1032 real tmp0 = x[i][0] - x[j][0];
1033 real tmp1 = x[i][1] - x[j][1];
1034 real tmp2 = x[i][2] - x[j][2];
1035 real rlen = gmx::invsqrt(tmp0 * tmp0 + tmp1 * tmp1 + tmp2 * tmp2);
1036 r[b][0] = rlen * tmp0;
1037 r[b][1] = rlen * tmp1;
1038 r[b][2] = rlen * tmp2;
1039 /* 16 ncons flops */
1042 * (r[b][0] * (xp[i][0] - xp[j][0]) + r[b][1] * (xp[i][1] - xp[j][1])
1043 + r[b][2] * (xp[i][2] - xp[j][2]) - bllen[b]);
1048 /* Together: 26*ncons + 6*nrtot flops */
1051 #endif // GMX_SIMD_HAVE_REAL
1053 if (lincsd->bTaskDep)
1055 /* We need a barrier, since the matrix construction below
1056 * can access entries in r of other threads.
1061 /* Construct the (sparse) LINCS matrix */
1062 for (int b = b0; b < b1; b++)
1064 for (int n = blnr[b]; n < blnr[b + 1]; n++)
1066 blcc[n] = blmf[n] * gmx::dot(r[b], r[blbnb[n]]);
1069 /* Together: 26*ncons + 6*nrtot flops */
1071 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
1072 /* nrec*(ncons+2*nrtot) flops */
1074 #if GMX_SIMD_HAVE_REAL
1075 for (int b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1077 SimdReal t1 = load<SimdReal>(blc.data() + b);
1078 SimdReal t2 = load<SimdReal>(sol.data() + b);
1079 store(mlambda.data() + b, t1 * t2);
1082 for (int b = b0; b < b1; b++)
1084 mlambda[b] = blc[b] * sol[b];
1086 #endif // GMX_SIMD_HAVE_REAL
1088 /* Update the coordinates */
1089 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1092 ******** Correction for centripetal effects ********
1097 wfac = std::cos(gmx::c_deg2Rad * wangle);
1100 for (int iter = 0; iter < lincsd->nIter; iter++)
1102 if ((lincsd->bCommIter && haveDDAtomOrdering(*cr) && cr->dd->constraints))
1107 /* Communicate the corrected non-local coordinates */
1108 if (haveDDAtomOrdering(*cr))
1110 dd_move_x_constraints(cr->dd, box, xpPadded.unpaddedArrayRef(), ArrayRef<RVec>(), FALSE);
1115 else if (lincsd->bTaskDep)
1120 #if GMX_SIMD_HAVE_REAL
1121 calc_dist_iter_simd(
1122 b0, b1, atoms, xp, bllen.data(), blc.data(), pbc_simd, wfac, rhs1.data(), sol.data(), bWarn);
1124 calc_dist_iter(b0, b1, atoms, xp, bllen.data(), blc.data(), pbc, wfac, rhs1.data(), sol.data(), bWarn);
1125 /* 20*ncons flops */
1126 #endif // GMX_SIMD_HAVE_REAL
1128 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
1129 /* nrec*(ncons+2*nrtot) flops */
1131 #if GMX_SIMD_HAVE_REAL
1132 for (int b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1134 SimdReal t1 = load<SimdReal>(blc.data() + b);
1135 SimdReal t2 = load<SimdReal>(sol.data() + b);
1136 SimdReal mvb = t1 * t2;
1137 store(blc_sol.data() + b, mvb);
1138 store(mlambda.data() + b, load<SimdReal>(mlambda.data() + b) + mvb);
1141 for (int b = b0; b < b1; b++)
1143 real mvb = blc[b] * sol[b];
1147 #endif // GMX_SIMD_HAVE_REAL
1149 /* Update the coordinates */
1150 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1152 /* nit*ncons*(37+9*nrec) flops */
1156 /* Update the velocities */
1157 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1158 /* 16 ncons flops */
1161 if (!nlocat.empty() && (bCalcDHDL || bCalcVir))
1163 if (lincsd->bTaskDep)
1165 /* In lincs_update_atoms threads might cross-read mlambda */
1169 /* Only account for local atoms */
1170 for (int b = b0; b < b1; b++)
1172 mlambda[b] *= 0.5 * nlocat[b];
1179 for (int b = b0; b < b1; b++)
1181 /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1182 * later after the contributions are reduced over the threads.
1184 dhdl -= lincsd->mlambda[b] * lincsd->ddist[b];
1186 lincsd->task[th].dhdlambda = dhdl;
1191 /* Constraint virial */
1192 for (int b = b0; b < b1; b++)
1194 real tmp0 = -bllen[b] * mlambda[b];
1195 for (int i = 0; i < DIM; i++)
1197 real tmp1 = tmp0 * r[b][i];
1198 for (int j = 0; j < DIM; j++)
1200 vir_r_m_dr[i][j] -= tmp1 * r[b][j];
1203 } /* 22 ncons flops */
1207 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1208 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1210 * (26+nrec)*ncons + (6+2*nrec)*nrtot
1211 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1213 * (63+nrec)*ncons + (6+4*nrec)*nrtot
1217 /*! \brief Sets the elements in the LINCS matrix for task task. */
1218 static void set_lincs_matrix_task(Lincs* li,
1220 ArrayRef<const real> invmass,
1222 int* nCrossTaskTriangles)
1224 /* Construct the coupling coefficient matrix blmf */
1225 li_task->ntriangle = 0;
1227 *nCrossTaskTriangles = 0;
1228 for (int i = li_task->b0; i < li_task->b1; i++)
1230 const int a1 = li->atoms[i].index1;
1231 const int a2 = li->atoms[i].index2;
1232 for (int n = li->blnr[i]; n < li->blnr[i + 1]; n++)
1234 const int k = li->blbnb[n];
1236 /* If we are using multiple, independent tasks for LINCS,
1237 * the calls to check_assign_connected should have
1238 * put all connected constraints in our task.
1240 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1243 if (a1 == li->atoms[k].index1 || a2 == li->atoms[k].index2)
1253 if (a1 == li->atoms[k].index1 || a1 == li->atoms[k].index2)
1263 li->blmf[n] = sign * invmass[center] * li->blc[i] * li->blc[k];
1264 li->blmf1[n] = sign * 0.5;
1265 if (li->ncg_triangle > 0)
1267 /* Look for constraint triangles */
1268 for (int nk = li->blnr[k]; nk < li->blnr[k + 1]; nk++)
1270 const int kk = li->blbnb[nk];
1271 if (kk != i && kk != k && (li->atoms[kk].index1 == end || li->atoms[kk].index2 == end))
1273 /* Check if the constraints in this triangle actually
1274 * belong to a different task. We still assign them
1275 * here, since it's convenient for the triangle
1276 * iterations, but we then need an extra barrier.
1278 if (k < li_task->b0 || k >= li_task->b1 || kk < li_task->b0 || kk >= li_task->b1)
1280 (*nCrossTaskTriangles)++;
1283 if (li_task->ntriangle == 0 || li_task->triangle[li_task->ntriangle - 1] < i)
1285 /* Add this constraint to the triangle list */
1286 li_task->triangle[li_task->ntriangle] = i;
1287 li_task->tri_bits[li_task->ntriangle] = 0;
1288 li_task->ntriangle++;
1289 if (li->blnr[i + 1] - li->blnr[i]
1290 > static_cast<int>(sizeof(li_task->tri_bits[0]) * 8 - 1))
1293 "A constraint is connected to %d constraints, this is "
1294 "more than the %zu allowed for constraints participating "
1296 li->blnr[i + 1] - li->blnr[i],
1297 sizeof(li_task->tri_bits[0]) * 8 - 1);
1300 li_task->tri_bits[li_task->ntriangle - 1] |= (1 << (n - li->blnr[i]));
1309 /*! \brief Sets the elements in the LINCS matrix. */
1310 static void set_lincs_matrix(Lincs* li, ArrayRef<const real> invmass, real lambda)
1312 const real invsqrt2 = 0.7071067811865475244;
1314 for (int i = 0; (i < li->nc); i++)
1316 const int a1 = li->atoms[i].index1;
1317 const int a2 = li->atoms[i].index2;
1318 li->blc[i] = gmx::invsqrt(invmass[a1] + invmass[a2]);
1319 li->blc1[i] = invsqrt2;
1322 /* Construct the coupling coefficient matrix blmf */
1323 int ntriangle = 0, ncc_triangle = 0, nCrossTaskTriangles = 0;
1324 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle, nCrossTaskTriangles) num_threads(li->ntask) schedule(static)
1325 for (int th = 0; th < li->ntask; th++)
1329 set_lincs_matrix_task(li, &li->task[th], invmass, &ncc_triangle, &nCrossTaskTriangles);
1330 ntriangle += li->task[th].ntriangle;
1332 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1334 li->ntriangle = ntriangle;
1335 li->ncc_triangle = ncc_triangle;
1336 li->bTaskDepTri = (nCrossTaskTriangles > 0);
1340 fprintf(debug, "The %d constraints participate in %d triangles\n", li->nc, li->ntriangle);
1341 fprintf(debug, "There are %d constraint couplings, of which %d in triangles\n", li->ncc, li->ncc_triangle);
1342 if (li->ntriangle > 0 && li->ntask > 1)
1345 "%d constraint triangles contain constraints assigned to different tasks\n",
1346 nCrossTaskTriangles);
1351 * so we know with which lambda value the masses have been set.
1353 li->matlam = lambda;
1356 //! Finds all triangles of atoms that share constraints to a central atom.
1357 static int count_triangle_constraints(const InteractionLists& ilist, const ListOfLists<int>& at2con)
1359 const int ncon1 = ilist[F_CONSTR].size() / 3;
1360 const int ncon_tot = ncon1 + ilist[F_CONSTRNC].size() / 3;
1362 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1363 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1365 int ncon_triangle = 0;
1366 for (int c0 = 0; c0 < ncon_tot; c0++)
1368 bool bTriangle = FALSE;
1369 const int* iap = constr_iatomptr(ia1, ia2, c0);
1370 const int a00 = iap[1];
1371 const int a01 = iap[2];
1372 for (const int c1 : at2con[a01])
1376 const int* iap = constr_iatomptr(ia1, ia2, c1);
1377 const int a10 = iap[1];
1378 const int a11 = iap[2];
1388 for (const int c2 : at2con[ac1])
1390 if (c2 != c0 && c2 != c1)
1392 const int* iap = constr_iatomptr(ia1, ia2, c2);
1393 const int a20 = iap[1];
1394 const int a21 = iap[2];
1395 if (a20 == a00 || a21 == a00)
1409 return ncon_triangle;
1412 //! Finds sequences of sequential constraints.
1413 static bool more_than_two_sequential_constraints(const InteractionLists& ilist, const ListOfLists<int>& at2con)
1415 const int ncon1 = ilist[F_CONSTR].size() / 3;
1416 const int ncon_tot = ncon1 + ilist[F_CONSTRNC].size() / 3;
1418 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1419 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1421 for (int c = 0; c < ncon_tot; c++)
1423 const int* iap = constr_iatomptr(ia1, ia2, c);
1424 const int a1 = iap[1];
1425 const int a2 = iap[2];
1426 /* Check if this constraint has constraints connected at both atoms */
1427 if (at2con[a1].ssize() > 1 && at2con[a2].ssize() > 1)
1436 Lincs* init_lincs(FILE* fplog,
1437 const gmx_mtop_t& mtop,
1438 int nflexcon_global,
1439 ArrayRef<const ListOfLists<int>> atomToConstraintsPerMolType,
1444 // TODO this should become a unique_ptr
1446 bool bMoreThanTwoSeq;
1450 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n", bPLINCS ? " Parallel" : "");
1455 li->ncg = gmx_mtop_ftype_count(mtop, F_CONSTR) + gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1456 li->ncg_flex = nflexcon_global;
1459 li->nOrder = nProjOrder;
1461 li->max_connect = 0;
1462 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1464 const auto& at2con = atomToConstraintsPerMolType[mt];
1465 for (int a = 0; a < mtop.moltype[mt].atoms.nr; a++)
1467 li->max_connect = std::max(li->max_connect, int(at2con[a].ssize()));
1471 li->ncg_triangle = 0;
1472 bMoreThanTwoSeq = FALSE;
1473 for (const gmx_molblock_t& molb : mtop.molblock)
1475 const gmx_moltype_t& molt = mtop.moltype[molb.type];
1476 const auto& at2con = atomToConstraintsPerMolType[molb.type];
1478 li->ncg_triangle += molb.nmol * count_triangle_constraints(molt.ilist, at2con);
1480 if (!bMoreThanTwoSeq && more_than_two_sequential_constraints(molt.ilist, at2con))
1482 bMoreThanTwoSeq = TRUE;
1486 /* Check if we need to communicate not only before LINCS,
1487 * but also before each iteration.
1488 * The check for only two sequential constraints is only
1489 * useful for the common case of H-bond only constraints.
1490 * With more effort we could also make it useful for small
1491 * molecules with nr. sequential constraints <= nOrder-1.
1493 li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1495 if (debug && bPLINCS)
1497 fprintf(debug, "PLINCS communication before each iteration: %d\n", static_cast<int>(li->bCommIter));
1500 /* LINCS can run on any number of threads.
1501 * Currently the number is fixed for the whole simulation,
1502 * but it could be set in set_lincs().
1503 * The current constraint to task assignment code can create independent
1504 * tasks only when not more than two constraints are connected sequentially.
1506 li->ntask = gmx_omp_nthreads_get(ModuleMultiThread::Lincs);
1507 li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1510 fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n", li->ntask, li->bTaskDep ? "" : "in");
1518 /* Allocate an extra elements for "task-overlap" constraints */
1519 li->task.resize(li->ntask + 1);
1522 if (bPLINCS || li->ncg_triangle > 0)
1524 please_cite(fplog, "Hess2008a");
1528 please_cite(fplog, "Hess97a");
1533 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1537 "There are constraints between atoms in different decomposition domains,\n"
1538 "will communicate selected coordinates each lincs iteration\n");
1540 if (li->ncg_triangle > 0)
1543 "%d constraints are involved in constraint triangles,\n"
1544 "will apply an additional matrix expansion of order %d for couplings\n"
1545 "between constraints inside triangles\n",
1554 void done_lincs(Lincs* li)
1559 /*! \brief Sets up the work division over the threads. */
1560 static void lincs_thread_setup(Lincs* li, int natoms)
1562 li->atf.resize(natoms);
1564 gmx::ArrayRef<gmx_bitmask_t> atf = li->atf;
1566 /* Clear the atom flags */
1567 for (gmx_bitmask_t& mask : atf)
1569 bitmask_clear(&mask);
1572 if (li->ntask > BITMASK_SIZE)
1574 gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1577 for (int th = 0; th < li->ntask; th++)
1579 const Task* li_task = &li->task[th];
1581 /* For each atom set a flag for constraints from each */
1582 for (int b = li_task->b0; b < li_task->b1; b++)
1584 bitmask_set_bit(&atf[li->atoms[b].index1], th);
1585 bitmask_set_bit(&atf[li->atoms[b].index2], th);
1589 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1590 for (int th = 0; th < li->ntask; th++)
1598 li_task = &li->task[th];
1600 bitmask_init_low_bits(&mask, th);
1602 li_task->ind.clear();
1603 li_task->ind_r.clear();
1604 for (b = li_task->b0; b < li_task->b1; b++)
1606 /* We let the constraint with the lowest thread index
1607 * operate on atoms with constraints from multiple threads.
1609 if (bitmask_is_disjoint(atf[li->atoms[b].index1], mask)
1610 && bitmask_is_disjoint(atf[li->atoms[b].index2], mask))
1612 /* Add the constraint to the local atom update index */
1613 li_task->ind.push_back(b);
1617 /* Add the constraint to the rest block */
1618 li_task->ind_r.push_back(b);
1622 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1625 /* We need to copy all constraints which have not be assigned
1626 * to a thread to a separate list which will be handled by one thread.
1628 Task* li_m = &li->task[li->ntask];
1631 for (int th = 0; th < li->ntask; th++)
1633 const Task& li_task = li->task[th];
1635 for (int ind_r : li_task.ind_r)
1637 li_m->ind.push_back(ind_r);
1642 fprintf(debug, "LINCS thread %d: %zu constraints\n", th, li_task.ind.size());
1648 fprintf(debug, "LINCS thread r: %zu constraints\n", li_m->ind.size());
1652 //! Assign a constraint.
1653 static void assign_constraint(Lincs* li,
1654 int constraint_index,
1659 const ListOfLists<int>& at2con)
1665 /* Make an mapping of local topology constraint index to LINCS index */
1666 li->con_index[constraint_index] = con;
1668 li->bllen0[con] = lenA;
1669 li->ddist[con] = lenB - lenA;
1670 /* Set the length to the topology A length */
1671 li->bllen[con] = lenA;
1672 li->atoms[con].index1 = a1;
1673 li->atoms[con].index2 = a2;
1675 /* Make space in the constraint connection matrix for constraints
1676 * connected to both end of the current constraint.
1678 li->ncc += at2con[a1].ssize() - 1 + at2con[a2].ssize() - 1;
1680 li->blnr[con + 1] = li->ncc;
1682 /* Increase the constraint count */
1686 /*! \brief Check if constraint with topology index constraint_index is connected
1687 * to other constraints, and if so add those connected constraints to our task. */
1688 static void check_assign_connected(Lincs* li,
1689 gmx::ArrayRef<const int> iatom,
1690 const InteractionDefinitions& idef,
1694 const ListOfLists<int>& at2con)
1696 /* Currently this function only supports constraint groups
1697 * in which all constraints share at least one atom
1698 * (e.g. H-bond constraints).
1699 * Check both ends of the current constraint for
1700 * connected constraints. We need to assign those
1703 for (int end = 0; end < 2; end++)
1705 const int a = (end == 0 ? a1 : a2);
1707 for (const int cc : at2con[a])
1709 /* Check if constraint cc has not yet been assigned */
1710 if (li->con_index[cc] == -1)
1712 const int type = iatom[cc * 3];
1713 const real lenA = idef.iparams[type].constr.dA;
1714 const real lenB = idef.iparams[type].constr.dB;
1716 if (bDynamics || lenA != 0 || lenB != 0)
1718 assign_constraint(li, cc, iatom[3 * cc + 1], iatom[3 * cc + 2], lenA, lenB, at2con);
1725 /*! \brief Check if constraint with topology index constraint_index is involved
1726 * in a constraint triangle, and if so add the other two constraints
1727 * in the triangle to our task. */
1728 static void check_assign_triangle(Lincs* li,
1729 gmx::ArrayRef<const int> iatom,
1730 const InteractionDefinitions& idef,
1732 int constraint_index,
1735 const ListOfLists<int>& at2con)
1737 int nca, cc[32], ca[32];
1738 int c_triangle[2] = { -1, -1 };
1741 for (const int c : at2con[a1])
1743 if (c != constraint_index)
1747 aa1 = iatom[c * 3 + 1];
1748 aa2 = iatom[c * 3 + 2];
1764 for (const int c : at2con[a2])
1766 if (c != constraint_index)
1770 aa1 = iatom[c * 3 + 1];
1771 aa2 = iatom[c * 3 + 2];
1774 for (i = 0; i < nca; i++)
1778 c_triangle[0] = cc[i];
1785 for (i = 0; i < nca; i++)
1789 c_triangle[0] = cc[i];
1797 if (c_triangle[0] >= 0)
1801 for (end = 0; end < 2; end++)
1803 /* Check if constraint c_triangle[end] has not yet been assigned */
1804 if (li->con_index[c_triangle[end]] == -1)
1809 i = c_triangle[end] * 3;
1811 lenA = idef.iparams[type].constr.dA;
1812 lenB = idef.iparams[type].constr.dB;
1814 if (bDynamics || lenA != 0 || lenB != 0)
1816 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1823 //! Sets matrix indices.
1824 static void set_matrix_indices(Lincs* li, const Task& li_task, const ListOfLists<int>& at2con, bool bSortMatrix)
1826 for (int b = li_task.b0; b < li_task.b1; b++)
1828 const int a1 = li->atoms[b].index1;
1829 const int a2 = li->atoms[b].index2;
1831 int i = li->blnr[b];
1832 for (const int constraint : at2con[a1])
1834 const int concon = li->con_index[constraint];
1837 li->blbnb[i++] = concon;
1840 for (const int constraint : at2con[a2])
1842 const int concon = li->con_index[constraint];
1845 li->blbnb[i++] = concon;
1851 /* Order the blbnb matrix to optimize memory access */
1852 std::sort(li->blbnb.begin() + li->blnr[b], li->blbnb.begin() + li->blnr[b + 1]);
1857 void set_lincs(const InteractionDefinitions& idef,
1859 ArrayRef<const real> invmass,
1862 const t_commrec* cr,
1868 /* Zero the thread index ranges.
1869 * Otherwise without local constraints we could return with old ranges.
1871 for (int i = 0; i < li->ntask; i++)
1875 li->task[i].ind.clear();
1879 li->task[li->ntask].ind.clear();
1882 /* This is the local topology, so there are only F_CONSTR constraints */
1883 if (idef.il[F_CONSTR].empty())
1885 /* There are no constraints,
1886 * we do not need to fill any data structures.
1893 fprintf(debug, "Building the LINCS connectivity\n");
1897 if (haveDDAtomOrdering(*cr))
1899 if (cr->dd->constraints)
1903 dd_get_constraint_range(*cr->dd, &start, &natoms);
1907 natoms = dd_numHomeAtoms(*cr->dd);
1915 const ListOfLists<int> at2con =
1916 make_at2con(natoms, idef.il, idef.iparams, flexibleConstraintTreatment(bDynamics));
1918 const int ncon_tot = idef.il[F_CONSTR].size() / 3;
1920 /* Ensure we have enough padding for aligned loads for each thread */
1921 const int numEntries = ncon_tot + li->ntask * simd_width;
1922 li->con_index.resize(numEntries);
1923 li->bllen0.resize(numEntries);
1924 li->ddist.resize(numEntries);
1925 li->atoms.resize(numEntries);
1926 li->blc.resize(numEntries);
1927 li->blc1.resize(numEntries);
1928 li->blnr.resize(numEntries + 1);
1929 li->bllen.resize(numEntries);
1930 li->tmpv.resizeWithPadding(numEntries);
1931 if (haveDDAtomOrdering(*cr))
1933 li->nlocat.resize(numEntries);
1935 li->tmp1.resize(numEntries);
1936 li->tmp2.resize(numEntries);
1937 li->tmp3.resize(numEntries);
1938 li->tmp4.resize(numEntries);
1939 li->mlambda.resize(numEntries);
1941 gmx::ArrayRef<const int> iatom = idef.il[F_CONSTR].iatoms;
1943 li->blnr[0] = li->ncc;
1945 /* Assign the constraints for li->ntask LINCS tasks.
1946 * We target a uniform distribution of constraints over the tasks.
1947 * Note that when flexible constraints are present, but are removed here
1948 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
1949 * happen during normal MD, that's ok.
1952 /* Determine the number of constraints we need to assign here */
1953 int ncon_assign = ncon_tot;
1956 /* With energy minimization, flexible constraints are ignored
1957 * (and thus minimized, as they should be).
1959 ncon_assign -= countFlexibleConstraints(idef.il, idef.iparams);
1962 /* Set the target constraint count per task to exactly uniform,
1963 * this might be overridden below.
1965 int ncon_target = (ncon_assign + li->ntask - 1) / li->ntask;
1967 /* Mark all constraints as unassigned by setting their index to -1 */
1968 for (int con = 0; con < ncon_tot; con++)
1970 li->con_index[con] = -1;
1974 for (int th = 0; th < li->ntask; th++)
1978 li_task = &li->task[th];
1980 #if GMX_SIMD_HAVE_REAL
1981 /* With indepedent tasks we likely have H-bond constraints or constraint
1982 * pairs. The connected constraints will be pulled into the task, so the
1983 * constraints per task will often exceed ncon_target.
1984 * Triangle constraints can also increase the count, but there are
1985 * relatively few of those, so we usually expect to get ncon_target.
1989 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
1990 * since otherwise a lot of operations can be wasted.
1991 * There are several ways to round here, we choose the one
1992 * that alternates block sizes, which helps with Intel HT.
1994 ncon_target = ((ncon_assign * (th + 1)) / li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1)
1995 & ~(GMX_SIMD_REAL_WIDTH - 1);
1997 #endif // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
1999 /* Continue filling the arrays where we left off with the previous task,
2000 * including padding for SIMD.
2002 li_task->b0 = li->nc;
2004 gmx::ArrayRef<const t_iparams> iparams = idef.iparams;
2006 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2008 if (li->con_index[con] == -1)
2013 type = iatom[3 * con];
2014 a1 = iatom[3 * con + 1];
2015 a2 = iatom[3 * con + 2];
2016 lenA = iparams[type].constr.dA;
2017 lenB = iparams[type].constr.dB;
2018 /* Skip the flexible constraints when not doing dynamics */
2019 if (bDynamics || lenA != 0 || lenB != 0)
2021 assign_constraint(li, con, a1, a2, lenA, lenB, at2con);
2023 if (li->ntask > 1 && !li->bTaskDep)
2025 /* We can generate independent tasks. Check if we
2026 * need to assign connected constraints to our task.
2028 check_assign_connected(li, iatom, idef, bDynamics, a1, a2, at2con);
2030 if (li->ntask > 1 && li->ncg_triangle > 0)
2032 /* Ensure constraints in one triangle are assigned
2035 check_assign_triangle(li, iatom, idef, bDynamics, con, a1, a2, at2con);
2043 li_task->b1 = li->nc;
2047 /* Copy the last atom pair indices and lengths for constraints
2048 * up to a multiple of simd_width, such that we can do all
2049 * SIMD operations without having to worry about end effects.
2053 li->nc = ((li_task->b1 + simd_width - 1) / simd_width) * simd_width;
2054 last = li_task->b1 - 1;
2055 for (i = li_task->b1; i < li->nc; i++)
2057 li->atoms[i] = li->atoms[last];
2058 li->bllen0[i] = li->bllen0[last];
2059 li->ddist[i] = li->ddist[last];
2060 li->bllen[i] = li->bllen[last];
2061 li->blnr[i + 1] = li->blnr[last + 1];
2065 /* Keep track of how many constraints we assigned */
2066 li->nc_real += li_task->b1 - li_task->b0;
2070 fprintf(debug, "LINCS task %d constraints %d - %d\n", th, li_task->b0, li_task->b1);
2074 assert(li->nc_real == ncon_assign);
2078 /* Without DD we order the blbnb matrix to optimize memory access.
2079 * With DD the overhead of sorting is more than the gain during access.
2081 bSortMatrix = !haveDDAtomOrdering(*cr);
2083 li->blbnb.resize(li->ncc);
2085 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2086 for (int th = 0; th < li->ntask; th++)
2090 Task& li_task = li->task[th];
2092 if (li->ncg_triangle > 0)
2094 /* This is allocating too much, but it is difficult to improve */
2095 li_task.triangle.resize(li_task.b1 - li_task.b0);
2096 li_task.tri_bits.resize(li_task.b1 - li_task.b0);
2099 set_matrix_indices(li, li_task, at2con, bSortMatrix);
2101 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2104 if (cr->dd == nullptr)
2106 /* Since the matrix is static, we should free some memory */
2107 li->blbnb.resize(li->ncc);
2110 li->blmf.resize(li->ncc);
2111 li->blmf1.resize(li->ncc);
2112 li->tmpncc.resize(li->ncc);
2114 gmx::ArrayRef<const int> nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2115 if (!nlocat_dd.empty())
2117 /* Convert nlocat from local topology to LINCS constraint indexing */
2118 for (con = 0; con < ncon_tot; con++)
2120 li->nlocat[li->con_index[con]] = nlocat_dd[con];
2130 fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n", li->nc_real, li->nc, li->ncc);
2135 lincs_thread_setup(li, numAtoms);
2138 set_lincs_matrix(li, invmass, lambda);
2140 li->rmsdData[0] = 0.0;
2141 li->rmsdData[1] = 0.0;
2144 //! Issues a warning when LINCS constraints cannot be satisfied.
2145 static void lincs_warning(gmx_domdec_t* dd,
2146 ArrayRef<const RVec> x,
2147 ArrayRef<const RVec> xprime,
2150 gmx::ArrayRef<const AtomPair> atoms,
2151 gmx::ArrayRef<const real> bllen,
2156 real wfac = std::cos(gmx::c_deg2Rad * wangle);
2159 "bonds that rotated more than %g degrees:\n"
2160 " atom 1 atom 2 angle previous, current, constraint length\n",
2163 for (int b = 0; b < ncons; b++)
2165 const int i = atoms[b].index1;
2166 const int j = atoms[b].index2;
2171 pbc_dx_aiuc(pbc, x[i], x[j], v0);
2172 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2176 rvec_sub(x[i], x[j], v0);
2177 rvec_sub(xprime[i], xprime[j], v1);
2181 real cosine = ::iprod(v0, v1) / (d0 * d1);
2185 " %6d %6d %5.1f %8.4f %8.4f %8.4f\n",
2188 gmx::c_rad2Deg * std::acos(cosine),
2192 if (!std::isfinite(d1))
2194 gmx_fatal(FARGS, "Bond length not finite.");
2200 if (*warncount > maxwarn)
2202 too_many_constraint_warnings(ConstraintAlgorithm::Lincs, *warncount);
2206 //! Status information about how well LINCS satisified the constraints in this domain
2207 struct LincsDeviations
2209 //! The maximum over all bonds in this domain of the relative deviation in bond lengths
2210 real maxDeviation = 0;
2211 //! Sum over all bonds in this domain of the squared relative deviation
2212 real sumSquaredDeviation = 0;
2213 //! Index of bond with max deviation
2214 int indexOfMaxDeviation = -1;
2215 //! Number of bonds constrained in this domain
2216 int numConstraints = 0;
2219 //! Determine how well the constraints have been satisfied.
2220 static LincsDeviations makeLincsDeviations(const Lincs& lincsd, ArrayRef<const RVec> x, const t_pbc* pbc)
2222 LincsDeviations result;
2223 const ArrayRef<const AtomPair> atoms = lincsd.atoms;
2224 const ArrayRef<const real> bllen = lincsd.bllen;
2225 const ArrayRef<const int> nlocat = lincsd.nlocat;
2227 for (int task = 0; task < lincsd.ntask; task++)
2229 for (int b = lincsd.task[task].b0; b < lincsd.task[task].b1; b++)
2234 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
2238 rvec_sub(x[atoms[b].index1], x[atoms[b].index2], dx);
2240 real r2 = ::norm2(dx);
2241 real len = r2 * gmx::invsqrt(r2);
2242 real d = std::abs(len / bllen[b] - 1.0_real);
2243 if (d > result.maxDeviation && (nlocat.empty() || nlocat[b]))
2245 result.maxDeviation = d;
2246 result.indexOfMaxDeviation = b;
2250 result.sumSquaredDeviation += d * d;
2251 result.numConstraints++;
2255 result.sumSquaredDeviation += nlocat[b] * d * d;
2256 result.numConstraints += nlocat[b];
2261 if (!nlocat.empty())
2263 result.numConstraints /= 2;
2264 result.sumSquaredDeviation *= 0.5;
2269 bool constrain_lincs(bool computeRmsd,
2270 const t_inputrec& ir,
2273 ArrayRef<const real> invmass,
2274 const t_commrec* cr,
2275 const gmx_multisim_t* ms,
2276 ArrayRefWithPadding<const RVec> xPadded,
2277 ArrayRefWithPadding<RVec> xprimePadded,
2278 ArrayRef<RVec> min_proj,
2281 const bool hasMassPerturbed,
2288 ConstraintVariable econq,
2295 /* This boolean should be set by a flag passed to this routine.
2296 * We can also easily check if any constraint length is changed,
2297 * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2299 bool bCalcDHDL = (ir.efep != FreeEnergyPerturbationType::No && dvdlambda != nullptr);
2301 if (lincsd->nc == 0 && cr->dd == nullptr)
2305 lincsd->rmsdData = { { 0 } };
2311 ArrayRef<const RVec> x = xPadded.unpaddedArrayRef();
2312 ArrayRef<RVec> xprime = xprimePadded.unpaddedArrayRef();
2314 if (econq == ConstraintVariable::Positions)
2316 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2317 * also with efep!=fepNO.
2319 if (ir.efep != FreeEnergyPerturbationType::No)
2321 if (hasMassPerturbed && lincsd->matlam != lambda)
2323 set_lincs_matrix(lincsd, invmass, lambda);
2326 for (int i = 0; i < lincsd->nc; i++)
2328 lincsd->bllen[i] = lincsd->bllen0[i] + lambda * lincsd->ddist[i];
2332 if (lincsd->ncg_flex)
2334 /* Set the flexible constraint lengths to the old lengths */
2337 for (int i = 0; i < lincsd->nc; i++)
2339 if (lincsd->bllen[i] == 0)
2342 pbc_dx_aiuc(pbc, x[lincsd->atoms[i].index1], x[lincsd->atoms[i].index2], dx);
2343 lincsd->bllen[i] = norm(dx);
2349 for (int i = 0; i < lincsd->nc; i++)
2351 if (lincsd->bllen[i] == 0)
2353 lincsd->bllen[i] = std::sqrt(
2354 distance2(x[lincsd->atoms[i].index1], x[lincsd->atoms[i].index2]));
2360 const bool printDebugOutput = ((debug != nullptr) && lincsd->nc > 0);
2361 if (printDebugOutput)
2363 LincsDeviations deviations = makeLincsDeviations(*lincsd, xprime, pbc);
2364 fprintf(debug, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2366 " Before LINCS %.6f %.6f %6d %6d\n",
2367 std::sqrt(deviations.sumSquaredDeviation / deviations.numConstraints),
2368 deviations.maxDeviation,
2369 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index1),
2370 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index2));
2373 /* This bWarn var can be updated by multiple threads
2374 * at the same time. But as we only need to detect
2375 * if a warning occurred or not, this is not an issue.
2379 /* The OpenMP parallel region of constrain_lincs for coords */
2380 #pragma omp parallel num_threads(lincsd->ntask)
2384 int th = gmx_omp_get_thread_num();
2386 clear_mat(lincsd->task[th].vir_r_m_dr);
2402 th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2404 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2407 if (computeRmsd || printDebugOutput || bWarn)
2409 LincsDeviations deviations = makeLincsDeviations(*lincsd, xprime, pbc);
2413 // This is reduced across domains in compute_globals and
2414 // reported to the log file.
2415 lincsd->rmsdData[0] = deviations.numConstraints;
2416 lincsd->rmsdData[1] = deviations.sumSquaredDeviation;
2420 // This is never read
2421 lincsd->rmsdData = { { 0 } };
2423 if (printDebugOutput)
2426 " After LINCS %.6f %.6f %6d %6d\n\n",
2427 std::sqrt(deviations.sumSquaredDeviation / deviations.numConstraints),
2428 deviations.maxDeviation,
2429 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index1),
2430 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index2));
2435 if (maxwarn < INT_MAX)
2437 std::string simMesg;
2440 simMesg += gmx::formatString(" in simulation %d", ms->simulationIndex_);
2444 ", time %g (ps) LINCS WARNING%s\n"
2445 "relative constraint deviation after LINCS:\n"
2446 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2448 ir.init_t + step * ir.delta_t,
2450 std::sqrt(deviations.sumSquaredDeviation / deviations.numConstraints),
2451 deviations.maxDeviation,
2452 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index1),
2453 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index2));
2456 cr->dd, x, xprime, pbc, lincsd->nc, lincsd->atoms, lincsd->bllen, ir.LincsWarnAngle, maxwarn, warncount);
2458 bOK = (deviations.maxDeviation < 0.5);
2462 if (lincsd->ncg_flex)
2464 for (int i = 0; (i < lincsd->nc); i++)
2466 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2468 lincsd->bllen[i] = 0;
2475 /* The OpenMP parallel region of constrain_lincs for derivatives */
2476 #pragma omp parallel num_threads(lincsd->ntask)
2480 int th = gmx_omp_get_thread_num();
2492 th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2494 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2500 /* Reduce the dH/dlambda contributions over the threads */
2505 for (th = 0; th < lincsd->ntask; th++)
2507 dhdlambda += lincsd->task[th].dhdlambda;
2509 if (econq == ConstraintVariable::Positions)
2511 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2512 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2513 dhdlambda /= ir.delta_t * ir.delta_t;
2515 *dvdlambda += dhdlambda;
2518 if (bCalcVir && lincsd->ntask > 1)
2520 for (int i = 1; i < lincsd->ntask; i++)
2522 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2526 /* count assuming nit=1 */
2527 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2528 inc_nrnb(nrnb, eNR_LINCSMAT, (2 + lincsd->nOrder) * lincsd->ncc);
2529 if (lincsd->ntriangle > 0)
2531 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder * lincsd->ncc_triangle);
2535 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real * 2);
2539 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);