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, 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"
88 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
93 //! Indices of the two atoms involved in a single constraint
96 //! \brief Constructor, does not initialize to catch bugs and faster construction
105 //! Unit of work within LINCS.
108 //! First constraint for this task.
110 //! b1-1 is the last constraint for this task.
112 //! The number of constraints in triangles.
114 //! The list of triangle constraints.
115 std::vector<int> triangle;
116 //! The bits tell if the matrix element should be used.
117 std::vector<int> tri_bits;
118 //! Constraint index for updating atom data.
119 std::vector<int> ind;
120 //! Constraint index for updating atom data.
121 std::vector<int> ind_r;
122 //! Temporary variable for virial calculation.
123 tensor vir_r_m_dr = { { 0 } };
124 //! Temporary variable for lambda derivative.
128 /*! \brief Data for LINCS algorithm.
133 //! The global number of constraints.
135 //! The global number of flexible constraints.
137 //! The global number of constraints in triangles.
138 int ncg_triangle = 0;
139 //! The number of iterations.
141 //! The order of the matrix expansion.
143 //! The maximum number of constraints connected to a single atom.
146 //! The number of real constraints.
148 //! The number of constraints including padding for SIMD.
150 //! The number of constraint connections.
152 //! The FE lambda value used for filling blc and blmf.
154 //! mapping from topology to LINCS constraints.
155 std::vector<int> con_index;
156 //! The reference distance in topology A.
157 std::vector<real, AlignedAllocator<real>> bllen0;
158 //! The reference distance in top B - the r.d. in top A.
159 std::vector<real, AlignedAllocator<real>> ddist;
160 //! The atom pairs involved in the constraints.
161 std::vector<AtomPair> atoms;
162 //! 1/sqrt(invmass1 invmass2).
163 std::vector<real, AlignedAllocator<real>> blc;
164 //! As blc, but with all masses 1.
165 std::vector<real, AlignedAllocator<real>> blc1;
166 //! Index into blbnb and blmf.
167 std::vector<int> blnr;
168 //! List of constraint connections.
169 std::vector<int> blbnb;
170 //! The local number of constraints in triangles.
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 std::vector<real> blmf;
178 //! As blmf, but with all masses 1.
179 std::vector<real> blmf1;
180 //! The reference bond length.
181 std::vector<real, AlignedAllocator<real>> bllen;
182 //! The local atom count per constraint, can be NULL.
183 std::vector<int> nlocat;
185 /*! \brief The number of tasks used for LINCS work.
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,
193 /*! \brief LINCS thread division */
194 std::vector<Task> task;
195 //! Atom flags for thread parallelization.
196 std::vector<gmx_bitmask_t> atf;
197 //! Are the LINCS tasks interdependent?
198 bool bTaskDep = false;
199 //! Are there triangle constraints that cross task borders?
200 bool bTaskDepTri = false;
201 //! Arrays for temporary storage in the LINCS algorithm.
203 PaddedVector<gmx::RVec> tmpv;
204 std::vector<real> tmpncc;
205 std::vector<real, AlignedAllocator<real>> tmp1;
206 std::vector<real, AlignedAllocator<real>> tmp2;
207 std::vector<real, AlignedAllocator<real>> tmp3;
208 std::vector<real, AlignedAllocator<real>> tmp4;
210 //! The Lagrange multipliers times -1.
211 std::vector<real, AlignedAllocator<real>> mlambda;
212 //! Storage for the constraint RMS relative deviation output.
213 std::array<real, 2> rmsdData = { { 0 } };
216 /*! \brief Define simd_width for memory allocation used for SIMD code */
217 #if GMX_SIMD_HAVE_REAL
218 static const int simd_width = GMX_SIMD_REAL_WIDTH;
220 static const int simd_width = 1;
223 ArrayRef<real> lincs_rmsdData(Lincs* lincsd)
225 return lincsd->rmsdData;
228 real lincs_rmsd(const Lincs* lincsd)
230 if (lincsd->rmsdData[0] > 0)
232 return std::sqrt(lincsd->rmsdData[1] / lincsd->rmsdData[0]);
240 /*! \brief Do a set of nrec LINCS matrix multiplications.
242 * This function will return with up to date thread-local
243 * constraint data, without an OpenMP barrier.
245 static void lincs_matrix_expand(const Lincs& lincsd,
247 gmx::ArrayRef<const real> blcc,
248 gmx::ArrayRef<real> rhs1,
249 gmx::ArrayRef<real> rhs2,
250 gmx::ArrayRef<real> sol)
252 gmx::ArrayRef<const int> blnr = lincsd.blnr;
253 gmx::ArrayRef<const int> blbnb = lincsd.blbnb;
255 const int b0 = li_task.b0;
256 const int b1 = li_task.b1;
257 const int nrec = lincsd.nOrder;
259 for (int rec = 0; rec < nrec; rec++)
265 for (int b = b0; b < b1; b++)
271 for (n = blnr[b]; n < blnr[b + 1]; n++)
273 mvb = mvb + blcc[n] * rhs1[blbnb[n]];
276 sol[b] = sol[b] + mvb;
279 std::swap(rhs1, rhs2);
280 } /* nrec*(ncons+2*nrtot) flops */
282 if (lincsd.ntriangle > 0)
284 /* Perform an extra nrec recursions for only the constraints
285 * involved in rigid triangles.
286 * In this way their accuracy should come close to those of the other
287 * constraints, since traingles of constraints can produce eigenvalues
288 * around 0.7, while the effective eigenvalue for bond constraints
289 * is around 0.4 (and 0.7*0.7=0.5).
294 /* We need a barrier here, since other threads might still be
295 * reading the contents of rhs1 and/o rhs2.
296 * We could avoid this barrier by introducing two extra rhs
297 * arrays for the triangle constraints only.
302 /* Constraints involved in a triangle are ensured to be in the same
303 * LINCS task. This means no barriers are required during the extra
304 * iterations for the triangle constraints.
306 gmx::ArrayRef<const int> triangle = li_task.triangle;
307 gmx::ArrayRef<const int> tri_bits = li_task.tri_bits;
309 for (int rec = 0; rec < nrec; rec++)
311 for (int tb = 0; tb < li_task.ntriangle; tb++)
313 int b, bits, nr0, nr1, n;
321 for (n = nr0; n < nr1; n++)
323 if (bits & (1 << (n - nr0)))
325 mvb = mvb + blcc[n] * rhs1[blbnb[n]];
329 sol[b] = sol[b] + mvb;
332 std::swap(rhs1, rhs2);
333 } /* nrec*(ntriangle + ncc_triangle*2) flops */
335 if (lincsd.bTaskDepTri)
337 /* The constraints triangles are decoupled from each other,
338 * but constraints in one triangle cross thread task borders.
339 * We could probably avoid this with more advanced setup code.
346 //! Update atomic coordinates when an index is not required.
347 static void lincs_update_atoms_noind(int ncons,
348 gmx::ArrayRef<const AtomPair> atoms,
350 gmx::ArrayRef<const real> fac,
351 gmx::ArrayRef<const gmx::RVec> r,
355 if (invmass != nullptr)
357 for (int b = 0; b < ncons; b++)
359 int i = atoms[b].index1;
360 int j = atoms[b].index2;
361 real mvb = preFactor * fac[b];
362 real im1 = invmass[i];
363 real im2 = invmass[j];
364 real tmp0 = r[b][0] * mvb;
365 real tmp1 = r[b][1] * mvb;
366 real tmp2 = r[b][2] * mvb;
367 x[i][0] -= tmp0 * im1;
368 x[i][1] -= tmp1 * im1;
369 x[i][2] -= tmp2 * im1;
370 x[j][0] += tmp0 * im2;
371 x[j][1] += tmp1 * im2;
372 x[j][2] += tmp2 * im2;
373 } /* 16 ncons flops */
377 for (int b = 0; b < ncons; b++)
379 int i = atoms[b].index1;
380 int j = atoms[b].index2;
381 real mvb = preFactor * fac[b];
382 real tmp0 = r[b][0] * mvb;
383 real tmp1 = r[b][1] * mvb;
384 real tmp2 = r[b][2] * mvb;
395 //! Update atomic coordinates when an index is required.
396 static void lincs_update_atoms_ind(gmx::ArrayRef<const int> ind,
397 gmx::ArrayRef<const AtomPair> atoms,
399 gmx::ArrayRef<const real> fac,
400 gmx::ArrayRef<const gmx::RVec> r,
404 if (invmass != nullptr)
408 int i = atoms[b].index1;
409 int j = atoms[b].index2;
410 real mvb = preFactor * fac[b];
411 real im1 = invmass[i];
412 real im2 = invmass[j];
413 real tmp0 = r[b][0] * mvb;
414 real tmp1 = r[b][1] * mvb;
415 real tmp2 = r[b][2] * mvb;
416 x[i][0] -= tmp0 * im1;
417 x[i][1] -= tmp1 * im1;
418 x[i][2] -= tmp2 * im1;
419 x[j][0] += tmp0 * im2;
420 x[j][1] += tmp1 * im2;
421 x[j][2] += tmp2 * im2;
422 } /* 16 ncons flops */
428 int i = atoms[b].index1;
429 int j = atoms[b].index2;
430 real mvb = preFactor * fac[b];
431 real tmp0 = r[b][0] * mvb;
432 real tmp1 = r[b][1] * mvb;
433 real tmp2 = r[b][2] * mvb;
440 } /* 16 ncons flops */
444 //! Update coordinates for atoms.
445 static void lincs_update_atoms(Lincs* li,
448 gmx::ArrayRef<const real> fac,
449 gmx::ArrayRef<const gmx::RVec> r,
455 /* Single thread, we simply update for all constraints */
456 lincs_update_atoms_noind(li->nc_real, li->atoms, preFactor, fac, r, invmass, x);
460 /* Update the atom vector components for our thread local
461 * constraints that only access our local atom range.
462 * This can be done without a barrier.
464 lincs_update_atoms_ind(li->task[th].ind, li->atoms, preFactor, fac, r, invmass, x);
466 if (!li->task[li->ntask].ind.empty())
468 /* Update the constraints that operate on atoms
469 * in multiple thread atom blocks on the master thread.
474 lincs_update_atoms_ind(li->task[li->ntask].ind, li->atoms, preFactor, fac, r, invmass, x);
480 #if GMX_SIMD_HAVE_REAL
481 //! Helper function so that we can run TSAN with SIMD support (where implemented).
483 static inline void gmx_simdcall gatherLoadUTransposeTSANSafe(const real* base,
484 const std::int32_t* offset,
489 # if (CMAKE_BUILD_TYPE == CMAKE_BUILD_TYPE_TSAN) && GMX_SIMD_X86_AVX2_256
490 // This function is only implemented in this case
491 gatherLoadUTransposeSafe<align>(base, offset, v0, v1, v2);
493 gatherLoadUTranspose<align>(base, offset, v0, v1, v2);
497 /*! \brief Calculate the constraint distance vectors r to project on from x.
499 * Determine the right-hand side of the matrix equation using quantity f.
500 * This function only differs from calc_dr_x_xp_simd below in that
501 * no constraint length is subtracted and no PBC is used for f. */
502 static void gmx_simdcall calc_dr_x_f_simd(int b0,
504 gmx::ArrayRef<const AtomPair> atoms,
505 const rvec* gmx_restrict x,
506 const rvec* gmx_restrict f,
507 const real* gmx_restrict blc,
508 const real* pbc_simd,
509 rvec* gmx_restrict r,
510 real* gmx_restrict rhs,
511 real* gmx_restrict sol)
513 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
515 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
517 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
522 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
524 SimdReal x0_S, y0_S, z0_S;
525 SimdReal x1_S, y1_S, z1_S;
526 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
527 SimdReal fx_S, fy_S, fz_S, ip_S, rhs_S;
528 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
529 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
531 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
533 offset0[i] = atoms[bs + i].index1;
534 offset1[i] = atoms[bs + i].index2;
537 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset0, &x0_S, &y0_S, &z0_S);
538 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset1, &x1_S, &y1_S, &z1_S);
543 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
545 n2_S = norm2(rx_S, ry_S, rz_S);
546 il_S = invsqrt(n2_S);
552 transposeScatterStoreU<3>(reinterpret_cast<real*>(r + bs), offset2, rx_S, ry_S, rz_S);
554 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(f), offset0, &x0_S, &y0_S, &z0_S);
555 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(f), offset1, &x1_S, &y1_S, &z1_S);
560 ip_S = iprod(rx_S, ry_S, rz_S, fx_S, fy_S, fz_S);
562 rhs_S = load<SimdReal>(blc + bs) * ip_S;
564 store(rhs + bs, rhs_S);
565 store(sol + bs, rhs_S);
568 #endif // GMX_SIMD_HAVE_REAL
570 /*! \brief LINCS projection, works on derivatives of the coordinates. */
571 static void do_lincsp(ArrayRefWithPadding<const RVec> xPadded,
572 ArrayRefWithPadding<RVec> fPadded,
578 ConstraintVariable econq,
583 const int b0 = lincsd->task[th].b0;
584 const int b1 = lincsd->task[th].b1;
586 gmx::ArrayRef<const AtomPair> atoms = lincsd->atoms;
587 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
588 gmx::ArrayRef<const int> blnr = lincsd->blnr;
589 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
591 gmx::ArrayRef<const real> blc;
592 gmx::ArrayRef<const real> blmf;
593 if (econq != ConstraintVariable::Force)
595 /* Use mass-weighted parameters */
601 /* Use non mass-weighted parameters */
603 blmf = lincsd->blmf1;
605 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
606 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
607 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
608 gmx::ArrayRef<real> sol = lincsd->tmp3;
610 const rvec* x = as_rvec_array(xPadded.paddedArrayRef().data());
611 rvec* f = as_rvec_array(fPadded.paddedArrayRef().data());
613 #if GMX_SIMD_HAVE_REAL
614 /* This SIMD code does the same as the plain-C code after the #else.
615 * The only difference is that we always call pbc code, as with SIMD
616 * the overhead of pbc computation (when not needed) is small.
618 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9 * GMX_SIMD_REAL_WIDTH];
620 /* Convert the pbc struct for SIMD */
621 set_pbc_simd(pbc, pbc_simd);
623 /* Compute normalized x i-j vectors, store in r.
624 * Compute the inner product of r and xp i-j and store in rhs1.
626 calc_dr_x_f_simd(b0, b1, atoms, x, f, blc.data(), pbc_simd, as_rvec_array(r.data()),
627 rhs1.data(), sol.data());
629 #else // GMX_SIMD_HAVE_REAL
631 /* Compute normalized i-j vectors */
634 for (int b = b0; b < b1; b++)
638 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
644 for (int b = b0; b < b1; b++)
648 rvec_sub(x[atoms[b].index1], x[atoms[b].index2], dx);
650 } /* 16 ncons flops */
653 for (int b = b0; b < b1; b++)
655 int i = atoms[b].index1;
656 int j = atoms[b].index2;
658 * (r[b][0] * (f[i][0] - f[j][0]) + r[b][1] * (f[i][1] - f[j][1])
659 + r[b][2] * (f[i][2] - f[j][2]));
665 #endif // GMX_SIMD_HAVE_REAL
667 if (lincsd->bTaskDep)
669 /* We need a barrier, since the matrix construction below
670 * can access entries in r of other threads.
675 /* Construct the (sparse) LINCS matrix */
676 for (int b = b0; b < b1; b++)
678 for (int n = blnr[b]; n < blnr[b + 1]; n++)
680 blcc[n] = blmf[n] * ::iprod(r[b], r[blbnb[n]]);
683 /* Together: 23*ncons + 6*nrtot flops */
685 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
686 /* nrec*(ncons+2*nrtot) flops */
688 if (econq == ConstraintVariable::Deriv_FlexCon)
690 /* We only want to constraint the flexible constraints,
691 * so we mask out the normal ones by setting sol to 0.
693 for (int b = b0; b < b1; b++)
695 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
702 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
703 for (int b = b0; b < b1; b++)
708 /* When constraining forces, we should not use mass weighting,
709 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
711 lincs_update_atoms(lincsd, th, 1.0, sol, r, (econq != ConstraintVariable::Force) ? invmass : nullptr,
712 as_rvec_array(fp.data()));
717 for (int b = b0; b < b1; b++)
719 dhdlambda -= sol[b] * lincsd->ddist[b];
722 lincsd->task[th].dhdlambda = dhdlambda;
727 /* Constraint virial,
728 * determines sum r_bond x delta f,
729 * where delta f is the constraint correction
730 * of the quantity that is being constrained.
732 for (int b = b0; b < b1; b++)
734 const real mvb = lincsd->bllen[b] * sol[b];
735 for (int i = 0; i < DIM; i++)
737 const real tmp1 = mvb * r[b][i];
738 for (int j = 0; j < DIM; j++)
740 rmdf[i][j] += tmp1 * r[b][j];
743 } /* 23 ncons flops */
747 #if GMX_SIMD_HAVE_REAL
749 /*! \brief Calculate the constraint distance vectors r to project on from x.
751 * Determine the right-hand side of the matrix equation using coordinates xp. */
752 static void gmx_simdcall calc_dr_x_xp_simd(int b0,
754 gmx::ArrayRef<const AtomPair> atoms,
755 const rvec* gmx_restrict x,
756 const rvec* gmx_restrict xp,
757 const real* gmx_restrict bllen,
758 const real* gmx_restrict blc,
759 const real* pbc_simd,
760 rvec* gmx_restrict r,
761 real* gmx_restrict rhs,
762 real* gmx_restrict sol)
764 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
765 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
767 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
772 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
774 SimdReal x0_S, y0_S, z0_S;
775 SimdReal x1_S, y1_S, z1_S;
776 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
777 SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
778 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
779 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
781 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
783 offset0[i] = atoms[bs + i].index1;
784 offset1[i] = atoms[bs + i].index2;
787 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset0, &x0_S, &y0_S, &z0_S);
788 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset1, &x1_S, &y1_S, &z1_S);
793 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
795 n2_S = norm2(rx_S, ry_S, rz_S);
796 il_S = invsqrt(n2_S);
802 transposeScatterStoreU<3>(reinterpret_cast<real*>(r + bs), offset2, rx_S, ry_S, rz_S);
804 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(xp), offset0, &x0_S, &y0_S, &z0_S);
805 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(xp), offset1, &x1_S, &y1_S, &z1_S);
811 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
813 ip_S = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
815 rhs_S = load<SimdReal>(blc + bs) * (ip_S - load<SimdReal>(bllen + bs));
817 store(rhs + bs, rhs_S);
818 store(sol + bs, rhs_S);
821 #endif // GMX_SIMD_HAVE_REAL
823 /*! \brief Determine the distances and right-hand side for the next iteration. */
824 gmx_unused static void calc_dist_iter(int b0,
826 gmx::ArrayRef<const AtomPair> atoms,
827 const rvec* gmx_restrict xp,
828 const real* gmx_restrict bllen,
829 const real* gmx_restrict blc,
832 real* gmx_restrict rhs,
833 real* gmx_restrict sol,
836 for (int b = b0; b < b1; b++)
842 pbc_dx_aiuc(pbc, xp[atoms[b].index1], xp[atoms[b].index2], dx);
846 rvec_sub(xp[atoms[b].index1], xp[atoms[b].index2], dx);
848 real len2 = len * len;
849 real dlen2 = 2 * len2 - ::norm2(dx);
850 if (dlen2 < wfac * len2)
852 /* not race free - see detailed comment in caller */
858 mvb = blc[b] * (len - dlen2 * gmx::invsqrt(dlen2));
866 } /* 20*ncons flops */
869 #if GMX_SIMD_HAVE_REAL
870 /*! \brief As calc_dist_iter(), but using SIMD intrinsics. */
871 static void gmx_simdcall calc_dist_iter_simd(int b0,
873 gmx::ArrayRef<const AtomPair> atoms,
874 const rvec* gmx_restrict x,
875 const real* gmx_restrict bllen,
876 const real* gmx_restrict blc,
877 const real* pbc_simd,
879 real* gmx_restrict rhs,
880 real* gmx_restrict sol,
883 SimdReal min_S(GMX_REAL_MIN);
885 SimdReal wfac_S(wfac);
888 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
890 /* Initialize all to FALSE */
891 warn_S = (two_S < setZero());
893 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
895 SimdReal x0_S, y0_S, z0_S;
896 SimdReal x1_S, y1_S, z1_S;
897 SimdReal rx_S, ry_S, rz_S, n2_S;
898 SimdReal len_S, len2_S, dlen2_S, lc_S, blc_S;
899 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
900 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
902 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
904 offset0[i] = atoms[bs + i].index1;
905 offset1[i] = atoms[bs + i].index2;
908 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset0, &x0_S, &y0_S, &z0_S);
909 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset1, &x1_S, &y1_S, &z1_S);
915 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
917 n2_S = norm2(rx_S, ry_S, rz_S);
919 len_S = load<SimdReal>(bllen + bs);
920 len2_S = len_S * len_S;
922 dlen2_S = fms(two_S, len2_S, n2_S);
924 warn_S = warn_S || (dlen2_S < (wfac_S * len2_S));
926 /* Avoid 1/0 by taking the max with REAL_MIN.
927 * Note: when dlen2 is close to zero (90 degree constraint rotation),
928 * the accuracy of the algorithm is no longer relevant.
930 dlen2_S = max(dlen2_S, min_S);
932 lc_S = fnma(dlen2_S, invsqrt(dlen2_S), len_S);
934 blc_S = load<SimdReal>(blc + bs);
938 store(rhs + bs, lc_S);
939 store(sol + bs, lc_S);
947 #endif // GMX_SIMD_HAVE_REAL
949 //! Implements LINCS constraining.
950 static void do_lincs(ArrayRefWithPadding<const RVec> xPadded,
951 ArrayRefWithPadding<RVec> xpPadded,
966 const rvec* x = as_rvec_array(xPadded.paddedArrayRef().data());
967 rvec* xp = as_rvec_array(xpPadded.paddedArrayRef().data());
968 rvec* gmx_restrict v = as_rvec_array(vRef.data());
970 const int b0 = lincsd->task[th].b0;
971 const int b1 = lincsd->task[th].b1;
973 gmx::ArrayRef<const AtomPair> atoms = lincsd->atoms;
974 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
975 gmx::ArrayRef<const int> blnr = lincsd->blnr;
976 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
977 gmx::ArrayRef<const real> blc = lincsd->blc;
978 gmx::ArrayRef<const real> blmf = lincsd->blmf;
979 gmx::ArrayRef<const real> bllen = lincsd->bllen;
980 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
981 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
982 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
983 gmx::ArrayRef<real> sol = lincsd->tmp3;
984 gmx::ArrayRef<real> blc_sol = lincsd->tmp4;
985 gmx::ArrayRef<real> mlambda = lincsd->mlambda;
986 gmx::ArrayRef<const int> nlocat = lincsd->nlocat;
988 #if GMX_SIMD_HAVE_REAL
990 /* This SIMD code does the same as the plain-C code after the #else.
991 * The only difference is that we always call pbc code, as with SIMD
992 * the overhead of pbc computation (when not needed) is small.
994 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9 * GMX_SIMD_REAL_WIDTH];
996 /* Convert the pbc struct for SIMD */
997 set_pbc_simd(pbc, pbc_simd);
999 /* Compute normalized x i-j vectors, store in r.
1000 * Compute the inner product of r and xp i-j and store in rhs1.
1002 calc_dr_x_xp_simd(b0, b1, atoms, x, xp, bllen.data(), blc.data(), pbc_simd,
1003 as_rvec_array(r.data()), rhs1.data(), sol.data());
1005 #else // GMX_SIMD_HAVE_REAL
1009 /* Compute normalized i-j vectors */
1010 for (int b = b0; b < b1; b++)
1013 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
1016 pbc_dx_aiuc(pbc, xp[atoms[b].index1], xp[atoms[b].index2], dx);
1017 real mvb = blc[b] * (::iprod(r[b], dx) - bllen[b]);
1024 /* Compute normalized i-j vectors */
1025 for (int b = b0; b < b1; b++)
1027 int i = atoms[b].index1;
1028 int j = atoms[b].index2;
1029 real tmp0 = x[i][0] - x[j][0];
1030 real tmp1 = x[i][1] - x[j][1];
1031 real tmp2 = x[i][2] - x[j][2];
1032 real rlen = gmx::invsqrt(tmp0 * tmp0 + tmp1 * tmp1 + tmp2 * tmp2);
1033 r[b][0] = rlen * tmp0;
1034 r[b][1] = rlen * tmp1;
1035 r[b][2] = rlen * tmp2;
1036 /* 16 ncons flops */
1039 * (r[b][0] * (xp[i][0] - xp[j][0]) + r[b][1] * (xp[i][1] - xp[j][1])
1040 + r[b][2] * (xp[i][2] - xp[j][2]) - bllen[b]);
1045 /* Together: 26*ncons + 6*nrtot flops */
1048 #endif // GMX_SIMD_HAVE_REAL
1050 if (lincsd->bTaskDep)
1052 /* We need a barrier, since the matrix construction below
1053 * can access entries in r of other threads.
1058 /* Construct the (sparse) LINCS matrix */
1059 for (int b = b0; b < b1; b++)
1061 for (int n = blnr[b]; n < blnr[b + 1]; n++)
1063 blcc[n] = blmf[n] * gmx::dot(r[b], r[blbnb[n]]);
1066 /* Together: 26*ncons + 6*nrtot flops */
1068 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
1069 /* nrec*(ncons+2*nrtot) flops */
1071 #if GMX_SIMD_HAVE_REAL
1072 for (int b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1074 SimdReal t1 = load<SimdReal>(blc.data() + b);
1075 SimdReal t2 = load<SimdReal>(sol.data() + b);
1076 store(mlambda.data() + b, t1 * t2);
1079 for (int b = b0; b < b1; b++)
1081 mlambda[b] = blc[b] * sol[b];
1083 #endif // GMX_SIMD_HAVE_REAL
1085 /* Update the coordinates */
1086 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1089 ******** Correction for centripetal effects ********
1094 wfac = std::cos(DEG2RAD * wangle);
1097 for (int iter = 0; iter < lincsd->nIter; iter++)
1099 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
1104 /* Communicate the corrected non-local coordinates */
1105 if (DOMAINDECOMP(cr))
1107 dd_move_x_constraints(cr->dd, box, xpPadded.unpaddedArrayRef(), ArrayRef<RVec>(), FALSE);
1112 else if (lincsd->bTaskDep)
1117 #if GMX_SIMD_HAVE_REAL
1118 calc_dist_iter_simd(b0, b1, atoms, xp, bllen.data(), blc.data(), pbc_simd, wfac,
1119 rhs1.data(), sol.data(), bWarn);
1121 calc_dist_iter(b0, b1, atoms, xp, bllen.data(), blc.data(), pbc, wfac, rhs1.data(),
1123 /* 20*ncons flops */
1124 #endif // GMX_SIMD_HAVE_REAL
1126 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
1127 /* nrec*(ncons+2*nrtot) flops */
1129 #if GMX_SIMD_HAVE_REAL
1130 for (int b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1132 SimdReal t1 = load<SimdReal>(blc.data() + b);
1133 SimdReal t2 = load<SimdReal>(sol.data() + b);
1134 SimdReal mvb = t1 * t2;
1135 store(blc_sol.data() + b, mvb);
1136 store(mlambda.data() + b, load<SimdReal>(mlambda.data() + b) + mvb);
1139 for (int b = b0; b < b1; b++)
1141 real mvb = blc[b] * sol[b];
1145 #endif // GMX_SIMD_HAVE_REAL
1147 /* Update the coordinates */
1148 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1150 /* nit*ncons*(37+9*nrec) flops */
1154 /* Update the velocities */
1155 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1156 /* 16 ncons flops */
1159 if (!nlocat.empty() && (bCalcDHDL || bCalcVir))
1161 if (lincsd->bTaskDep)
1163 /* In lincs_update_atoms threads might cross-read mlambda */
1167 /* Only account for local atoms */
1168 for (int b = b0; b < b1; b++)
1170 mlambda[b] *= 0.5 * nlocat[b];
1177 for (int b = b0; b < b1; b++)
1179 /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1180 * later after the contributions are reduced over the threads.
1182 dhdl -= lincsd->mlambda[b] * lincsd->ddist[b];
1184 lincsd->task[th].dhdlambda = dhdl;
1189 /* Constraint virial */
1190 for (int b = b0; b < b1; b++)
1192 real tmp0 = -bllen[b] * mlambda[b];
1193 for (int i = 0; i < DIM; i++)
1195 real tmp1 = tmp0 * r[b][i];
1196 for (int j = 0; j < DIM; j++)
1198 vir_r_m_dr[i][j] -= tmp1 * r[b][j];
1201 } /* 22 ncons flops */
1205 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1206 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1208 * (26+nrec)*ncons + (6+2*nrec)*nrtot
1209 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1211 * (63+nrec)*ncons + (6+4*nrec)*nrtot
1215 /*! \brief Sets the elements in the LINCS matrix for task task. */
1216 static void set_lincs_matrix_task(Lincs* li, Task* li_task, const real* invmass, int* ncc_triangle, int* nCrossTaskTriangles)
1218 /* Construct the coupling coefficient matrix blmf */
1219 li_task->ntriangle = 0;
1221 *nCrossTaskTriangles = 0;
1222 for (int i = li_task->b0; i < li_task->b1; i++)
1224 const int a1 = li->atoms[i].index1;
1225 const int a2 = li->atoms[i].index2;
1226 for (int n = li->blnr[i]; n < li->blnr[i + 1]; n++)
1228 const int k = li->blbnb[n];
1230 /* If we are using multiple, independent tasks for LINCS,
1231 * the calls to check_assign_connected should have
1232 * put all connected constraints in our task.
1234 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1237 if (a1 == li->atoms[k].index1 || a2 == li->atoms[k].index2)
1247 if (a1 == li->atoms[k].index1 || a1 == li->atoms[k].index2)
1257 li->blmf[n] = sign * invmass[center] * li->blc[i] * li->blc[k];
1258 li->blmf1[n] = sign * 0.5;
1259 if (li->ncg_triangle > 0)
1261 /* Look for constraint triangles */
1262 for (int nk = li->blnr[k]; nk < li->blnr[k + 1]; nk++)
1264 const int kk = li->blbnb[nk];
1265 if (kk != i && kk != k && (li->atoms[kk].index1 == end || li->atoms[kk].index2 == end))
1267 /* Check if the constraints in this triangle actually
1268 * belong to a different task. We still assign them
1269 * here, since it's convenient for the triangle
1270 * iterations, but we then need an extra barrier.
1272 if (k < li_task->b0 || k >= li_task->b1 || kk < li_task->b0 || kk >= li_task->b1)
1274 (*nCrossTaskTriangles)++;
1277 if (li_task->ntriangle == 0 || li_task->triangle[li_task->ntriangle - 1] < i)
1279 /* Add this constraint to the triangle list */
1280 li_task->triangle[li_task->ntriangle] = i;
1281 li_task->tri_bits[li_task->ntriangle] = 0;
1282 li_task->ntriangle++;
1283 if (li->blnr[i + 1] - li->blnr[i]
1284 > static_cast<int>(sizeof(li_task->tri_bits[0]) * 8 - 1))
1287 "A constraint is connected to %d constraints, this is "
1288 "more than the %zu allowed for constraints participating "
1290 li->blnr[i + 1] - li->blnr[i],
1291 sizeof(li_task->tri_bits[0]) * 8 - 1);
1294 li_task->tri_bits[li_task->ntriangle - 1] |= (1 << (n - li->blnr[i]));
1303 /*! \brief Sets the elements in the LINCS matrix. */
1304 static void set_lincs_matrix(Lincs* li, const real* invmass, real lambda)
1306 const real invsqrt2 = 0.7071067811865475244;
1308 for (int i = 0; (i < li->nc); i++)
1310 const int a1 = li->atoms[i].index1;
1311 const int a2 = li->atoms[i].index2;
1312 li->blc[i] = gmx::invsqrt(invmass[a1] + invmass[a2]);
1313 li->blc1[i] = invsqrt2;
1316 /* Construct the coupling coefficient matrix blmf */
1317 int ntriangle = 0, ncc_triangle = 0, nCrossTaskTriangles = 0;
1318 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle, nCrossTaskTriangles) num_threads(li->ntask) schedule(static)
1319 for (int th = 0; th < li->ntask; th++)
1323 set_lincs_matrix_task(li, &li->task[th], invmass, &ncc_triangle, &nCrossTaskTriangles);
1324 ntriangle += li->task[th].ntriangle;
1326 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1328 li->ntriangle = ntriangle;
1329 li->ncc_triangle = ncc_triangle;
1330 li->bTaskDepTri = (nCrossTaskTriangles > 0);
1334 fprintf(debug, "The %d constraints participate in %d triangles\n", li->nc, li->ntriangle);
1335 fprintf(debug, "There are %d constraint couplings, of which %d in triangles\n", li->ncc,
1337 if (li->ntriangle > 0 && li->ntask > 1)
1340 "%d constraint triangles contain constraints assigned to different tasks\n",
1341 nCrossTaskTriangles);
1346 * so we know with which lambda value the masses have been set.
1348 li->matlam = lambda;
1351 //! Finds all triangles of atoms that share constraints to a central atom.
1352 static int count_triangle_constraints(const InteractionLists& ilist, const ListOfLists<int>& at2con)
1354 const int ncon1 = ilist[F_CONSTR].size() / 3;
1355 const int ncon_tot = ncon1 + ilist[F_CONSTRNC].size() / 3;
1357 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1358 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1360 int ncon_triangle = 0;
1361 for (int c0 = 0; c0 < ncon_tot; c0++)
1363 bool bTriangle = FALSE;
1364 const int* iap = constr_iatomptr(ia1, ia2, c0);
1365 const int a00 = iap[1];
1366 const int a01 = iap[2];
1367 for (const int c1 : at2con[a01])
1371 const int* iap = constr_iatomptr(ia1, ia2, c1);
1372 const int a10 = iap[1];
1373 const int a11 = iap[2];
1383 for (const int c2 : at2con[ac1])
1385 if (c2 != c0 && c2 != c1)
1387 const int* iap = constr_iatomptr(ia1, ia2, c2);
1388 const int a20 = iap[1];
1389 const int a21 = iap[2];
1390 if (a20 == a00 || a21 == a00)
1404 return ncon_triangle;
1407 //! Finds sequences of sequential constraints.
1408 static bool more_than_two_sequential_constraints(const InteractionLists& ilist, const ListOfLists<int>& at2con)
1410 const int ncon1 = ilist[F_CONSTR].size() / 3;
1411 const int ncon_tot = ncon1 + ilist[F_CONSTRNC].size() / 3;
1413 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1414 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1416 for (int c = 0; c < ncon_tot; c++)
1418 const int* iap = constr_iatomptr(ia1, ia2, c);
1419 const int a1 = iap[1];
1420 const int a2 = iap[2];
1421 /* Check if this constraint has constraints connected at both atoms */
1422 if (at2con[a1].ssize() > 1 && at2con[a2].ssize() > 1)
1431 Lincs* init_lincs(FILE* fplog,
1432 const gmx_mtop_t& mtop,
1433 int nflexcon_global,
1434 ArrayRef<const ListOfLists<int>> atomToConstraintsPerMolType,
1439 // TODO this should become a unique_ptr
1441 bool bMoreThanTwoSeq;
1445 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n", bPLINCS ? " Parallel" : "");
1450 li->ncg = gmx_mtop_ftype_count(mtop, F_CONSTR) + gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1451 li->ncg_flex = nflexcon_global;
1454 li->nOrder = nProjOrder;
1456 li->max_connect = 0;
1457 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1459 const auto& at2con = atomToConstraintsPerMolType[mt];
1460 for (int a = 0; a < mtop.moltype[mt].atoms.nr; a++)
1462 li->max_connect = std::max(li->max_connect, int(at2con[a].ssize()));
1466 li->ncg_triangle = 0;
1467 bMoreThanTwoSeq = FALSE;
1468 for (const gmx_molblock_t& molb : mtop.molblock)
1470 const gmx_moltype_t& molt = mtop.moltype[molb.type];
1471 const auto& at2con = atomToConstraintsPerMolType[molb.type];
1473 li->ncg_triangle += molb.nmol * count_triangle_constraints(molt.ilist, at2con);
1475 if (!bMoreThanTwoSeq && more_than_two_sequential_constraints(molt.ilist, at2con))
1477 bMoreThanTwoSeq = TRUE;
1481 /* Check if we need to communicate not only before LINCS,
1482 * but also before each iteration.
1483 * The check for only two sequential constraints is only
1484 * useful for the common case of H-bond only constraints.
1485 * With more effort we could also make it useful for small
1486 * molecules with nr. sequential constraints <= nOrder-1.
1488 li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1490 if (debug && bPLINCS)
1492 fprintf(debug, "PLINCS communication before each iteration: %d\n", static_cast<int>(li->bCommIter));
1495 /* LINCS can run on any number of threads.
1496 * Currently the number is fixed for the whole simulation,
1497 * but it could be set in set_lincs().
1498 * The current constraint to task assignment code can create independent
1499 * tasks only when not more than two constraints are connected sequentially.
1501 li->ntask = gmx_omp_nthreads_get(emntLINCS);
1502 li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1505 fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n", li->ntask,
1506 li->bTaskDep ? "" : "in");
1514 /* Allocate an extra elements for "task-overlap" constraints */
1515 li->task.resize(li->ntask + 1);
1518 if (bPLINCS || li->ncg_triangle > 0)
1520 please_cite(fplog, "Hess2008a");
1524 please_cite(fplog, "Hess97a");
1529 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1533 "There are constraints between atoms in different decomposition domains,\n"
1534 "will communicate selected coordinates each lincs iteration\n");
1536 if (li->ncg_triangle > 0)
1539 "%d constraints are involved in constraint triangles,\n"
1540 "will apply an additional matrix expansion of order %d for couplings\n"
1541 "between constraints inside triangles\n",
1542 li->ncg_triangle, li->nOrder);
1549 void done_lincs(Lincs* li)
1554 /*! \brief Sets up the work division over the threads. */
1555 static void lincs_thread_setup(Lincs* li, int natoms)
1557 li->atf.resize(natoms);
1559 gmx::ArrayRef<gmx_bitmask_t> atf = li->atf;
1561 /* Clear the atom flags */
1562 for (gmx_bitmask_t& mask : atf)
1564 bitmask_clear(&mask);
1567 if (li->ntask > BITMASK_SIZE)
1569 gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1572 for (int th = 0; th < li->ntask; th++)
1574 const Task* li_task = &li->task[th];
1576 /* For each atom set a flag for constraints from each */
1577 for (int b = li_task->b0; b < li_task->b1; b++)
1579 bitmask_set_bit(&atf[li->atoms[b].index1], th);
1580 bitmask_set_bit(&atf[li->atoms[b].index2], th);
1584 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1585 for (int th = 0; th < li->ntask; th++)
1593 li_task = &li->task[th];
1595 bitmask_init_low_bits(&mask, th);
1597 li_task->ind.clear();
1598 li_task->ind_r.clear();
1599 for (b = li_task->b0; b < li_task->b1; b++)
1601 /* We let the constraint with the lowest thread index
1602 * operate on atoms with constraints from multiple threads.
1604 if (bitmask_is_disjoint(atf[li->atoms[b].index1], mask)
1605 && bitmask_is_disjoint(atf[li->atoms[b].index2], mask))
1607 /* Add the constraint to the local atom update index */
1608 li_task->ind.push_back(b);
1612 /* Add the constraint to the rest block */
1613 li_task->ind_r.push_back(b);
1617 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1620 /* We need to copy all constraints which have not be assigned
1621 * to a thread to a separate list which will be handled by one thread.
1623 Task* li_m = &li->task[li->ntask];
1626 for (int th = 0; th < li->ntask; th++)
1628 const Task& li_task = li->task[th];
1630 for (int ind_r : li_task.ind_r)
1632 li_m->ind.push_back(ind_r);
1637 fprintf(debug, "LINCS thread %d: %zu constraints\n", th, li_task.ind.size());
1643 fprintf(debug, "LINCS thread r: %zu constraints\n", li_m->ind.size());
1647 //! Assign a constraint.
1648 static void assign_constraint(Lincs* li,
1649 int constraint_index,
1654 const ListOfLists<int>& at2con)
1660 /* Make an mapping of local topology constraint index to LINCS index */
1661 li->con_index[constraint_index] = con;
1663 li->bllen0[con] = lenA;
1664 li->ddist[con] = lenB - lenA;
1665 /* Set the length to the topology A length */
1666 li->bllen[con] = lenA;
1667 li->atoms[con].index1 = a1;
1668 li->atoms[con].index2 = a2;
1670 /* Make space in the constraint connection matrix for constraints
1671 * connected to both end of the current constraint.
1673 li->ncc += at2con[a1].ssize() - 1 + at2con[a2].ssize() - 1;
1675 li->blnr[con + 1] = li->ncc;
1677 /* Increase the constraint count */
1681 /*! \brief Check if constraint with topology index constraint_index is connected
1682 * to other constraints, and if so add those connected constraints to our task. */
1683 static void check_assign_connected(Lincs* li,
1684 gmx::ArrayRef<const int> iatom,
1685 const InteractionDefinitions& idef,
1689 const ListOfLists<int>& at2con)
1691 /* Currently this function only supports constraint groups
1692 * in which all constraints share at least one atom
1693 * (e.g. H-bond constraints).
1694 * Check both ends of the current constraint for
1695 * connected constraints. We need to assign those
1698 for (int end = 0; end < 2; end++)
1700 const int a = (end == 0 ? a1 : a2);
1702 for (const int cc : at2con[a])
1704 /* Check if constraint cc has not yet been assigned */
1705 if (li->con_index[cc] == -1)
1707 const int type = iatom[cc * 3];
1708 const real lenA = idef.iparams[type].constr.dA;
1709 const real lenB = idef.iparams[type].constr.dB;
1711 if (bDynamics || lenA != 0 || lenB != 0)
1713 assign_constraint(li, cc, iatom[3 * cc + 1], iatom[3 * cc + 2], lenA, lenB, at2con);
1720 /*! \brief Check if constraint with topology index constraint_index is involved
1721 * in a constraint triangle, and if so add the other two constraints
1722 * in the triangle to our task. */
1723 static void check_assign_triangle(Lincs* li,
1724 gmx::ArrayRef<const int> iatom,
1725 const InteractionDefinitions& idef,
1727 int constraint_index,
1730 const ListOfLists<int>& at2con)
1732 int nca, cc[32], ca[32];
1733 int c_triangle[2] = { -1, -1 };
1736 for (const int c : at2con[a1])
1738 if (c != constraint_index)
1742 aa1 = iatom[c * 3 + 1];
1743 aa2 = iatom[c * 3 + 2];
1759 for (const int c : at2con[a2])
1761 if (c != constraint_index)
1765 aa1 = iatom[c * 3 + 1];
1766 aa2 = iatom[c * 3 + 2];
1769 for (i = 0; i < nca; i++)
1773 c_triangle[0] = cc[i];
1780 for (i = 0; i < nca; i++)
1784 c_triangle[0] = cc[i];
1792 if (c_triangle[0] >= 0)
1796 for (end = 0; end < 2; end++)
1798 /* Check if constraint c_triangle[end] has not yet been assigned */
1799 if (li->con_index[c_triangle[end]] == -1)
1804 i = c_triangle[end] * 3;
1806 lenA = idef.iparams[type].constr.dA;
1807 lenB = idef.iparams[type].constr.dB;
1809 if (bDynamics || lenA != 0 || lenB != 0)
1811 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1818 //! Sets matrix indices.
1819 static void set_matrix_indices(Lincs* li, const Task& li_task, const ListOfLists<int>& at2con, bool bSortMatrix)
1821 for (int b = li_task.b0; b < li_task.b1; b++)
1823 const int a1 = li->atoms[b].index1;
1824 const int a2 = li->atoms[b].index2;
1826 int i = li->blnr[b];
1827 for (const int constraint : at2con[a1])
1829 const int concon = li->con_index[constraint];
1832 li->blbnb[i++] = concon;
1835 for (const int constraint : at2con[a2])
1837 const int concon = li->con_index[constraint];
1840 li->blbnb[i++] = concon;
1846 /* Order the blbnb matrix to optimize memory access */
1847 std::sort(li->blbnb.begin() + li->blnr[b], li->blbnb.begin() + li->blnr[b + 1]);
1852 void set_lincs(const InteractionDefinitions& idef,
1854 const real* invmass,
1857 const t_commrec* cr,
1863 /* Zero the thread index ranges.
1864 * Otherwise without local constraints we could return with old ranges.
1866 for (int i = 0; i < li->ntask; i++)
1870 li->task[i].ind.clear();
1874 li->task[li->ntask].ind.clear();
1877 /* This is the local topology, so there are only F_CONSTR constraints */
1878 if (idef.il[F_CONSTR].empty())
1880 /* There are no constraints,
1881 * we do not need to fill any data structures.
1888 fprintf(debug, "Building the LINCS connectivity\n");
1892 if (DOMAINDECOMP(cr))
1894 if (cr->dd->constraints)
1898 dd_get_constraint_range(cr->dd, &start, &natoms);
1902 natoms = dd_numHomeAtoms(*cr->dd);
1910 const ListOfLists<int> at2con =
1911 make_at2con(natoms, idef.il, idef.iparams, flexibleConstraintTreatment(bDynamics));
1913 const int ncon_tot = idef.il[F_CONSTR].size() / 3;
1915 /* Ensure we have enough padding for aligned loads for each thread */
1916 const int numEntries = ncon_tot + li->ntask * simd_width;
1917 li->con_index.resize(numEntries);
1918 li->bllen0.resize(numEntries);
1919 li->ddist.resize(numEntries);
1920 li->atoms.resize(numEntries);
1921 li->blc.resize(numEntries);
1922 li->blc1.resize(numEntries);
1923 li->blnr.resize(numEntries + 1);
1924 li->bllen.resize(numEntries);
1925 li->tmpv.resizeWithPadding(numEntries);
1926 if (DOMAINDECOMP(cr))
1928 li->nlocat.resize(numEntries);
1930 li->tmp1.resize(numEntries);
1931 li->tmp2.resize(numEntries);
1932 li->tmp3.resize(numEntries);
1933 li->tmp4.resize(numEntries);
1934 li->mlambda.resize(numEntries);
1936 gmx::ArrayRef<const int> iatom = idef.il[F_CONSTR].iatoms;
1938 li->blnr[0] = li->ncc;
1940 /* Assign the constraints for li->ntask LINCS tasks.
1941 * We target a uniform distribution of constraints over the tasks.
1942 * Note that when flexible constraints are present, but are removed here
1943 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
1944 * happen during normal MD, that's ok.
1947 /* Determine the number of constraints we need to assign here */
1948 int ncon_assign = ncon_tot;
1951 /* With energy minimization, flexible constraints are ignored
1952 * (and thus minimized, as they should be).
1954 ncon_assign -= countFlexibleConstraints(idef.il, idef.iparams);
1957 /* Set the target constraint count per task to exactly uniform,
1958 * this might be overridden below.
1960 int ncon_target = (ncon_assign + li->ntask - 1) / li->ntask;
1962 /* Mark all constraints as unassigned by setting their index to -1 */
1963 for (int con = 0; con < ncon_tot; con++)
1965 li->con_index[con] = -1;
1969 for (int th = 0; th < li->ntask; th++)
1973 li_task = &li->task[th];
1975 #if GMX_SIMD_HAVE_REAL
1976 /* With indepedent tasks we likely have H-bond constraints or constraint
1977 * pairs. The connected constraints will be pulled into the task, so the
1978 * constraints per task will often exceed ncon_target.
1979 * Triangle constraints can also increase the count, but there are
1980 * relatively few of those, so we usually expect to get ncon_target.
1984 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
1985 * since otherwise a lot of operations can be wasted.
1986 * There are several ways to round here, we choose the one
1987 * that alternates block sizes, which helps with Intel HT.
1989 ncon_target = ((ncon_assign * (th + 1)) / li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1)
1990 & ~(GMX_SIMD_REAL_WIDTH - 1);
1992 #endif // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
1994 /* Continue filling the arrays where we left off with the previous task,
1995 * including padding for SIMD.
1997 li_task->b0 = li->nc;
1999 gmx::ArrayRef<const t_iparams> iparams = idef.iparams;
2001 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2003 if (li->con_index[con] == -1)
2008 type = iatom[3 * con];
2009 a1 = iatom[3 * con + 1];
2010 a2 = iatom[3 * con + 2];
2011 lenA = iparams[type].constr.dA;
2012 lenB = iparams[type].constr.dB;
2013 /* Skip the flexible constraints when not doing dynamics */
2014 if (bDynamics || lenA != 0 || lenB != 0)
2016 assign_constraint(li, con, a1, a2, lenA, lenB, at2con);
2018 if (li->ntask > 1 && !li->bTaskDep)
2020 /* We can generate independent tasks. Check if we
2021 * need to assign connected constraints to our task.
2023 check_assign_connected(li, iatom, idef, bDynamics, a1, a2, at2con);
2025 if (li->ntask > 1 && li->ncg_triangle > 0)
2027 /* Ensure constraints in one triangle are assigned
2030 check_assign_triangle(li, iatom, idef, bDynamics, con, a1, a2, at2con);
2038 li_task->b1 = li->nc;
2042 /* Copy the last atom pair indices and lengths for constraints
2043 * up to a multiple of simd_width, such that we can do all
2044 * SIMD operations without having to worry about end effects.
2048 li->nc = ((li_task->b1 + simd_width - 1) / simd_width) * simd_width;
2049 last = li_task->b1 - 1;
2050 for (i = li_task->b1; i < li->nc; i++)
2052 li->atoms[i] = li->atoms[last];
2053 li->bllen0[i] = li->bllen0[last];
2054 li->ddist[i] = li->ddist[last];
2055 li->bllen[i] = li->bllen[last];
2056 li->blnr[i + 1] = li->blnr[last + 1];
2060 /* Keep track of how many constraints we assigned */
2061 li->nc_real += li_task->b1 - li_task->b0;
2065 fprintf(debug, "LINCS task %d constraints %d - %d\n", th, li_task->b0, li_task->b1);
2069 assert(li->nc_real == ncon_assign);
2073 /* Without DD we order the blbnb matrix to optimize memory access.
2074 * With DD the overhead of sorting is more than the gain during access.
2076 bSortMatrix = !DOMAINDECOMP(cr);
2078 li->blbnb.resize(li->ncc);
2080 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2081 for (int th = 0; th < li->ntask; th++)
2085 Task& li_task = li->task[th];
2087 if (li->ncg_triangle > 0)
2089 /* This is allocating too much, but it is difficult to improve */
2090 li_task.triangle.resize(li_task.b1 - li_task.b0);
2091 li_task.tri_bits.resize(li_task.b1 - li_task.b0);
2094 set_matrix_indices(li, li_task, at2con, bSortMatrix);
2096 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2099 if (cr->dd == nullptr)
2101 /* Since the matrix is static, we should free some memory */
2102 li->blbnb.resize(li->ncc);
2105 li->blmf.resize(li->ncc);
2106 li->blmf1.resize(li->ncc);
2107 li->tmpncc.resize(li->ncc);
2109 gmx::ArrayRef<const int> nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2110 if (!nlocat_dd.empty())
2112 /* Convert nlocat from local topology to LINCS constraint indexing */
2113 for (con = 0; con < ncon_tot; con++)
2115 li->nlocat[li->con_index[con]] = nlocat_dd[con];
2125 fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n", li->nc_real,
2131 lincs_thread_setup(li, numAtoms);
2134 set_lincs_matrix(li, invmass, lambda);
2137 //! Issues a warning when LINCS constraints cannot be satisfied.
2138 static void lincs_warning(gmx_domdec_t* dd,
2139 ArrayRef<const RVec> x,
2140 ArrayRef<const RVec> xprime,
2143 gmx::ArrayRef<const AtomPair> atoms,
2144 gmx::ArrayRef<const real> bllen,
2149 real wfac = std::cos(DEG2RAD * wangle);
2152 "bonds that rotated more than %g degrees:\n"
2153 " atom 1 atom 2 angle previous, current, constraint length\n",
2156 for (int b = 0; b < ncons; b++)
2158 const int i = atoms[b].index1;
2159 const int j = atoms[b].index2;
2164 pbc_dx_aiuc(pbc, x[i], x[j], v0);
2165 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2169 rvec_sub(x[i], x[j], v0);
2170 rvec_sub(xprime[i], xprime[j], v1);
2174 real cosine = ::iprod(v0, v1) / (d0 * d1);
2177 fprintf(stderr, " %6d %6d %5.1f %8.4f %8.4f %8.4f\n", ddglatnr(dd, i),
2178 ddglatnr(dd, j), RAD2DEG * std::acos(cosine), d0, d1, bllen[b]);
2179 if (!std::isfinite(d1))
2181 gmx_fatal(FARGS, "Bond length not finite.");
2187 if (*warncount > maxwarn)
2189 too_many_constraint_warnings(econtLINCS, *warncount);
2193 //! Status information about how well LINCS satisified the constraints in this domain
2194 struct LincsDeviations
2196 //! The maximum over all bonds in this domain of the relative deviation in bond lengths
2197 real maxDeviation = 0;
2198 //! Sum over all bonds in this domain of the squared relative deviation
2199 real sumSquaredDeviation = 0;
2200 //! Index of bond with max deviation
2201 int indexOfMaxDeviation = -1;
2202 //! Number of bonds constrained in this domain
2203 int numConstraints = 0;
2206 //! Determine how well the constraints have been satisfied.
2207 static LincsDeviations makeLincsDeviations(const Lincs& lincsd, ArrayRef<const RVec> x, const t_pbc* pbc)
2209 LincsDeviations result;
2210 const ArrayRef<const AtomPair> atoms = lincsd.atoms;
2211 const ArrayRef<const real> bllen = lincsd.bllen;
2212 const ArrayRef<const int> nlocat = lincsd.nlocat;
2214 for (int task = 0; task < lincsd.ntask; task++)
2216 for (int b = lincsd.task[task].b0; b < lincsd.task[task].b1; b++)
2221 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
2225 rvec_sub(x[atoms[b].index1], x[atoms[b].index2], dx);
2227 real r2 = ::norm2(dx);
2228 real len = r2 * gmx::invsqrt(r2);
2229 real d = std::abs(len / bllen[b] - 1.0_real);
2230 if (d > result.maxDeviation && (nlocat.empty() || nlocat[b]))
2232 result.maxDeviation = d;
2233 result.indexOfMaxDeviation = b;
2237 result.sumSquaredDeviation += d * d;
2238 result.numConstraints++;
2242 result.sumSquaredDeviation += nlocat[b] * d * d;
2243 result.numConstraints += nlocat[b];
2248 if (!nlocat.empty())
2250 result.numConstraints /= 2;
2251 result.sumSquaredDeviation *= 0.5;
2256 bool constrain_lincs(bool computeRmsd,
2257 const t_inputrec& ir,
2260 const real* invmass,
2261 const t_commrec* cr,
2262 const gmx_multisim_t* ms,
2263 ArrayRefWithPadding<const RVec> xPadded,
2264 ArrayRefWithPadding<RVec> xprimePadded,
2265 ArrayRef<RVec> min_proj,
2268 const bool hasMassPerturbed,
2275 ConstraintVariable econq,
2282 /* This boolean should be set by a flag passed to this routine.
2283 * We can also easily check if any constraint length is changed,
2284 * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2286 bool bCalcDHDL = (ir.efep != efepNO && dvdlambda != nullptr);
2288 if (lincsd->nc == 0 && cr->dd == nullptr)
2292 lincsd->rmsdData = { { 0 } };
2298 ArrayRef<const RVec> x = xPadded.unpaddedArrayRef();
2299 ArrayRef<RVec> xprime = xprimePadded.unpaddedArrayRef();
2301 if (econq == ConstraintVariable::Positions)
2303 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2304 * also with efep!=fepNO.
2306 if (ir.efep != efepNO)
2308 if (hasMassPerturbed && lincsd->matlam != lambda)
2310 set_lincs_matrix(lincsd, invmass, lambda);
2313 for (int i = 0; i < lincsd->nc; i++)
2315 lincsd->bllen[i] = lincsd->bllen0[i] + lambda * lincsd->ddist[i];
2319 if (lincsd->ncg_flex)
2321 /* Set the flexible constraint lengths to the old lengths */
2324 for (int i = 0; i < lincsd->nc; i++)
2326 if (lincsd->bllen[i] == 0)
2329 pbc_dx_aiuc(pbc, x[lincsd->atoms[i].index1], x[lincsd->atoms[i].index2], dx);
2330 lincsd->bllen[i] = norm(dx);
2336 for (int i = 0; i < lincsd->nc; i++)
2338 if (lincsd->bllen[i] == 0)
2340 lincsd->bllen[i] = std::sqrt(
2341 distance2(x[lincsd->atoms[i].index1], x[lincsd->atoms[i].index2]));
2347 const bool printDebugOutput = ((debug != nullptr) && lincsd->nc > 0);
2348 if (printDebugOutput)
2350 LincsDeviations deviations = makeLincsDeviations(*lincsd, xprime, pbc);
2351 fprintf(debug, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2352 fprintf(debug, " Before LINCS %.6f %.6f %6d %6d\n",
2353 std::sqrt(deviations.sumSquaredDeviation / deviations.numConstraints),
2354 deviations.maxDeviation,
2355 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index1),
2356 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index2));
2359 /* This bWarn var can be updated by multiple threads
2360 * at the same time. But as we only need to detect
2361 * if a warning occurred or not, this is not an issue.
2365 /* The OpenMP parallel region of constrain_lincs for coords */
2366 #pragma omp parallel num_threads(lincsd->ntask)
2370 int th = gmx_omp_get_thread_num();
2372 clear_mat(lincsd->task[th].vir_r_m_dr);
2374 do_lincs(xPadded, xprimePadded, box, pbc, lincsd, th, invmass, cr, bCalcDHDL,
2375 ir.LincsWarnAngle, &bWarn, invdt, v, bCalcVir,
2376 th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2378 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2381 if (computeRmsd || printDebugOutput || bWarn)
2383 LincsDeviations deviations = makeLincsDeviations(*lincsd, xprime, pbc);
2387 // This is reduced across domains in compute_globals and
2388 // reported to the log file.
2389 lincsd->rmsdData[0] = deviations.numConstraints;
2390 lincsd->rmsdData[1] = deviations.sumSquaredDeviation;
2394 // This is never read
2395 lincsd->rmsdData = { { 0 } };
2397 if (printDebugOutput)
2399 fprintf(debug, " After LINCS %.6f %.6f %6d %6d\n\n",
2400 std::sqrt(deviations.sumSquaredDeviation / deviations.numConstraints),
2401 deviations.maxDeviation,
2402 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index1),
2403 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index2));
2408 if (maxwarn < INT_MAX)
2410 std::string simMesg;
2413 simMesg += gmx::formatString(" in simulation %d", ms->sim);
2417 ", time %g (ps) LINCS WARNING%s\n"
2418 "relative constraint deviation after LINCS:\n"
2419 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2420 step, ir.init_t + step * ir.delta_t, simMesg.c_str(),
2421 std::sqrt(deviations.sumSquaredDeviation / deviations.numConstraints),
2422 deviations.maxDeviation,
2423 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index1),
2424 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index2));
2426 lincs_warning(cr->dd, x, xprime, pbc, lincsd->nc, lincsd->atoms, lincsd->bllen,
2427 ir.LincsWarnAngle, maxwarn, warncount);
2429 bOK = (deviations.maxDeviation < 0.5);
2433 if (lincsd->ncg_flex)
2435 for (int i = 0; (i < lincsd->nc); i++)
2437 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2439 lincsd->bllen[i] = 0;
2446 /* The OpenMP parallel region of constrain_lincs for derivatives */
2447 #pragma omp parallel num_threads(lincsd->ntask)
2451 int th = gmx_omp_get_thread_num();
2453 do_lincsp(xPadded, xprimePadded, min_proj, pbc, lincsd, th, invmass, econq,
2454 bCalcDHDL, bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2456 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2462 /* Reduce the dH/dlambda contributions over the threads */
2467 for (th = 0; th < lincsd->ntask; th++)
2469 dhdlambda += lincsd->task[th].dhdlambda;
2471 if (econq == ConstraintVariable::Positions)
2473 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2474 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2475 dhdlambda /= ir.delta_t * ir.delta_t;
2477 *dvdlambda += dhdlambda;
2480 if (bCalcVir && lincsd->ntask > 1)
2482 for (int i = 1; i < lincsd->ntask; i++)
2484 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2488 /* count assuming nit=1 */
2489 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2490 inc_nrnb(nrnb, eNR_LINCSMAT, (2 + lincsd->nOrder) * lincsd->ncc);
2491 if (lincsd->ntriangle > 0)
2493 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder * lincsd->ncc_triangle);
2497 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real * 2);
2501 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);