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/mdtypes/mdatom.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/pbcutil/pbc_simd.h"
74 #include "gromacs/simd/simd.h"
75 #include "gromacs/simd/simd_math.h"
76 #include "gromacs/simd/vector_operations.h"
77 #include "gromacs/topology/mtop_util.h"
78 #include "gromacs/utility/alignedallocator.h"
79 #include "gromacs/utility/arrayref.h"
80 #include "gromacs/utility/basedefinitions.h"
81 #include "gromacs/utility/bitmask.h"
82 #include "gromacs/utility/cstringutil.h"
83 #include "gromacs/utility/exceptions.h"
84 #include "gromacs/utility/fatalerror.h"
85 #include "gromacs/utility/gmxomp.h"
86 #include "gromacs/utility/listoflists.h"
87 #include "gromacs/utility/pleasecite.h"
89 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
94 //! Indices of the two atoms involved in a single constraint
97 //! \brief Constructor, does not initialize to catch bugs and faster construction
106 //! Unit of work within LINCS.
109 //! First constraint for this task.
111 //! b1-1 is the last constraint for this task.
113 //! The number of constraints in triangles.
115 //! The list of triangle constraints.
116 std::vector<int> triangle;
117 //! The bits tell if the matrix element should be used.
118 std::vector<int> tri_bits;
119 //! Constraint index for updating atom data.
120 std::vector<int> ind;
121 //! Constraint index for updating atom data.
122 std::vector<int> ind_r;
123 //! Temporary variable for virial calculation.
124 tensor vir_r_m_dr = { { 0 } };
125 //! Temporary variable for lambda derivative.
129 /*! \brief Data for LINCS algorithm.
134 //! The global number of constraints.
136 //! The global number of flexible constraints.
138 //! The global number of constraints in triangles.
139 int ncg_triangle = 0;
140 //! The number of iterations.
142 //! The order of the matrix expansion.
144 //! The maximum number of constraints connected to a single atom.
147 //! The number of real constraints.
149 //! The number of constraints including padding for SIMD.
151 //! The number of constraint connections.
153 //! The FE lambda value used for filling blc and blmf.
155 //! mapping from topology to LINCS constraints.
156 std::vector<int> con_index;
157 //! The reference distance in topology A.
158 std::vector<real, AlignedAllocator<real>> bllen0;
159 //! The reference distance in top B - the r.d. in top A.
160 std::vector<real, AlignedAllocator<real>> ddist;
161 //! The atom pairs involved in the constraints.
162 std::vector<AtomPair> atoms;
163 //! 1/sqrt(invmass1 invmass2).
164 std::vector<real, AlignedAllocator<real>> blc;
165 //! As blc, but with all masses 1.
166 std::vector<real, AlignedAllocator<real>> blc1;
167 //! Index into blbnb and blmf.
168 std::vector<int> blnr;
169 //! List of constraint connections.
170 std::vector<int> blbnb;
171 //! The local number of constraints in triangles.
173 //! The number of constraint connections in triangles.
174 int ncc_triangle = 0;
175 //! Communicate before each LINCS interation.
176 bool bCommIter = false;
177 //! Matrix of mass factors for constraint connections.
178 std::vector<real> blmf;
179 //! As blmf, but with all masses 1.
180 std::vector<real> blmf1;
181 //! The reference bond length.
182 std::vector<real, AlignedAllocator<real>> bllen;
183 //! The local atom count per constraint, can be NULL.
184 std::vector<int> nlocat;
186 /*! \brief The number of tasks used for LINCS work.
188 * \todo This is mostly used to loop over \c task, which would
189 * be nicer to do with range-based for loops, but the thread
190 * index is used for constructing bit masks and organizing the
191 * virial output buffer, so other things need to change,
194 /*! \brief LINCS thread division */
195 std::vector<Task> task;
196 //! Atom flags for thread parallelization.
197 std::vector<gmx_bitmask_t> atf;
198 //! Are the LINCS tasks interdependent?
199 bool bTaskDep = false;
200 //! Are there triangle constraints that cross task borders?
201 bool bTaskDepTri = false;
202 //! Arrays for temporary storage in the LINCS algorithm.
204 PaddedVector<gmx::RVec> tmpv;
205 std::vector<real> tmpncc;
206 std::vector<real, AlignedAllocator<real>> tmp1;
207 std::vector<real, AlignedAllocator<real>> tmp2;
208 std::vector<real, AlignedAllocator<real>> tmp3;
209 std::vector<real, AlignedAllocator<real>> tmp4;
211 //! The Lagrange multipliers times -1.
212 std::vector<real, AlignedAllocator<real>> mlambda;
213 //! Storage for the constraint RMS relative deviation output.
214 std::array<real, 2> rmsdData = { { 0 } };
217 /*! \brief Define simd_width for memory allocation used for SIMD code */
218 #if GMX_SIMD_HAVE_REAL
219 static const int simd_width = GMX_SIMD_REAL_WIDTH;
221 static const int simd_width = 1;
224 ArrayRef<real> lincs_rmsdData(Lincs* lincsd)
226 return lincsd->rmsdData;
229 real lincs_rmsd(const Lincs* lincsd)
231 if (lincsd->rmsdData[0] > 0)
233 return std::sqrt(lincsd->rmsdData[1] / lincsd->rmsdData[0]);
241 /*! \brief Do a set of nrec LINCS matrix multiplications.
243 * This function will return with up to date thread-local
244 * constraint data, without an OpenMP barrier.
246 static void lincs_matrix_expand(const Lincs& lincsd,
248 gmx::ArrayRef<const real> blcc,
249 gmx::ArrayRef<real> rhs1,
250 gmx::ArrayRef<real> rhs2,
251 gmx::ArrayRef<real> sol)
253 gmx::ArrayRef<const int> blnr = lincsd.blnr;
254 gmx::ArrayRef<const int> blbnb = lincsd.blbnb;
256 const int b0 = li_task.b0;
257 const int b1 = li_task.b1;
258 const int nrec = lincsd.nOrder;
260 for (int rec = 0; rec < nrec; rec++)
266 for (int b = b0; b < b1; b++)
272 for (n = blnr[b]; n < blnr[b + 1]; n++)
274 mvb = mvb + blcc[n] * rhs1[blbnb[n]];
277 sol[b] = sol[b] + mvb;
280 std::swap(rhs1, rhs2);
281 } /* nrec*(ncons+2*nrtot) flops */
283 if (lincsd.ntriangle > 0)
285 /* Perform an extra nrec recursions for only the constraints
286 * involved in rigid triangles.
287 * In this way their accuracy should come close to those of the other
288 * constraints, since traingles of constraints can produce eigenvalues
289 * around 0.7, while the effective eigenvalue for bond constraints
290 * is around 0.4 (and 0.7*0.7=0.5).
295 /* We need a barrier here, since other threads might still be
296 * reading the contents of rhs1 and/o rhs2.
297 * We could avoid this barrier by introducing two extra rhs
298 * arrays for the triangle constraints only.
303 /* Constraints involved in a triangle are ensured to be in the same
304 * LINCS task. This means no barriers are required during the extra
305 * iterations for the triangle constraints.
307 gmx::ArrayRef<const int> triangle = li_task.triangle;
308 gmx::ArrayRef<const int> tri_bits = li_task.tri_bits;
310 for (int rec = 0; rec < nrec; rec++)
312 for (int tb = 0; tb < li_task.ntriangle; tb++)
314 int b, bits, nr0, nr1, n;
322 for (n = nr0; n < nr1; n++)
324 if (bits & (1 << (n - nr0)))
326 mvb = mvb + blcc[n] * rhs1[blbnb[n]];
330 sol[b] = sol[b] + mvb;
333 std::swap(rhs1, rhs2);
334 } /* nrec*(ntriangle + ncc_triangle*2) flops */
336 if (lincsd.bTaskDepTri)
338 /* The constraints triangles are decoupled from each other,
339 * but constraints in one triangle cross thread task borders.
340 * We could probably avoid this with more advanced setup code.
347 //! Update atomic coordinates when an index is not required.
348 static void lincs_update_atoms_noind(int ncons,
349 gmx::ArrayRef<const AtomPair> atoms,
351 gmx::ArrayRef<const real> fac,
352 gmx::ArrayRef<const gmx::RVec> r,
356 if (invmass != nullptr)
358 for (int b = 0; b < ncons; b++)
360 int i = atoms[b].index1;
361 int j = atoms[b].index2;
362 real mvb = preFactor * fac[b];
363 real im1 = invmass[i];
364 real im2 = invmass[j];
365 real tmp0 = r[b][0] * mvb;
366 real tmp1 = r[b][1] * mvb;
367 real tmp2 = r[b][2] * mvb;
368 x[i][0] -= tmp0 * im1;
369 x[i][1] -= tmp1 * im1;
370 x[i][2] -= tmp2 * im1;
371 x[j][0] += tmp0 * im2;
372 x[j][1] += tmp1 * im2;
373 x[j][2] += tmp2 * im2;
374 } /* 16 ncons flops */
378 for (int b = 0; b < ncons; b++)
380 int i = atoms[b].index1;
381 int j = atoms[b].index2;
382 real mvb = preFactor * fac[b];
383 real tmp0 = r[b][0] * mvb;
384 real tmp1 = r[b][1] * mvb;
385 real tmp2 = r[b][2] * mvb;
396 //! Update atomic coordinates when an index is required.
397 static void lincs_update_atoms_ind(gmx::ArrayRef<const int> ind,
398 gmx::ArrayRef<const AtomPair> atoms,
400 gmx::ArrayRef<const real> fac,
401 gmx::ArrayRef<const gmx::RVec> r,
405 if (invmass != nullptr)
409 int i = atoms[b].index1;
410 int j = atoms[b].index2;
411 real mvb = preFactor * fac[b];
412 real im1 = invmass[i];
413 real im2 = invmass[j];
414 real tmp0 = r[b][0] * mvb;
415 real tmp1 = r[b][1] * mvb;
416 real tmp2 = r[b][2] * mvb;
417 x[i][0] -= tmp0 * im1;
418 x[i][1] -= tmp1 * im1;
419 x[i][2] -= tmp2 * im1;
420 x[j][0] += tmp0 * im2;
421 x[j][1] += tmp1 * im2;
422 x[j][2] += tmp2 * im2;
423 } /* 16 ncons flops */
429 int i = atoms[b].index1;
430 int j = atoms[b].index2;
431 real mvb = preFactor * fac[b];
432 real tmp0 = r[b][0] * mvb;
433 real tmp1 = r[b][1] * mvb;
434 real tmp2 = r[b][2] * mvb;
441 } /* 16 ncons flops */
445 //! Update coordinates for atoms.
446 static void lincs_update_atoms(Lincs* li,
449 gmx::ArrayRef<const real> fac,
450 gmx::ArrayRef<const gmx::RVec> r,
456 /* Single thread, we simply update for all constraints */
457 lincs_update_atoms_noind(li->nc_real, li->atoms, preFactor, fac, r, invmass, x);
461 /* Update the atom vector components for our thread local
462 * constraints that only access our local atom range.
463 * This can be done without a barrier.
465 lincs_update_atoms_ind(li->task[th].ind, li->atoms, preFactor, fac, r, invmass, x);
467 if (!li->task[li->ntask].ind.empty())
469 /* Update the constraints that operate on atoms
470 * in multiple thread atom blocks on the master thread.
475 lincs_update_atoms_ind(li->task[li->ntask].ind, li->atoms, preFactor, fac, r, invmass, x);
481 #if GMX_SIMD_HAVE_REAL
482 //! Helper function so that we can run TSAN with SIMD support (where implemented).
484 static inline void gmx_simdcall gatherLoadUTransposeTSANSafe(const real* base,
485 const std::int32_t* offset,
490 # if (CMAKE_BUILD_TYPE == CMAKE_BUILD_TYPE_TSAN) && GMX_SIMD_X86_AVX2_256
491 // This function is only implemented in this case
492 gatherLoadUTransposeSafe<align>(base, offset, v0, v1, v2);
494 gatherLoadUTranspose<align>(base, offset, v0, v1, v2);
498 /*! \brief Calculate the constraint distance vectors r to project on from x.
500 * Determine the right-hand side of the matrix equation using quantity f.
501 * This function only differs from calc_dr_x_xp_simd below in that
502 * no constraint length is subtracted and no PBC is used for f. */
503 static void gmx_simdcall calc_dr_x_f_simd(int b0,
505 gmx::ArrayRef<const AtomPair> atoms,
506 const rvec* gmx_restrict x,
507 const rvec* gmx_restrict f,
508 const real* gmx_restrict blc,
509 const real* pbc_simd,
510 rvec* gmx_restrict r,
511 real* gmx_restrict rhs,
512 real* gmx_restrict sol)
514 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
516 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
518 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
523 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
525 SimdReal x0_S, y0_S, z0_S;
526 SimdReal x1_S, y1_S, z1_S;
527 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
528 SimdReal fx_S, fy_S, fz_S, ip_S, rhs_S;
529 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
530 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
532 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
534 offset0[i] = atoms[bs + i].index1;
535 offset1[i] = atoms[bs + i].index2;
538 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset0, &x0_S, &y0_S, &z0_S);
539 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset1, &x1_S, &y1_S, &z1_S);
544 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
546 n2_S = norm2(rx_S, ry_S, rz_S);
547 il_S = invsqrt(n2_S);
553 transposeScatterStoreU<3>(reinterpret_cast<real*>(r + bs), offset2, rx_S, ry_S, rz_S);
555 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(f), offset0, &x0_S, &y0_S, &z0_S);
556 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(f), offset1, &x1_S, &y1_S, &z1_S);
561 ip_S = iprod(rx_S, ry_S, rz_S, fx_S, fy_S, fz_S);
563 rhs_S = load<SimdReal>(blc + bs) * ip_S;
565 store(rhs + bs, rhs_S);
566 store(sol + bs, rhs_S);
569 #endif // GMX_SIMD_HAVE_REAL
571 /*! \brief LINCS projection, works on derivatives of the coordinates. */
572 static void do_lincsp(const rvec* x,
579 ConstraintVariable econq,
584 const int b0 = lincsd->task[th].b0;
585 const int b1 = lincsd->task[th].b1;
587 gmx::ArrayRef<const AtomPair> atoms = lincsd->atoms;
588 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
589 gmx::ArrayRef<const int> blnr = lincsd->blnr;
590 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
592 gmx::ArrayRef<const real> blc;
593 gmx::ArrayRef<const real> blmf;
594 if (econq != ConstraintVariable::Force)
596 /* Use mass-weighted parameters */
602 /* Use non mass-weighted parameters */
604 blmf = lincsd->blmf1;
606 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
607 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
608 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
609 gmx::ArrayRef<real> sol = lincsd->tmp3;
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.
624 calc_dr_x_f_simd(b0, b1, atoms, x, f, blc.data(), pbc_simd, as_rvec_array(r.data()),
625 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, th, 1.0, sol, r,
710 (econq != ConstraintVariable::Force) ? invmass : nullptr, fp);
715 for (int b = b0; b < b1; b++)
717 dhdlambda -= sol[b] * lincsd->ddist[b];
720 lincsd->task[th].dhdlambda = dhdlambda;
725 /* Constraint virial,
726 * determines sum r_bond x delta f,
727 * where delta f is the constraint correction
728 * of the quantity that is being constrained.
730 for (int b = b0; b < b1; b++)
732 const real mvb = lincsd->bllen[b] * sol[b];
733 for (int i = 0; i < DIM; i++)
735 const real tmp1 = mvb * r[b][i];
736 for (int j = 0; j < DIM; j++)
738 rmdf[i][j] += tmp1 * r[b][j];
741 } /* 23 ncons flops */
745 #if GMX_SIMD_HAVE_REAL
747 /*! \brief Calculate the constraint distance vectors r to project on from x.
749 * Determine the right-hand side of the matrix equation using coordinates xp. */
750 static void gmx_simdcall calc_dr_x_xp_simd(int b0,
752 gmx::ArrayRef<const AtomPair> atoms,
753 const rvec* gmx_restrict x,
754 const rvec* gmx_restrict xp,
755 const real* gmx_restrict bllen,
756 const real* gmx_restrict blc,
757 const real* pbc_simd,
758 rvec* gmx_restrict r,
759 real* gmx_restrict rhs,
760 real* gmx_restrict sol)
762 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
763 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
765 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
770 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
772 SimdReal x0_S, y0_S, z0_S;
773 SimdReal x1_S, y1_S, z1_S;
774 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
775 SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
776 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
777 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
779 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
781 offset0[i] = atoms[bs + i].index1;
782 offset1[i] = atoms[bs + i].index2;
785 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset0, &x0_S, &y0_S, &z0_S);
786 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset1, &x1_S, &y1_S, &z1_S);
791 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
793 n2_S = norm2(rx_S, ry_S, rz_S);
794 il_S = invsqrt(n2_S);
800 transposeScatterStoreU<3>(reinterpret_cast<real*>(r + bs), offset2, rx_S, ry_S, rz_S);
802 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(xp), offset0, &x0_S, &y0_S, &z0_S);
803 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(xp), offset1, &x1_S, &y1_S, &z1_S);
809 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
811 ip_S = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
813 rhs_S = load<SimdReal>(blc + bs) * (ip_S - load<SimdReal>(bllen + bs));
815 store(rhs + bs, rhs_S);
816 store(sol + bs, rhs_S);
819 #endif // GMX_SIMD_HAVE_REAL
821 /*! \brief Determine the distances and right-hand side for the next iteration. */
822 gmx_unused static void calc_dist_iter(int b0,
824 gmx::ArrayRef<const AtomPair> atoms,
825 const rvec* gmx_restrict xp,
826 const real* gmx_restrict bllen,
827 const real* gmx_restrict blc,
830 real* gmx_restrict rhs,
831 real* gmx_restrict sol,
834 for (int b = b0; b < b1; b++)
840 pbc_dx_aiuc(pbc, xp[atoms[b].index1], xp[atoms[b].index2], dx);
844 rvec_sub(xp[atoms[b].index1], xp[atoms[b].index2], dx);
846 real len2 = len * len;
847 real dlen2 = 2 * len2 - ::norm2(dx);
848 if (dlen2 < wfac * len2)
850 /* not race free - see detailed comment in caller */
856 mvb = blc[b] * (len - dlen2 * gmx::invsqrt(dlen2));
864 } /* 20*ncons flops */
867 #if GMX_SIMD_HAVE_REAL
868 /*! \brief As calc_dist_iter(), but using SIMD intrinsics. */
869 static void gmx_simdcall calc_dist_iter_simd(int b0,
871 gmx::ArrayRef<const AtomPair> atoms,
872 const rvec* gmx_restrict x,
873 const real* gmx_restrict bllen,
874 const real* gmx_restrict blc,
875 const real* pbc_simd,
877 real* gmx_restrict rhs,
878 real* gmx_restrict sol,
881 SimdReal min_S(GMX_REAL_MIN);
883 SimdReal wfac_S(wfac);
886 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
888 /* Initialize all to FALSE */
889 warn_S = (two_S < setZero());
891 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
893 SimdReal x0_S, y0_S, z0_S;
894 SimdReal x1_S, y1_S, z1_S;
895 SimdReal rx_S, ry_S, rz_S, n2_S;
896 SimdReal len_S, len2_S, dlen2_S, lc_S, blc_S;
897 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
898 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
900 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
902 offset0[i] = atoms[bs + i].index1;
903 offset1[i] = atoms[bs + i].index2;
906 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset0, &x0_S, &y0_S, &z0_S);
907 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset1, &x1_S, &y1_S, &z1_S);
913 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
915 n2_S = norm2(rx_S, ry_S, rz_S);
917 len_S = load<SimdReal>(bllen + bs);
918 len2_S = len_S * len_S;
920 dlen2_S = fms(two_S, len2_S, n2_S);
922 warn_S = warn_S || (dlen2_S < (wfac_S * len2_S));
924 /* Avoid 1/0 by taking the max with REAL_MIN.
925 * Note: when dlen2 is close to zero (90 degree constraint rotation),
926 * the accuracy of the algorithm is no longer relevant.
928 dlen2_S = max(dlen2_S, min_S);
930 lc_S = fnma(dlen2_S, invsqrt(dlen2_S), len_S);
932 blc_S = load<SimdReal>(blc + bs);
936 store(rhs + bs, lc_S);
937 store(sol + bs, lc_S);
945 #endif // GMX_SIMD_HAVE_REAL
947 //! Implements LINCS constraining.
948 static void do_lincs(const rvec* x,
960 rvec* gmx_restrict v,
964 const int b0 = lincsd->task[th].b0;
965 const int b1 = lincsd->task[th].b1;
967 gmx::ArrayRef<const AtomPair> atoms = lincsd->atoms;
968 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
969 gmx::ArrayRef<const int> blnr = lincsd->blnr;
970 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
971 gmx::ArrayRef<const real> blc = lincsd->blc;
972 gmx::ArrayRef<const real> blmf = lincsd->blmf;
973 gmx::ArrayRef<const real> bllen = lincsd->bllen;
974 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
975 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
976 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
977 gmx::ArrayRef<real> sol = lincsd->tmp3;
978 gmx::ArrayRef<real> blc_sol = lincsd->tmp4;
979 gmx::ArrayRef<real> mlambda = lincsd->mlambda;
980 gmx::ArrayRef<const int> nlocat = lincsd->nlocat;
982 #if GMX_SIMD_HAVE_REAL
984 /* This SIMD code does the same as the plain-C code after the #else.
985 * The only difference is that we always call pbc code, as with SIMD
986 * the overhead of pbc computation (when not needed) is small.
988 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9 * GMX_SIMD_REAL_WIDTH];
990 /* Convert the pbc struct for SIMD */
991 set_pbc_simd(pbc, pbc_simd);
993 /* Compute normalized x i-j vectors, store in r.
994 * Compute the inner product of r and xp i-j and store in rhs1.
996 calc_dr_x_xp_simd(b0, b1, atoms, x, xp, bllen.data(), blc.data(), pbc_simd,
997 as_rvec_array(r.data()), rhs1.data(), sol.data());
999 #else // GMX_SIMD_HAVE_REAL
1003 /* Compute normalized i-j vectors */
1004 for (int b = b0; b < b1; b++)
1007 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
1010 pbc_dx_aiuc(pbc, xp[atoms[b].index1], xp[atoms[b].index2], dx);
1011 real mvb = blc[b] * (::iprod(r[b], dx) - bllen[b]);
1018 /* Compute normalized i-j vectors */
1019 for (int b = b0; b < b1; b++)
1021 int i = atoms[b].index1;
1022 int j = atoms[b].index2;
1023 real tmp0 = x[i][0] - x[j][0];
1024 real tmp1 = x[i][1] - x[j][1];
1025 real tmp2 = x[i][2] - x[j][2];
1026 real rlen = gmx::invsqrt(tmp0 * tmp0 + tmp1 * tmp1 + tmp2 * tmp2);
1027 r[b][0] = rlen * tmp0;
1028 r[b][1] = rlen * tmp1;
1029 r[b][2] = rlen * tmp2;
1030 /* 16 ncons flops */
1033 * (r[b][0] * (xp[i][0] - xp[j][0]) + r[b][1] * (xp[i][1] - xp[j][1])
1034 + r[b][2] * (xp[i][2] - xp[j][2]) - bllen[b]);
1039 /* Together: 26*ncons + 6*nrtot flops */
1042 #endif // GMX_SIMD_HAVE_REAL
1044 if (lincsd->bTaskDep)
1046 /* We need a barrier, since the matrix construction below
1047 * can access entries in r of other threads.
1052 /* Construct the (sparse) LINCS matrix */
1053 for (int b = b0; b < b1; b++)
1055 for (int n = blnr[b]; n < blnr[b + 1]; n++)
1057 blcc[n] = blmf[n] * gmx::dot(r[b], r[blbnb[n]]);
1060 /* Together: 26*ncons + 6*nrtot flops */
1062 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
1063 /* nrec*(ncons+2*nrtot) flops */
1065 #if GMX_SIMD_HAVE_REAL
1066 for (int b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1068 SimdReal t1 = load<SimdReal>(blc.data() + b);
1069 SimdReal t2 = load<SimdReal>(sol.data() + b);
1070 store(mlambda.data() + b, t1 * t2);
1073 for (int b = b0; b < b1; b++)
1075 mlambda[b] = blc[b] * sol[b];
1077 #endif // GMX_SIMD_HAVE_REAL
1079 /* Update the coordinates */
1080 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1083 ******** Correction for centripetal effects ********
1088 wfac = std::cos(DEG2RAD * wangle);
1091 for (int iter = 0; iter < lincsd->nIter; iter++)
1093 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
1098 /* Communicate the corrected non-local coordinates */
1099 if (DOMAINDECOMP(cr))
1101 dd_move_x_constraints(cr->dd, box, xp, nullptr, FALSE);
1106 else if (lincsd->bTaskDep)
1111 #if GMX_SIMD_HAVE_REAL
1112 calc_dist_iter_simd(b0, b1, atoms, xp, bllen.data(), blc.data(), pbc_simd, wfac,
1113 rhs1.data(), sol.data(), bWarn);
1115 calc_dist_iter(b0, b1, atoms, xp, bllen.data(), blc.data(), pbc, wfac, rhs1.data(),
1117 /* 20*ncons flops */
1118 #endif // GMX_SIMD_HAVE_REAL
1120 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
1121 /* nrec*(ncons+2*nrtot) flops */
1123 #if GMX_SIMD_HAVE_REAL
1124 for (int b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1126 SimdReal t1 = load<SimdReal>(blc.data() + b);
1127 SimdReal t2 = load<SimdReal>(sol.data() + b);
1128 SimdReal mvb = t1 * t2;
1129 store(blc_sol.data() + b, mvb);
1130 store(mlambda.data() + b, load<SimdReal>(mlambda.data() + b) + mvb);
1133 for (int b = b0; b < b1; b++)
1135 real mvb = blc[b] * sol[b];
1139 #endif // GMX_SIMD_HAVE_REAL
1141 /* Update the coordinates */
1142 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1144 /* nit*ncons*(37+9*nrec) flops */
1148 /* Update the velocities */
1149 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1150 /* 16 ncons flops */
1153 if (!nlocat.empty() && (bCalcDHDL || bCalcVir))
1155 if (lincsd->bTaskDep)
1157 /* In lincs_update_atoms threads might cross-read mlambda */
1161 /* Only account for local atoms */
1162 for (int b = b0; b < b1; b++)
1164 mlambda[b] *= 0.5 * nlocat[b];
1171 for (int b = b0; b < b1; b++)
1173 /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1174 * later after the contributions are reduced over the threads.
1176 dhdl -= lincsd->mlambda[b] * lincsd->ddist[b];
1178 lincsd->task[th].dhdlambda = dhdl;
1183 /* Constraint virial */
1184 for (int b = b0; b < b1; b++)
1186 real tmp0 = -bllen[b] * mlambda[b];
1187 for (int i = 0; i < DIM; i++)
1189 real tmp1 = tmp0 * r[b][i];
1190 for (int j = 0; j < DIM; j++)
1192 vir_r_m_dr[i][j] -= tmp1 * r[b][j];
1195 } /* 22 ncons flops */
1199 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1200 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1202 * (26+nrec)*ncons + (6+2*nrec)*nrtot
1203 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1205 * (63+nrec)*ncons + (6+4*nrec)*nrtot
1209 /*! \brief Sets the elements in the LINCS matrix for task task. */
1210 static void set_lincs_matrix_task(Lincs* li, Task* li_task, const real* invmass, int* ncc_triangle, int* nCrossTaskTriangles)
1212 /* Construct the coupling coefficient matrix blmf */
1213 li_task->ntriangle = 0;
1215 *nCrossTaskTriangles = 0;
1216 for (int i = li_task->b0; i < li_task->b1; i++)
1218 const int a1 = li->atoms[i].index1;
1219 const int a2 = li->atoms[i].index2;
1220 for (int n = li->blnr[i]; n < li->blnr[i + 1]; n++)
1222 const int k = li->blbnb[n];
1224 /* If we are using multiple, independent tasks for LINCS,
1225 * the calls to check_assign_connected should have
1226 * put all connected constraints in our task.
1228 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1231 if (a1 == li->atoms[k].index1 || a2 == li->atoms[k].index2)
1241 if (a1 == li->atoms[k].index1 || a1 == li->atoms[k].index2)
1251 li->blmf[n] = sign * invmass[center] * li->blc[i] * li->blc[k];
1252 li->blmf1[n] = sign * 0.5;
1253 if (li->ncg_triangle > 0)
1255 /* Look for constraint triangles */
1256 for (int nk = li->blnr[k]; nk < li->blnr[k + 1]; nk++)
1258 const int kk = li->blbnb[nk];
1259 if (kk != i && kk != k && (li->atoms[kk].index1 == end || li->atoms[kk].index2 == end))
1261 /* Check if the constraints in this triangle actually
1262 * belong to a different task. We still assign them
1263 * here, since it's convenient for the triangle
1264 * iterations, but we then need an extra barrier.
1266 if (k < li_task->b0 || k >= li_task->b1 || kk < li_task->b0 || kk >= li_task->b1)
1268 (*nCrossTaskTriangles)++;
1271 if (li_task->ntriangle == 0 || li_task->triangle[li_task->ntriangle - 1] < i)
1273 /* Add this constraint to the triangle list */
1274 li_task->triangle[li_task->ntriangle] = i;
1275 li_task->tri_bits[li_task->ntriangle] = 0;
1276 li_task->ntriangle++;
1277 if (li->blnr[i + 1] - li->blnr[i]
1278 > static_cast<int>(sizeof(li_task->tri_bits[0]) * 8 - 1))
1281 "A constraint is connected to %d constraints, this is "
1282 "more than the %zu allowed for constraints participating "
1284 li->blnr[i + 1] - li->blnr[i],
1285 sizeof(li_task->tri_bits[0]) * 8 - 1);
1288 li_task->tri_bits[li_task->ntriangle - 1] |= (1 << (n - li->blnr[i]));
1297 /*! \brief Sets the elements in the LINCS matrix. */
1298 static void set_lincs_matrix(Lincs* li, real* invmass, real lambda)
1300 const real invsqrt2 = 0.7071067811865475244;
1302 for (int i = 0; (i < li->nc); i++)
1304 const int a1 = li->atoms[i].index1;
1305 const int a2 = li->atoms[i].index2;
1306 li->blc[i] = gmx::invsqrt(invmass[a1] + invmass[a2]);
1307 li->blc1[i] = invsqrt2;
1310 /* Construct the coupling coefficient matrix blmf */
1311 int ntriangle = 0, ncc_triangle = 0, nCrossTaskTriangles = 0;
1312 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle, nCrossTaskTriangles) num_threads(li->ntask) schedule(static)
1313 for (int th = 0; th < li->ntask; th++)
1317 set_lincs_matrix_task(li, &li->task[th], invmass, &ncc_triangle, &nCrossTaskTriangles);
1318 ntriangle += li->task[th].ntriangle;
1320 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1322 li->ntriangle = ntriangle;
1323 li->ncc_triangle = ncc_triangle;
1324 li->bTaskDepTri = (nCrossTaskTriangles > 0);
1328 fprintf(debug, "The %d constraints participate in %d triangles\n", li->nc, li->ntriangle);
1329 fprintf(debug, "There are %d constraint couplings, of which %d in triangles\n", li->ncc,
1331 if (li->ntriangle > 0 && li->ntask > 1)
1334 "%d constraint triangles contain constraints assigned to different tasks\n",
1335 nCrossTaskTriangles);
1340 * so we know with which lambda value the masses have been set.
1342 li->matlam = lambda;
1345 //! Finds all triangles of atoms that share constraints to a central atom.
1346 static int count_triangle_constraints(const InteractionLists& ilist, const ListOfLists<int>& at2con)
1348 const int ncon1 = ilist[F_CONSTR].size() / 3;
1349 const int ncon_tot = ncon1 + ilist[F_CONSTRNC].size() / 3;
1351 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1352 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1354 int ncon_triangle = 0;
1355 for (int c0 = 0; c0 < ncon_tot; c0++)
1357 bool bTriangle = FALSE;
1358 const int* iap = constr_iatomptr(ia1, ia2, c0);
1359 const int a00 = iap[1];
1360 const int a01 = iap[2];
1361 for (const int c1 : at2con[a01])
1365 const int* iap = constr_iatomptr(ia1, ia2, c1);
1366 const int a10 = iap[1];
1367 const int a11 = iap[2];
1377 for (const int c2 : at2con[ac1])
1379 if (c2 != c0 && c2 != c1)
1381 const int* iap = constr_iatomptr(ia1, ia2, c2);
1382 const int a20 = iap[1];
1383 const int a21 = iap[2];
1384 if (a20 == a00 || a21 == a00)
1398 return ncon_triangle;
1401 //! Finds sequences of sequential constraints.
1402 static bool more_than_two_sequential_constraints(const InteractionLists& ilist, const ListOfLists<int>& at2con)
1404 const int ncon1 = ilist[F_CONSTR].size() / 3;
1405 const int ncon_tot = ncon1 + ilist[F_CONSTRNC].size() / 3;
1407 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1408 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1410 for (int c = 0; c < ncon_tot; c++)
1412 const int* iap = constr_iatomptr(ia1, ia2, c);
1413 const int a1 = iap[1];
1414 const int a2 = iap[2];
1415 /* Check if this constraint has constraints connected at both atoms */
1416 if (at2con[a1].ssize() > 1 && at2con[a2].ssize() > 1)
1425 Lincs* init_lincs(FILE* fplog,
1426 const gmx_mtop_t& mtop,
1427 int nflexcon_global,
1428 ArrayRef<const ListOfLists<int>> atomToConstraintsPerMolType,
1433 // TODO this should become a unique_ptr
1435 bool bMoreThanTwoSeq;
1439 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n", bPLINCS ? " Parallel" : "");
1444 li->ncg = gmx_mtop_ftype_count(mtop, F_CONSTR) + gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1445 li->ncg_flex = nflexcon_global;
1448 li->nOrder = nProjOrder;
1450 li->max_connect = 0;
1451 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1453 const auto& at2con = atomToConstraintsPerMolType[mt];
1454 for (int a = 0; a < mtop.moltype[mt].atoms.nr; a++)
1456 li->max_connect = std::max(li->max_connect, int(at2con[a].ssize()));
1460 li->ncg_triangle = 0;
1461 bMoreThanTwoSeq = FALSE;
1462 for (const gmx_molblock_t& molb : mtop.molblock)
1464 const gmx_moltype_t& molt = mtop.moltype[molb.type];
1465 const auto& at2con = atomToConstraintsPerMolType[molb.type];
1467 li->ncg_triangle += molb.nmol * count_triangle_constraints(molt.ilist, at2con);
1469 if (!bMoreThanTwoSeq && more_than_two_sequential_constraints(molt.ilist, at2con))
1471 bMoreThanTwoSeq = TRUE;
1475 /* Check if we need to communicate not only before LINCS,
1476 * but also before each iteration.
1477 * The check for only two sequential constraints is only
1478 * useful for the common case of H-bond only constraints.
1479 * With more effort we could also make it useful for small
1480 * molecules with nr. sequential constraints <= nOrder-1.
1482 li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1484 if (debug && bPLINCS)
1486 fprintf(debug, "PLINCS communication before each iteration: %d\n", static_cast<int>(li->bCommIter));
1489 /* LINCS can run on any number of threads.
1490 * Currently the number is fixed for the whole simulation,
1491 * but it could be set in set_lincs().
1492 * The current constraint to task assignment code can create independent
1493 * tasks only when not more than two constraints are connected sequentially.
1495 li->ntask = gmx_omp_nthreads_get(emntLINCS);
1496 li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1499 fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n", li->ntask,
1500 li->bTaskDep ? "" : "in");
1508 /* Allocate an extra elements for "task-overlap" constraints */
1509 li->task.resize(li->ntask + 1);
1512 if (bPLINCS || li->ncg_triangle > 0)
1514 please_cite(fplog, "Hess2008a");
1518 please_cite(fplog, "Hess97a");
1523 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1527 "There are constraints between atoms in different decomposition domains,\n"
1528 "will communicate selected coordinates each lincs iteration\n");
1530 if (li->ncg_triangle > 0)
1533 "%d constraints are involved in constraint triangles,\n"
1534 "will apply an additional matrix expansion of order %d for couplings\n"
1535 "between constraints inside triangles\n",
1536 li->ncg_triangle, li->nOrder);
1543 void done_lincs(Lincs* li)
1548 /*! \brief Sets up the work division over the threads. */
1549 static void lincs_thread_setup(Lincs* li, int natoms)
1551 li->atf.resize(natoms);
1553 gmx::ArrayRef<gmx_bitmask_t> atf = li->atf;
1555 /* Clear the atom flags */
1556 for (gmx_bitmask_t& mask : atf)
1558 bitmask_clear(&mask);
1561 if (li->ntask > BITMASK_SIZE)
1563 gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1566 for (int th = 0; th < li->ntask; th++)
1568 const Task* li_task = &li->task[th];
1570 /* For each atom set a flag for constraints from each */
1571 for (int b = li_task->b0; b < li_task->b1; b++)
1573 bitmask_set_bit(&atf[li->atoms[b].index1], th);
1574 bitmask_set_bit(&atf[li->atoms[b].index2], th);
1578 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1579 for (int th = 0; th < li->ntask; th++)
1587 li_task = &li->task[th];
1589 bitmask_init_low_bits(&mask, th);
1591 li_task->ind.clear();
1592 li_task->ind_r.clear();
1593 for (b = li_task->b0; b < li_task->b1; b++)
1595 /* We let the constraint with the lowest thread index
1596 * operate on atoms with constraints from multiple threads.
1598 if (bitmask_is_disjoint(atf[li->atoms[b].index1], mask)
1599 && bitmask_is_disjoint(atf[li->atoms[b].index2], mask))
1601 /* Add the constraint to the local atom update index */
1602 li_task->ind.push_back(b);
1606 /* Add the constraint to the rest block */
1607 li_task->ind_r.push_back(b);
1611 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1614 /* We need to copy all constraints which have not be assigned
1615 * to a thread to a separate list which will be handled by one thread.
1617 Task* li_m = &li->task[li->ntask];
1620 for (int th = 0; th < li->ntask; th++)
1622 const Task& li_task = li->task[th];
1624 for (int ind_r : li_task.ind_r)
1626 li_m->ind.push_back(ind_r);
1631 fprintf(debug, "LINCS thread %d: %zu constraints\n", th, li_task.ind.size());
1637 fprintf(debug, "LINCS thread r: %zu constraints\n", li_m->ind.size());
1641 //! Assign a constraint.
1642 static void assign_constraint(Lincs* li,
1643 int constraint_index,
1648 const ListOfLists<int>& at2con)
1654 /* Make an mapping of local topology constraint index to LINCS index */
1655 li->con_index[constraint_index] = con;
1657 li->bllen0[con] = lenA;
1658 li->ddist[con] = lenB - lenA;
1659 /* Set the length to the topology A length */
1660 li->bllen[con] = lenA;
1661 li->atoms[con].index1 = a1;
1662 li->atoms[con].index2 = a2;
1664 /* Make space in the constraint connection matrix for constraints
1665 * connected to both end of the current constraint.
1667 li->ncc += at2con[a1].ssize() - 1 + at2con[a2].ssize() - 1;
1669 li->blnr[con + 1] = li->ncc;
1671 /* Increase the constraint count */
1675 /*! \brief Check if constraint with topology index constraint_index is connected
1676 * to other constraints, and if so add those connected constraints to our task. */
1677 static void check_assign_connected(Lincs* li,
1678 const t_iatom* iatom,
1683 const ListOfLists<int>& at2con)
1685 /* Currently this function only supports constraint groups
1686 * in which all constraints share at least one atom
1687 * (e.g. H-bond constraints).
1688 * Check both ends of the current constraint for
1689 * connected constraints. We need to assign those
1692 for (int end = 0; end < 2; end++)
1694 const int a = (end == 0 ? a1 : a2);
1696 for (const int cc : at2con[a])
1698 /* Check if constraint cc has not yet been assigned */
1699 if (li->con_index[cc] == -1)
1701 const int type = iatom[cc * 3];
1702 const real lenA = idef.iparams[type].constr.dA;
1703 const real lenB = idef.iparams[type].constr.dB;
1705 if (bDynamics || lenA != 0 || lenB != 0)
1707 assign_constraint(li, cc, iatom[3 * cc + 1], iatom[3 * cc + 2], lenA, lenB, at2con);
1714 /*! \brief Check if constraint with topology index constraint_index is involved
1715 * in a constraint triangle, and if so add the other two constraints
1716 * in the triangle to our task. */
1717 static void check_assign_triangle(Lincs* li,
1718 const t_iatom* iatom,
1721 int constraint_index,
1724 const ListOfLists<int>& at2con)
1726 int nca, cc[32], ca[32];
1727 int c_triangle[2] = { -1, -1 };
1730 for (const int c : at2con[a1])
1732 if (c != constraint_index)
1736 aa1 = iatom[c * 3 + 1];
1737 aa2 = iatom[c * 3 + 2];
1753 for (const int c : at2con[a2])
1755 if (c != constraint_index)
1759 aa1 = iatom[c * 3 + 1];
1760 aa2 = iatom[c * 3 + 2];
1763 for (i = 0; i < nca; i++)
1767 c_triangle[0] = cc[i];
1774 for (i = 0; i < nca; i++)
1778 c_triangle[0] = cc[i];
1786 if (c_triangle[0] >= 0)
1790 for (end = 0; end < 2; end++)
1792 /* Check if constraint c_triangle[end] has not yet been assigned */
1793 if (li->con_index[c_triangle[end]] == -1)
1798 i = c_triangle[end] * 3;
1800 lenA = idef.iparams[type].constr.dA;
1801 lenB = idef.iparams[type].constr.dB;
1803 if (bDynamics || lenA != 0 || lenB != 0)
1805 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1812 //! Sets matrix indices.
1813 static void set_matrix_indices(Lincs* li, const Task& li_task, const ListOfLists<int>& at2con, bool bSortMatrix)
1815 for (int b = li_task.b0; b < li_task.b1; b++)
1817 const int a1 = li->atoms[b].index1;
1818 const int a2 = li->atoms[b].index2;
1820 int i = li->blnr[b];
1821 for (const int constraint : at2con[a1])
1823 const int concon = li->con_index[constraint];
1826 li->blbnb[i++] = concon;
1829 for (const int constraint : at2con[a2])
1831 const int concon = li->con_index[constraint];
1834 li->blbnb[i++] = concon;
1840 /* Order the blbnb matrix to optimize memory access */
1841 std::sort(li->blbnb.begin() + li->blnr[b], li->blbnb.begin() + li->blnr[b + 1]);
1846 void set_lincs(const t_idef& idef, const t_mdatoms& md, bool bDynamics, const t_commrec* cr, Lincs* li)
1854 /* Zero the thread index ranges.
1855 * Otherwise without local constraints we could return with old ranges.
1857 for (int i = 0; i < li->ntask; i++)
1861 li->task[i].ind.clear();
1865 li->task[li->ntask].ind.clear();
1868 /* This is the local topology, so there are only F_CONSTR constraints */
1869 if (idef.il[F_CONSTR].nr == 0)
1871 /* There are no constraints,
1872 * we do not need to fill any data structures.
1879 fprintf(debug, "Building the LINCS connectivity\n");
1882 if (DOMAINDECOMP(cr))
1884 if (cr->dd->constraints)
1888 dd_get_constraint_range(cr->dd, &start, &natoms);
1892 natoms = dd_numHomeAtoms(*cr->dd);
1900 const ListOfLists<int> at2con =
1901 make_at2con(natoms, idef.il, idef.iparams, flexibleConstraintTreatment(bDynamics));
1903 const int ncon_tot = idef.il[F_CONSTR].nr / 3;
1905 /* Ensure we have enough padding for aligned loads for each thread */
1906 const int numEntries = ncon_tot + li->ntask * simd_width;
1907 li->con_index.resize(numEntries);
1908 li->bllen0.resize(numEntries);
1909 li->ddist.resize(numEntries);
1910 li->atoms.resize(numEntries);
1911 li->blc.resize(numEntries);
1912 li->blc1.resize(numEntries);
1913 li->blnr.resize(numEntries + 1);
1914 li->bllen.resize(numEntries);
1915 li->tmpv.resizeWithPadding(numEntries);
1916 if (DOMAINDECOMP(cr))
1918 li->nlocat.resize(numEntries);
1920 li->tmp1.resize(numEntries);
1921 li->tmp2.resize(numEntries);
1922 li->tmp3.resize(numEntries);
1923 li->tmp4.resize(numEntries);
1924 li->mlambda.resize(numEntries);
1926 iatom = idef.il[F_CONSTR].iatoms;
1928 li->blnr[0] = li->ncc;
1930 /* Assign the constraints for li->ntask LINCS tasks.
1931 * We target a uniform distribution of constraints over the tasks.
1932 * Note that when flexible constraints are present, but are removed here
1933 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
1934 * happen during normal MD, that's ok.
1937 /* Determine the number of constraints we need to assign here */
1938 int ncon_assign = ncon_tot;
1941 /* With energy minimization, flexible constraints are ignored
1942 * (and thus minimized, as they should be).
1944 ncon_assign -= countFlexibleConstraints(idef.il, idef.iparams);
1947 /* Set the target constraint count per task to exactly uniform,
1948 * this might be overridden below.
1950 int ncon_target = (ncon_assign + li->ntask - 1) / li->ntask;
1952 /* Mark all constraints as unassigned by setting their index to -1 */
1953 for (int con = 0; con < ncon_tot; con++)
1955 li->con_index[con] = -1;
1959 for (int th = 0; th < li->ntask; th++)
1963 li_task = &li->task[th];
1965 #if GMX_SIMD_HAVE_REAL
1966 /* With indepedent tasks we likely have H-bond constraints or constraint
1967 * pairs. The connected constraints will be pulled into the task, so the
1968 * constraints per task will often exceed ncon_target.
1969 * Triangle constraints can also increase the count, but there are
1970 * relatively few of those, so we usually expect to get ncon_target.
1974 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
1975 * since otherwise a lot of operations can be wasted.
1976 * There are several ways to round here, we choose the one
1977 * that alternates block sizes, which helps with Intel HT.
1979 ncon_target = ((ncon_assign * (th + 1)) / li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1)
1980 & ~(GMX_SIMD_REAL_WIDTH - 1);
1982 #endif // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
1984 /* Continue filling the arrays where we left off with the previous task,
1985 * including padding for SIMD.
1987 li_task->b0 = li->nc;
1989 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
1991 if (li->con_index[con] == -1)
1996 type = iatom[3 * con];
1997 a1 = iatom[3 * con + 1];
1998 a2 = iatom[3 * con + 2];
1999 lenA = idef.iparams[type].constr.dA;
2000 lenB = idef.iparams[type].constr.dB;
2001 /* Skip the flexible constraints when not doing dynamics */
2002 if (bDynamics || lenA != 0 || lenB != 0)
2004 assign_constraint(li, con, a1, a2, lenA, lenB, at2con);
2006 if (li->ntask > 1 && !li->bTaskDep)
2008 /* We can generate independent tasks. Check if we
2009 * need to assign connected constraints to our task.
2011 check_assign_connected(li, iatom, idef, bDynamics, a1, a2, at2con);
2013 if (li->ntask > 1 && li->ncg_triangle > 0)
2015 /* Ensure constraints in one triangle are assigned
2018 check_assign_triangle(li, iatom, idef, bDynamics, con, a1, a2, at2con);
2026 li_task->b1 = li->nc;
2030 /* Copy the last atom pair indices and lengths for constraints
2031 * up to a multiple of simd_width, such that we can do all
2032 * SIMD operations without having to worry about end effects.
2036 li->nc = ((li_task->b1 + simd_width - 1) / simd_width) * simd_width;
2037 last = li_task->b1 - 1;
2038 for (i = li_task->b1; i < li->nc; i++)
2040 li->atoms[i] = li->atoms[last];
2041 li->bllen0[i] = li->bllen0[last];
2042 li->ddist[i] = li->ddist[last];
2043 li->bllen[i] = li->bllen[last];
2044 li->blnr[i + 1] = li->blnr[last + 1];
2048 /* Keep track of how many constraints we assigned */
2049 li->nc_real += li_task->b1 - li_task->b0;
2053 fprintf(debug, "LINCS task %d constraints %d - %d\n", th, li_task->b0, li_task->b1);
2057 assert(li->nc_real == ncon_assign);
2061 /* Without DD we order the blbnb matrix to optimize memory access.
2062 * With DD the overhead of sorting is more than the gain during access.
2064 bSortMatrix = !DOMAINDECOMP(cr);
2066 li->blbnb.resize(li->ncc);
2068 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2069 for (int th = 0; th < li->ntask; th++)
2073 Task& li_task = li->task[th];
2075 if (li->ncg_triangle > 0)
2077 /* This is allocating too much, but it is difficult to improve */
2078 li_task.triangle.resize(li_task.b1 - li_task.b0);
2079 li_task.tri_bits.resize(li_task.b1 - li_task.b0);
2082 set_matrix_indices(li, li_task, at2con, bSortMatrix);
2084 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2087 if (cr->dd == nullptr)
2089 /* Since the matrix is static, we should free some memory */
2090 li->blbnb.resize(li->ncc);
2093 li->blmf.resize(li->ncc);
2094 li->blmf1.resize(li->ncc);
2095 li->tmpncc.resize(li->ncc);
2097 gmx::ArrayRef<const int> nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2098 if (!nlocat_dd.empty())
2100 /* Convert nlocat from local topology to LINCS constraint indexing */
2101 for (con = 0; con < ncon_tot; con++)
2103 li->nlocat[li->con_index[con]] = nlocat_dd[con];
2113 fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n", li->nc_real,
2119 lincs_thread_setup(li, md.nr);
2122 set_lincs_matrix(li, md.invmass, md.lambda);
2125 //! Issues a warning when LINCS constraints cannot be satisfied.
2126 static void lincs_warning(gmx_domdec_t* dd,
2131 gmx::ArrayRef<const AtomPair> atoms,
2132 gmx::ArrayRef<const real> bllen,
2137 real wfac = std::cos(DEG2RAD * wangle);
2140 "bonds that rotated more than %g degrees:\n"
2141 " atom 1 atom 2 angle previous, current, constraint length\n",
2144 for (int b = 0; b < ncons; b++)
2146 const int i = atoms[b].index1;
2147 const int j = atoms[b].index2;
2152 pbc_dx_aiuc(pbc, x[i], x[j], v0);
2153 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2157 rvec_sub(x[i], x[j], v0);
2158 rvec_sub(xprime[i], xprime[j], v1);
2162 real cosine = ::iprod(v0, v1) / (d0 * d1);
2165 fprintf(stderr, " %6d %6d %5.1f %8.4f %8.4f %8.4f\n", ddglatnr(dd, i),
2166 ddglatnr(dd, j), RAD2DEG * std::acos(cosine), d0, d1, bllen[b]);
2167 if (!std::isfinite(d1))
2169 gmx_fatal(FARGS, "Bond length not finite.");
2175 if (*warncount > maxwarn)
2177 too_many_constraint_warnings(econtLINCS, *warncount);
2181 //! Status information about how well LINCS satisified the constraints in this domain
2182 struct LincsDeviations
2184 //! The maximum over all bonds in this domain of the relative deviation in bond lengths
2185 real maxDeviation = 0;
2186 //! Sum over all bonds in this domain of the squared relative deviation
2187 real sumSquaredDeviation = 0;
2188 //! Index of bond with max deviation
2189 int indexOfMaxDeviation = -1;
2190 //! Number of bonds constrained in this domain
2191 int numConstraints = 0;
2194 //! Determine how well the constraints have been satisfied.
2195 static LincsDeviations makeLincsDeviations(const Lincs& lincsd, const rvec* x, const t_pbc* pbc)
2197 LincsDeviations result;
2198 const ArrayRef<const AtomPair> atoms = lincsd.atoms;
2199 const ArrayRef<const real> bllen = lincsd.bllen;
2200 const ArrayRef<const int> nlocat = lincsd.nlocat;
2202 for (int task = 0; task < lincsd.ntask; task++)
2204 for (int b = lincsd.task[task].b0; b < lincsd.task[task].b1; b++)
2209 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
2213 rvec_sub(x[atoms[b].index1], x[atoms[b].index2], dx);
2215 real r2 = ::norm2(dx);
2216 real len = r2 * gmx::invsqrt(r2);
2217 real d = std::abs(len / bllen[b] - 1.0_real);
2218 if (d > result.maxDeviation && (nlocat.empty() || nlocat[b]))
2220 result.maxDeviation = d;
2221 result.indexOfMaxDeviation = b;
2225 result.sumSquaredDeviation += d * d;
2226 result.numConstraints++;
2230 result.sumSquaredDeviation += nlocat[b] * d * d;
2231 result.numConstraints += nlocat[b];
2236 if (!nlocat.empty())
2238 result.numConstraints /= 2;
2239 result.sumSquaredDeviation *= 0.5;
2244 bool constrain_lincs(bool computeRmsd,
2245 const t_inputrec& ir,
2248 const t_mdatoms& md,
2249 const t_commrec* cr,
2250 const gmx_multisim_t* ms,
2262 ConstraintVariable econq,
2269 /* This boolean should be set by a flag passed to this routine.
2270 * We can also easily check if any constraint length is changed,
2271 * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2273 bool bCalcDHDL = (ir.efep != efepNO && dvdlambda != nullptr);
2275 if (lincsd->nc == 0 && cr->dd == nullptr)
2279 lincsd->rmsdData = { { 0 } };
2285 if (econq == ConstraintVariable::Positions)
2287 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2288 * also with efep!=fepNO.
2290 if (ir.efep != efepNO)
2292 if (md.nMassPerturbed && lincsd->matlam != md.lambda)
2294 set_lincs_matrix(lincsd, md.invmass, md.lambda);
2297 for (int i = 0; i < lincsd->nc; i++)
2299 lincsd->bllen[i] = lincsd->bllen0[i] + lambda * lincsd->ddist[i];
2303 if (lincsd->ncg_flex)
2305 /* Set the flexible constraint lengths to the old lengths */
2308 for (int i = 0; i < lincsd->nc; i++)
2310 if (lincsd->bllen[i] == 0)
2313 pbc_dx_aiuc(pbc, x[lincsd->atoms[i].index1], x[lincsd->atoms[i].index2], dx);
2314 lincsd->bllen[i] = norm(dx);
2320 for (int i = 0; i < lincsd->nc; i++)
2322 if (lincsd->bllen[i] == 0)
2324 lincsd->bllen[i] = std::sqrt(
2325 distance2(x[lincsd->atoms[i].index1], x[lincsd->atoms[i].index2]));
2331 const bool printDebugOutput = ((debug != nullptr) && lincsd->nc > 0);
2332 if (printDebugOutput)
2334 LincsDeviations deviations = makeLincsDeviations(*lincsd, xprime, pbc);
2335 fprintf(debug, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2336 fprintf(debug, " Before LINCS %.6f %.6f %6d %6d\n",
2337 std::sqrt(deviations.sumSquaredDeviation / deviations.numConstraints),
2338 deviations.maxDeviation,
2339 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index1),
2340 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index2));
2343 /* This bWarn var can be updated by multiple threads
2344 * at the same time. But as we only need to detect
2345 * if a warning occurred or not, this is not an issue.
2349 /* The OpenMP parallel region of constrain_lincs for coords */
2350 #pragma omp parallel num_threads(lincsd->ntask)
2354 int th = gmx_omp_get_thread_num();
2356 clear_mat(lincsd->task[th].vir_r_m_dr);
2358 do_lincs(x, xprime, box, pbc, lincsd, th, md.invmass, cr, bCalcDHDL, ir.LincsWarnAngle,
2359 &bWarn, invdt, v, bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2361 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2364 if (computeRmsd || printDebugOutput || bWarn)
2366 LincsDeviations deviations = makeLincsDeviations(*lincsd, xprime, pbc);
2370 // This is reduced across domains in compute_globals and
2371 // reported to the log file.
2372 lincsd->rmsdData[0] = deviations.numConstraints;
2373 lincsd->rmsdData[1] = deviations.sumSquaredDeviation;
2377 // This is never read
2378 lincsd->rmsdData = { { 0 } };
2380 if (printDebugOutput)
2382 fprintf(debug, " After LINCS %.6f %.6f %6d %6d\n\n",
2383 std::sqrt(deviations.sumSquaredDeviation / deviations.numConstraints),
2384 deviations.maxDeviation,
2385 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index1),
2386 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index2));
2391 if (maxwarn < INT_MAX)
2393 std::string simMesg;
2396 simMesg += gmx::formatString(" in simulation %d", ms->sim);
2400 ", time %g (ps) LINCS WARNING%s\n"
2401 "relative constraint deviation after LINCS:\n"
2402 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2403 step, ir.init_t + step * ir.delta_t, simMesg.c_str(),
2404 std::sqrt(deviations.sumSquaredDeviation / deviations.numConstraints),
2405 deviations.maxDeviation,
2406 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index1),
2407 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index2));
2409 lincs_warning(cr->dd, x, xprime, pbc, lincsd->nc, lincsd->atoms, lincsd->bllen,
2410 ir.LincsWarnAngle, maxwarn, warncount);
2412 bOK = (deviations.maxDeviation < 0.5);
2416 if (lincsd->ncg_flex)
2418 for (int i = 0; (i < lincsd->nc); i++)
2420 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2422 lincsd->bllen[i] = 0;
2429 /* The OpenMP parallel region of constrain_lincs for derivatives */
2430 #pragma omp parallel num_threads(lincsd->ntask)
2434 int th = gmx_omp_get_thread_num();
2436 do_lincsp(x, xprime, min_proj, pbc, lincsd, th, md.invmass, econq, bCalcDHDL,
2437 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2439 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2445 /* Reduce the dH/dlambda contributions over the threads */
2450 for (th = 0; th < lincsd->ntask; th++)
2452 dhdlambda += lincsd->task[th].dhdlambda;
2454 if (econq == ConstraintVariable::Positions)
2456 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2457 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2458 dhdlambda /= ir.delta_t * ir.delta_t;
2460 *dvdlambda += dhdlambda;
2463 if (bCalcVir && lincsd->ntask > 1)
2465 for (int i = 1; i < lincsd->ntask; i++)
2467 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2471 /* count assuming nit=1 */
2472 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2473 inc_nrnb(nrnb, eNR_LINCSMAT, (2 + lincsd->nOrder) * lincsd->ncc);
2474 if (lincsd->ntriangle > 0)
2476 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder * lincsd->ncc_triangle);
2480 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real * 2);
2484 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);