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,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief Defines LINCS code.
40 * \author Berk Hess <hess@kth.se>
41 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \ingroup module_mdlib
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/gmxlib/nrnb.h"
60 #include "gromacs/math/functions.h"
61 #include "gromacs/math/paddedvector.h"
62 #include "gromacs/math/units.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/mdlib/constr.h"
65 #include "gromacs/mdlib/gmx_omp_nthreads.h"
66 #include "gromacs/mdrunutility/multisim.h"
67 #include "gromacs/mdtypes/commrec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/mdatom.h"
71 #include "gromacs/pbcutil/pbc.h"
72 #include "gromacs/pbcutil/pbc_simd.h"
73 #include "gromacs/simd/simd.h"
74 #include "gromacs/simd/simd_math.h"
75 #include "gromacs/simd/vector_operations.h"
76 #include "gromacs/topology/block.h"
77 #include "gromacs/topology/mtop_util.h"
78 #include "gromacs/utility/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/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(const rvec* x,
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 #if GMX_SIMD_HAVE_REAL
611 /* This SIMD code does the same as the plain-C code after the #else.
612 * The only difference is that we always call pbc code, as with SIMD
613 * the overhead of pbc computation (when not needed) is small.
615 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9 * GMX_SIMD_REAL_WIDTH];
617 /* Convert the pbc struct for SIMD */
618 set_pbc_simd(pbc, pbc_simd);
620 /* Compute normalized x i-j vectors, store in r.
621 * Compute the inner product of r and xp i-j and store in rhs1.
623 calc_dr_x_f_simd(b0, b1, atoms, x, f, blc.data(), pbc_simd, as_rvec_array(r.data()),
624 rhs1.data(), sol.data());
626 #else // GMX_SIMD_HAVE_REAL
628 /* Compute normalized i-j vectors */
631 for (int b = b0; b < b1; b++)
635 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
641 for (int b = b0; b < b1; b++)
645 rvec_sub(x[atoms[b].index1], x[atoms[b].index2], dx);
647 } /* 16 ncons flops */
650 for (int b = b0; b < b1; b++)
652 int i = atoms[b].index1;
653 int j = atoms[b].index2;
655 * (r[b][0] * (f[i][0] - f[j][0]) + r[b][1] * (f[i][1] - f[j][1])
656 + r[b][2] * (f[i][2] - f[j][2]));
662 #endif // GMX_SIMD_HAVE_REAL
664 if (lincsd->bTaskDep)
666 /* We need a barrier, since the matrix construction below
667 * can access entries in r of other threads.
672 /* Construct the (sparse) LINCS matrix */
673 for (int b = b0; b < b1; b++)
675 for (int n = blnr[b]; n < blnr[b + 1]; n++)
677 blcc[n] = blmf[n] * ::iprod(r[b], r[blbnb[n]]);
680 /* Together: 23*ncons + 6*nrtot flops */
682 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
683 /* nrec*(ncons+2*nrtot) flops */
685 if (econq == ConstraintVariable::Deriv_FlexCon)
687 /* We only want to constraint the flexible constraints,
688 * so we mask out the normal ones by setting sol to 0.
690 for (int b = b0; b < b1; b++)
692 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
699 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
700 for (int b = b0; b < b1; b++)
705 /* When constraining forces, we should not use mass weighting,
706 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
708 lincs_update_atoms(lincsd, th, 1.0, sol, r,
709 (econq != ConstraintVariable::Force) ? invmass : nullptr, fp);
714 for (int b = b0; b < b1; b++)
716 dhdlambda -= sol[b] * lincsd->ddist[b];
719 lincsd->task[th].dhdlambda = dhdlambda;
724 /* Constraint virial,
725 * determines sum r_bond x delta f,
726 * where delta f is the constraint correction
727 * of the quantity that is being constrained.
729 for (int b = b0; b < b1; b++)
731 const real mvb = lincsd->bllen[b] * sol[b];
732 for (int i = 0; i < DIM; i++)
734 const real tmp1 = mvb * r[b][i];
735 for (int j = 0; j < DIM; j++)
737 rmdf[i][j] += tmp1 * r[b][j];
740 } /* 23 ncons flops */
744 #if GMX_SIMD_HAVE_REAL
746 /*! \brief Calculate the constraint distance vectors r to project on from x.
748 * Determine the right-hand side of the matrix equation using coordinates xp. */
749 static void gmx_simdcall calc_dr_x_xp_simd(int b0,
751 gmx::ArrayRef<const AtomPair> atoms,
752 const rvec* gmx_restrict x,
753 const rvec* gmx_restrict xp,
754 const real* gmx_restrict bllen,
755 const real* gmx_restrict blc,
756 const real* pbc_simd,
757 rvec* gmx_restrict r,
758 real* gmx_restrict rhs,
759 real* gmx_restrict sol)
761 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
762 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
764 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
769 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
771 SimdReal x0_S, y0_S, z0_S;
772 SimdReal x1_S, y1_S, z1_S;
773 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
774 SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
775 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
776 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
778 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
780 offset0[i] = atoms[bs + i].index1;
781 offset1[i] = atoms[bs + i].index2;
784 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset0, &x0_S, &y0_S, &z0_S);
785 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset1, &x1_S, &y1_S, &z1_S);
790 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
792 n2_S = norm2(rx_S, ry_S, rz_S);
793 il_S = invsqrt(n2_S);
799 transposeScatterStoreU<3>(reinterpret_cast<real*>(r + bs), offset2, rx_S, ry_S, rz_S);
801 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(xp), offset0, &x0_S, &y0_S, &z0_S);
802 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(xp), offset1, &x1_S, &y1_S, &z1_S);
808 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
810 ip_S = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
812 rhs_S = load<SimdReal>(blc + bs) * (ip_S - load<SimdReal>(bllen + bs));
814 store(rhs + bs, rhs_S);
815 store(sol + bs, rhs_S);
818 #endif // GMX_SIMD_HAVE_REAL
820 /*! \brief Determine the distances and right-hand side for the next iteration. */
821 gmx_unused static void calc_dist_iter(int b0,
823 gmx::ArrayRef<const AtomPair> atoms,
824 const rvec* gmx_restrict xp,
825 const real* gmx_restrict bllen,
826 const real* gmx_restrict blc,
829 real* gmx_restrict rhs,
830 real* gmx_restrict sol,
833 for (int b = b0; b < b1; b++)
839 pbc_dx_aiuc(pbc, xp[atoms[b].index1], xp[atoms[b].index2], dx);
843 rvec_sub(xp[atoms[b].index1], xp[atoms[b].index2], dx);
845 real len2 = len * len;
846 real dlen2 = 2 * len2 - ::norm2(dx);
847 if (dlen2 < wfac * len2)
849 /* not race free - see detailed comment in caller */
855 mvb = blc[b] * (len - dlen2 * gmx::invsqrt(dlen2));
863 } /* 20*ncons flops */
866 #if GMX_SIMD_HAVE_REAL
867 /*! \brief As calc_dist_iter(), but using SIMD intrinsics. */
868 static void gmx_simdcall calc_dist_iter_simd(int b0,
870 gmx::ArrayRef<const AtomPair> atoms,
871 const rvec* gmx_restrict x,
872 const real* gmx_restrict bllen,
873 const real* gmx_restrict blc,
874 const real* pbc_simd,
876 real* gmx_restrict rhs,
877 real* gmx_restrict sol,
880 SimdReal min_S(GMX_REAL_MIN);
882 SimdReal wfac_S(wfac);
885 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
887 /* Initialize all to FALSE */
888 warn_S = (two_S < setZero());
890 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
892 SimdReal x0_S, y0_S, z0_S;
893 SimdReal x1_S, y1_S, z1_S;
894 SimdReal rx_S, ry_S, rz_S, n2_S;
895 SimdReal len_S, len2_S, dlen2_S, lc_S, blc_S;
896 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
897 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
899 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
901 offset0[i] = atoms[bs + i].index1;
902 offset1[i] = atoms[bs + i].index2;
905 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset0, &x0_S, &y0_S, &z0_S);
906 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real*>(x), offset1, &x1_S, &y1_S, &z1_S);
912 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
914 n2_S = norm2(rx_S, ry_S, rz_S);
916 len_S = load<SimdReal>(bllen + bs);
917 len2_S = len_S * len_S;
919 dlen2_S = fms(two_S, len2_S, n2_S);
921 warn_S = warn_S || (dlen2_S < (wfac_S * len2_S));
923 /* Avoid 1/0 by taking the max with REAL_MIN.
924 * Note: when dlen2 is close to zero (90 degree constraint rotation),
925 * the accuracy of the algorithm is no longer relevant.
927 dlen2_S = max(dlen2_S, min_S);
929 lc_S = fnma(dlen2_S, invsqrt(dlen2_S), len_S);
931 blc_S = load<SimdReal>(blc + bs);
935 store(rhs + bs, lc_S);
936 store(sol + bs, lc_S);
944 #endif // GMX_SIMD_HAVE_REAL
946 //! Implements LINCS constraining.
947 static void do_lincs(const rvec* x,
959 rvec* gmx_restrict v,
963 const int b0 = lincsd->task[th].b0;
964 const int b1 = lincsd->task[th].b1;
966 gmx::ArrayRef<const AtomPair> atoms = lincsd->atoms;
967 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
968 gmx::ArrayRef<const int> blnr = lincsd->blnr;
969 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
970 gmx::ArrayRef<const real> blc = lincsd->blc;
971 gmx::ArrayRef<const real> blmf = lincsd->blmf;
972 gmx::ArrayRef<const real> bllen = lincsd->bllen;
973 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
974 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
975 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
976 gmx::ArrayRef<real> sol = lincsd->tmp3;
977 gmx::ArrayRef<real> blc_sol = lincsd->tmp4;
978 gmx::ArrayRef<real> mlambda = lincsd->mlambda;
979 gmx::ArrayRef<const int> nlocat = lincsd->nlocat;
981 #if GMX_SIMD_HAVE_REAL
983 /* This SIMD code does the same as the plain-C code after the #else.
984 * The only difference is that we always call pbc code, as with SIMD
985 * the overhead of pbc computation (when not needed) is small.
987 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9 * GMX_SIMD_REAL_WIDTH];
989 /* Convert the pbc struct for SIMD */
990 set_pbc_simd(pbc, pbc_simd);
992 /* Compute normalized x i-j vectors, store in r.
993 * Compute the inner product of r and xp i-j and store in rhs1.
995 calc_dr_x_xp_simd(b0, b1, atoms, x, xp, bllen.data(), blc.data(), pbc_simd,
996 as_rvec_array(r.data()), rhs1.data(), sol.data());
998 #else // GMX_SIMD_HAVE_REAL
1002 /* Compute normalized i-j vectors */
1003 for (int b = b0; b < b1; b++)
1006 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
1009 pbc_dx_aiuc(pbc, xp[atoms[b].index1], xp[atoms[b].index2], dx);
1010 real mvb = blc[b] * (::iprod(r[b], dx) - bllen[b]);
1017 /* Compute normalized i-j vectors */
1018 for (int b = b0; b < b1; b++)
1020 int i = atoms[b].index1;
1021 int j = atoms[b].index2;
1022 real tmp0 = x[i][0] - x[j][0];
1023 real tmp1 = x[i][1] - x[j][1];
1024 real tmp2 = x[i][2] - x[j][2];
1025 real rlen = gmx::invsqrt(tmp0 * tmp0 + tmp1 * tmp1 + tmp2 * tmp2);
1026 r[b][0] = rlen * tmp0;
1027 r[b][1] = rlen * tmp1;
1028 r[b][2] = rlen * tmp2;
1029 /* 16 ncons flops */
1032 * (r[b][0] * (xp[i][0] - xp[j][0]) + r[b][1] * (xp[i][1] - xp[j][1])
1033 + r[b][2] * (xp[i][2] - xp[j][2]) - bllen[b]);
1038 /* Together: 26*ncons + 6*nrtot flops */
1041 #endif // GMX_SIMD_HAVE_REAL
1043 if (lincsd->bTaskDep)
1045 /* We need a barrier, since the matrix construction below
1046 * can access entries in r of other threads.
1051 /* Construct the (sparse) LINCS matrix */
1052 for (int b = b0; b < b1; b++)
1054 for (int n = blnr[b]; n < blnr[b + 1]; n++)
1056 blcc[n] = blmf[n] * gmx::dot(r[b], r[blbnb[n]]);
1059 /* Together: 26*ncons + 6*nrtot flops */
1061 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
1062 /* nrec*(ncons+2*nrtot) flops */
1064 #if GMX_SIMD_HAVE_REAL
1065 for (int b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1067 SimdReal t1 = load<SimdReal>(blc.data() + b);
1068 SimdReal t2 = load<SimdReal>(sol.data() + b);
1069 store(mlambda.data() + b, t1 * t2);
1072 for (int b = b0; b < b1; b++)
1074 mlambda[b] = blc[b] * sol[b];
1076 #endif // GMX_SIMD_HAVE_REAL
1078 /* Update the coordinates */
1079 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1082 ******** Correction for centripetal effects ********
1087 wfac = std::cos(DEG2RAD * wangle);
1090 for (int iter = 0; iter < lincsd->nIter; iter++)
1092 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
1097 /* Communicate the corrected non-local coordinates */
1098 if (DOMAINDECOMP(cr))
1100 dd_move_x_constraints(cr->dd, box, xp, nullptr, FALSE);
1105 else if (lincsd->bTaskDep)
1110 #if GMX_SIMD_HAVE_REAL
1111 calc_dist_iter_simd(b0, b1, atoms, xp, bllen.data(), blc.data(), pbc_simd, wfac,
1112 rhs1.data(), sol.data(), bWarn);
1114 calc_dist_iter(b0, b1, atoms, xp, bllen.data(), blc.data(), pbc, wfac, rhs1.data(),
1116 /* 20*ncons flops */
1117 #endif // GMX_SIMD_HAVE_REAL
1119 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
1120 /* nrec*(ncons+2*nrtot) flops */
1122 #if GMX_SIMD_HAVE_REAL
1123 for (int b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1125 SimdReal t1 = load<SimdReal>(blc.data() + b);
1126 SimdReal t2 = load<SimdReal>(sol.data() + b);
1127 SimdReal mvb = t1 * t2;
1128 store(blc_sol.data() + b, mvb);
1129 store(mlambda.data() + b, load<SimdReal>(mlambda.data() + b) + mvb);
1132 for (int b = b0; b < b1; b++)
1134 real mvb = blc[b] * sol[b];
1138 #endif // GMX_SIMD_HAVE_REAL
1140 /* Update the coordinates */
1141 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1143 /* nit*ncons*(37+9*nrec) flops */
1147 /* Update the velocities */
1148 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1149 /* 16 ncons flops */
1152 if (!nlocat.empty() && (bCalcDHDL || bCalcVir))
1154 if (lincsd->bTaskDep)
1156 /* In lincs_update_atoms threads might cross-read mlambda */
1160 /* Only account for local atoms */
1161 for (int b = b0; b < b1; b++)
1163 mlambda[b] *= 0.5 * nlocat[b];
1170 for (int b = b0; b < b1; b++)
1172 /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1173 * later after the contributions are reduced over the threads.
1175 dhdl -= lincsd->mlambda[b] * lincsd->ddist[b];
1177 lincsd->task[th].dhdlambda = dhdl;
1182 /* Constraint virial */
1183 for (int b = b0; b < b1; b++)
1185 real tmp0 = -bllen[b] * mlambda[b];
1186 for (int i = 0; i < DIM; i++)
1188 real tmp1 = tmp0 * r[b][i];
1189 for (int j = 0; j < DIM; j++)
1191 vir_r_m_dr[i][j] -= tmp1 * r[b][j];
1194 } /* 22 ncons flops */
1198 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1199 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1201 * (26+nrec)*ncons + (6+2*nrec)*nrtot
1202 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1204 * (63+nrec)*ncons + (6+4*nrec)*nrtot
1208 /*! \brief Sets the elements in the LINCS matrix for task task. */
1209 static void set_lincs_matrix_task(Lincs* li, Task* li_task, const real* invmass, int* ncc_triangle, int* nCrossTaskTriangles)
1211 /* Construct the coupling coefficient matrix blmf */
1212 li_task->ntriangle = 0;
1214 *nCrossTaskTriangles = 0;
1215 for (int i = li_task->b0; i < li_task->b1; i++)
1217 const int a1 = li->atoms[i].index1;
1218 const int a2 = li->atoms[i].index2;
1219 for (int n = li->blnr[i]; n < li->blnr[i + 1]; n++)
1221 const int k = li->blbnb[n];
1223 /* If we are using multiple, independent tasks for LINCS,
1224 * the calls to check_assign_connected should have
1225 * put all connected constraints in our task.
1227 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1230 if (a1 == li->atoms[k].index1 || a2 == li->atoms[k].index2)
1240 if (a1 == li->atoms[k].index1 || a1 == li->atoms[k].index2)
1250 li->blmf[n] = sign * invmass[center] * li->blc[i] * li->blc[k];
1251 li->blmf1[n] = sign * 0.5;
1252 if (li->ncg_triangle > 0)
1254 /* Look for constraint triangles */
1255 for (int nk = li->blnr[k]; nk < li->blnr[k + 1]; nk++)
1257 const int kk = li->blbnb[nk];
1258 if (kk != i && kk != k && (li->atoms[kk].index1 == end || li->atoms[kk].index2 == end))
1260 /* Check if the constraints in this triangle actually
1261 * belong to a different task. We still assign them
1262 * here, since it's convenient for the triangle
1263 * iterations, but we then need an extra barrier.
1265 if (k < li_task->b0 || k >= li_task->b1 || kk < li_task->b0 || kk >= li_task->b1)
1267 (*nCrossTaskTriangles)++;
1270 if (li_task->ntriangle == 0 || li_task->triangle[li_task->ntriangle - 1] < i)
1272 /* Add this constraint to the triangle list */
1273 li_task->triangle[li_task->ntriangle] = i;
1274 li_task->tri_bits[li_task->ntriangle] = 0;
1275 li_task->ntriangle++;
1276 if (li->blnr[i + 1] - li->blnr[i]
1277 > static_cast<int>(sizeof(li_task->tri_bits[0]) * 8 - 1))
1280 "A constraint is connected to %d constraints, this is "
1281 "more than the %zu allowed for constraints participating "
1283 li->blnr[i + 1] - li->blnr[i],
1284 sizeof(li_task->tri_bits[0]) * 8 - 1);
1287 li_task->tri_bits[li_task->ntriangle - 1] |= (1 << (n - li->blnr[i]));
1296 /*! \brief Sets the elements in the LINCS matrix. */
1297 static void set_lincs_matrix(Lincs* li, real* invmass, real lambda)
1299 const real invsqrt2 = 0.7071067811865475244;
1301 for (int i = 0; (i < li->nc); i++)
1303 const int a1 = li->atoms[i].index1;
1304 const int a2 = li->atoms[i].index2;
1305 li->blc[i] = gmx::invsqrt(invmass[a1] + invmass[a2]);
1306 li->blc1[i] = invsqrt2;
1309 /* Construct the coupling coefficient matrix blmf */
1310 int ntriangle = 0, ncc_triangle = 0, nCrossTaskTriangles = 0;
1311 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle, nCrossTaskTriangles) num_threads(li->ntask) schedule(static)
1312 for (int th = 0; th < li->ntask; th++)
1316 set_lincs_matrix_task(li, &li->task[th], invmass, &ncc_triangle, &nCrossTaskTriangles);
1317 ntriangle += li->task[th].ntriangle;
1319 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1321 li->ntriangle = ntriangle;
1322 li->ncc_triangle = ncc_triangle;
1323 li->bTaskDepTri = (nCrossTaskTriangles > 0);
1327 fprintf(debug, "The %d constraints participate in %d triangles\n", li->nc, li->ntriangle);
1328 fprintf(debug, "There are %d constraint couplings, of which %d in triangles\n", li->ncc,
1330 if (li->ntriangle > 0 && li->ntask > 1)
1333 "%d constraint triangles contain constraints assigned to different tasks\n",
1334 nCrossTaskTriangles);
1339 * so we know with which lambda value the masses have been set.
1341 li->matlam = lambda;
1344 //! Finds all triangles of atoms that share constraints to a central atom.
1345 static int count_triangle_constraints(const InteractionLists& ilist, const t_blocka& at2con)
1347 int ncon1, ncon_tot;
1348 int c0, n1, c1, ac1, n2, c2;
1351 ncon1 = ilist[F_CONSTR].size() / 3;
1352 ncon_tot = ncon1 + ilist[F_CONSTRNC].size() / 3;
1354 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1355 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1358 for (c0 = 0; c0 < ncon_tot; c0++)
1360 bool bTriangle = FALSE;
1361 const int* iap = constr_iatomptr(ia1, ia2, c0);
1362 const int a00 = iap[1];
1363 const int a01 = iap[2];
1364 for (n1 = at2con.index[a01]; n1 < at2con.index[a01 + 1]; n1++)
1369 const int* iap = constr_iatomptr(ia1, ia2, c1);
1370 const int a10 = iap[1];
1371 const int a11 = iap[2];
1380 for (n2 = at2con.index[ac1]; n2 < at2con.index[ac1 + 1]; n2++)
1383 if (c2 != c0 && c2 != c1)
1385 const int* iap = constr_iatomptr(ia1, ia2, c2);
1386 const int a20 = iap[1];
1387 const int a21 = iap[2];
1388 if (a20 == a00 || a21 == a00)
1402 return ncon_triangle;
1405 //! Finds sequences of sequential constraints.
1406 static bool more_than_two_sequential_constraints(const InteractionLists& ilist, const t_blocka& at2con)
1408 int ncon1, ncon_tot, c;
1409 bool bMoreThanTwoSequentialConstraints;
1411 ncon1 = ilist[F_CONSTR].size() / 3;
1412 ncon_tot = ncon1 + ilist[F_CONSTRNC].size() / 3;
1414 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1415 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1417 bMoreThanTwoSequentialConstraints = FALSE;
1418 for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
1420 const int* iap = constr_iatomptr(ia1, ia2, c);
1421 const int a1 = iap[1];
1422 const int a2 = iap[2];
1423 /* Check if this constraint has constraints connected at both atoms */
1424 if (at2con.index[a1 + 1] - at2con.index[a1] > 1 && at2con.index[a2 + 1] - at2con.index[a2] > 1)
1426 bMoreThanTwoSequentialConstraints = TRUE;
1430 return bMoreThanTwoSequentialConstraints;
1433 Lincs* init_lincs(FILE* fplog,
1434 const gmx_mtop_t& mtop,
1435 int nflexcon_global,
1436 ArrayRef<const t_blocka> at2con,
1441 // TODO this should become a unique_ptr
1443 bool bMoreThanTwoSeq;
1447 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n", bPLINCS ? " Parallel" : "");
1452 li->ncg = gmx_mtop_ftype_count(mtop, F_CONSTR) + gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1453 li->ncg_flex = nflexcon_global;
1456 li->nOrder = nProjOrder;
1458 li->max_connect = 0;
1459 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1461 for (int a = 0; a < mtop.moltype[mt].atoms.nr; a++)
1463 li->max_connect = std::max(li->max_connect, at2con[mt].index[a + 1] - at2con[mt].index[a]);
1467 li->ncg_triangle = 0;
1468 bMoreThanTwoSeq = FALSE;
1469 for (const gmx_molblock_t& molb : mtop.molblock)
1471 const gmx_moltype_t& molt = mtop.moltype[molb.type];
1473 li->ncg_triangle += molb.nmol * count_triangle_constraints(molt.ilist, at2con[molb.type]);
1475 if (!bMoreThanTwoSeq && more_than_two_sequential_constraints(molt.ilist, at2con[molb.type]))
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, int constraint_index, int a1, int a2, real lenA, real lenB, const t_blocka* 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->index[a1 + 1] - at2con->index[a1] - 1 + at2con->index[a2 + 1]
1668 - at2con->index[a2] - 1;
1670 li->blnr[con + 1] = li->ncc;
1672 /* Increase the constraint count */
1676 /*! \brief Check if constraint with topology index constraint_index is connected
1677 * to other constraints, and if so add those connected constraints to our task. */
1678 static void check_assign_connected(Lincs* li,
1679 const t_iatom* iatom,
1684 const t_blocka* at2con)
1686 /* Currently this function only supports constraint groups
1687 * in which all constraints share at least one atom
1688 * (e.g. H-bond constraints).
1689 * Check both ends of the current constraint for
1690 * connected constraints. We need to assign those
1695 for (end = 0; end < 2; end++)
1699 a = (end == 0 ? a1 : a2);
1701 for (k = at2con->index[a]; k < at2con->index[a + 1]; k++)
1706 /* Check if constraint cc has not yet been assigned */
1707 if (li->con_index[cc] == -1)
1712 type = iatom[cc * 3];
1713 lenA = idef.iparams[type].constr.dA;
1714 lenB = idef.iparams[type].constr.dB;
1716 if (bDynamics || lenA != 0 || lenB != 0)
1718 assign_constraint(li, cc, iatom[3 * cc + 1], iatom[3 * cc + 2], lenA, lenB, at2con);
1725 /*! \brief Check if constraint with topology index constraint_index is involved
1726 * in a constraint triangle, and if so add the other two constraints
1727 * in the triangle to our task. */
1728 static void check_assign_triangle(Lincs* li,
1729 const t_iatom* iatom,
1732 int constraint_index,
1735 const t_blocka* at2con)
1737 int nca, cc[32], ca[32], k;
1738 int c_triangle[2] = { -1, -1 };
1741 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1746 if (c != constraint_index)
1750 aa1 = iatom[c * 3 + 1];
1751 aa2 = iatom[c * 3 + 2];
1767 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1772 if (c != constraint_index)
1776 aa1 = iatom[c * 3 + 1];
1777 aa2 = iatom[c * 3 + 2];
1780 for (i = 0; i < nca; i++)
1784 c_triangle[0] = cc[i];
1791 for (i = 0; i < nca; i++)
1795 c_triangle[0] = cc[i];
1803 if (c_triangle[0] >= 0)
1807 for (end = 0; end < 2; end++)
1809 /* Check if constraint c_triangle[end] has not yet been assigned */
1810 if (li->con_index[c_triangle[end]] == -1)
1815 i = c_triangle[end] * 3;
1817 lenA = idef.iparams[type].constr.dA;
1818 lenB = idef.iparams[type].constr.dB;
1820 if (bDynamics || lenA != 0 || lenB != 0)
1822 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1829 //! Sets matrix indices.
1830 static void set_matrix_indices(Lincs* li, const Task& li_task, const t_blocka* at2con, bool bSortMatrix)
1832 for (int b = li_task.b0; b < li_task.b1; b++)
1834 const int a1 = li->atoms[b].index1;
1835 const int a2 = li->atoms[b].index2;
1837 int i = li->blnr[b];
1838 for (int k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1840 int concon = li->con_index[at2con->a[k]];
1843 li->blbnb[i++] = concon;
1846 for (int k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1848 int concon = li->con_index[at2con->a[k]];
1851 li->blbnb[i++] = concon;
1857 /* Order the blbnb matrix to optimize memory access */
1858 std::sort(li->blbnb.begin() + li->blnr[b], li->blbnb.begin() + li->blnr[b + 1]);
1863 void set_lincs(const t_idef& idef, const t_mdatoms& md, bool bDynamics, const t_commrec* cr, Lincs* li)
1872 /* Zero the thread index ranges.
1873 * Otherwise without local constraints we could return with old ranges.
1875 for (int i = 0; i < li->ntask; i++)
1879 li->task[i].ind.clear();
1883 li->task[li->ntask].ind.clear();
1886 /* This is the local topology, so there are only F_CONSTR constraints */
1887 if (idef.il[F_CONSTR].nr == 0)
1889 /* There are no constraints,
1890 * we do not need to fill any data structures.
1897 fprintf(debug, "Building the LINCS connectivity\n");
1900 if (DOMAINDECOMP(cr))
1902 if (cr->dd->constraints)
1906 dd_get_constraint_range(cr->dd, &start, &natoms);
1910 natoms = dd_numHomeAtoms(*cr->dd);
1918 at2con = make_at2con(natoms, idef.il, idef.iparams, flexibleConstraintTreatment(bDynamics));
1920 const int ncon_tot = idef.il[F_CONSTR].nr / 3;
1922 /* Ensure we have enough padding for aligned loads for each thread */
1923 const int numEntries = ncon_tot + li->ntask * simd_width;
1924 li->con_index.resize(numEntries);
1925 li->bllen0.resize(numEntries);
1926 li->ddist.resize(numEntries);
1927 li->atoms.resize(numEntries);
1928 li->blc.resize(numEntries);
1929 li->blc1.resize(numEntries);
1930 li->blnr.resize(numEntries + 1);
1931 li->bllen.resize(numEntries);
1932 li->tmpv.resizeWithPadding(numEntries);
1933 if (DOMAINDECOMP(cr))
1935 li->nlocat.resize(numEntries);
1937 li->tmp1.resize(numEntries);
1938 li->tmp2.resize(numEntries);
1939 li->tmp3.resize(numEntries);
1940 li->tmp4.resize(numEntries);
1941 li->mlambda.resize(numEntries);
1943 iatom = idef.il[F_CONSTR].iatoms;
1945 li->blnr[0] = li->ncc;
1947 /* Assign the constraints for li->ntask LINCS tasks.
1948 * We target a uniform distribution of constraints over the tasks.
1949 * Note that when flexible constraints are present, but are removed here
1950 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
1951 * happen during normal MD, that's ok.
1954 /* Determine the number of constraints we need to assign here */
1955 int ncon_assign = ncon_tot;
1958 /* With energy minimization, flexible constraints are ignored
1959 * (and thus minimized, as they should be).
1961 ncon_assign -= countFlexibleConstraints(idef.il, idef.iparams);
1964 /* Set the target constraint count per task to exactly uniform,
1965 * this might be overridden below.
1967 int ncon_target = (ncon_assign + li->ntask - 1) / li->ntask;
1969 /* Mark all constraints as unassigned by setting their index to -1 */
1970 for (int con = 0; con < ncon_tot; con++)
1972 li->con_index[con] = -1;
1976 for (int th = 0; th < li->ntask; th++)
1980 li_task = &li->task[th];
1982 #if GMX_SIMD_HAVE_REAL
1983 /* With indepedent tasks we likely have H-bond constraints or constraint
1984 * pairs. The connected constraints will be pulled into the task, so the
1985 * constraints per task will often exceed ncon_target.
1986 * Triangle constraints can also increase the count, but there are
1987 * relatively few of those, so we usually expect to get ncon_target.
1991 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
1992 * since otherwise a lot of operations can be wasted.
1993 * There are several ways to round here, we choose the one
1994 * that alternates block sizes, which helps with Intel HT.
1996 ncon_target = ((ncon_assign * (th + 1)) / li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1)
1997 & ~(GMX_SIMD_REAL_WIDTH - 1);
1999 #endif // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
2001 /* Continue filling the arrays where we left off with the previous task,
2002 * including padding for SIMD.
2004 li_task->b0 = li->nc;
2006 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2008 if (li->con_index[con] == -1)
2013 type = iatom[3 * con];
2014 a1 = iatom[3 * con + 1];
2015 a2 = iatom[3 * con + 2];
2016 lenA = idef.iparams[type].constr.dA;
2017 lenB = idef.iparams[type].constr.dB;
2018 /* Skip the flexible constraints when not doing dynamics */
2019 if (bDynamics || lenA != 0 || lenB != 0)
2021 assign_constraint(li, con, a1, a2, lenA, lenB, &at2con);
2023 if (li->ntask > 1 && !li->bTaskDep)
2025 /* We can generate independent tasks. Check if we
2026 * need to assign connected constraints to our task.
2028 check_assign_connected(li, iatom, idef, bDynamics, a1, a2, &at2con);
2030 if (li->ntask > 1 && li->ncg_triangle > 0)
2032 /* Ensure constraints in one triangle are assigned
2035 check_assign_triangle(li, iatom, idef, bDynamics, con, a1, a2, &at2con);
2043 li_task->b1 = li->nc;
2047 /* Copy the last atom pair indices and lengths for constraints
2048 * up to a multiple of simd_width, such that we can do all
2049 * SIMD operations without having to worry about end effects.
2053 li->nc = ((li_task->b1 + simd_width - 1) / simd_width) * simd_width;
2054 last = li_task->b1 - 1;
2055 for (i = li_task->b1; i < li->nc; i++)
2057 li->atoms[i] = li->atoms[last];
2058 li->bllen0[i] = li->bllen0[last];
2059 li->ddist[i] = li->ddist[last];
2060 li->bllen[i] = li->bllen[last];
2061 li->blnr[i + 1] = li->blnr[last + 1];
2065 /* Keep track of how many constraints we assigned */
2066 li->nc_real += li_task->b1 - li_task->b0;
2070 fprintf(debug, "LINCS task %d constraints %d - %d\n", th, li_task->b0, li_task->b1);
2074 assert(li->nc_real == ncon_assign);
2078 /* Without DD we order the blbnb matrix to optimize memory access.
2079 * With DD the overhead of sorting is more than the gain during access.
2081 bSortMatrix = !DOMAINDECOMP(cr);
2083 li->blbnb.resize(li->ncc);
2085 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2086 for (int th = 0; th < li->ntask; th++)
2090 Task& li_task = li->task[th];
2092 if (li->ncg_triangle > 0)
2094 /* This is allocating too much, but it is difficult to improve */
2095 li_task.triangle.resize(li_task.b1 - li_task.b0);
2096 li_task.tri_bits.resize(li_task.b1 - li_task.b0);
2099 set_matrix_indices(li, li_task, &at2con, bSortMatrix);
2101 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2104 done_blocka(&at2con);
2106 if (cr->dd == nullptr)
2108 /* Since the matrix is static, we should free some memory */
2109 li->blbnb.resize(li->ncc);
2112 li->blmf.resize(li->ncc);
2113 li->blmf1.resize(li->ncc);
2114 li->tmpncc.resize(li->ncc);
2116 gmx::ArrayRef<const int> nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2117 if (!nlocat_dd.empty())
2119 /* Convert nlocat from local topology to LINCS constraint indexing */
2120 for (con = 0; con < ncon_tot; con++)
2122 li->nlocat[li->con_index[con]] = nlocat_dd[con];
2132 fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n", li->nc_real,
2138 lincs_thread_setup(li, md.nr);
2141 set_lincs_matrix(li, md.invmass, md.lambda);
2144 //! Issues a warning when LINCS constraints cannot be satisfied.
2145 static void lincs_warning(gmx_domdec_t* dd,
2150 gmx::ArrayRef<const AtomPair> atoms,
2151 gmx::ArrayRef<const real> bllen,
2156 real wfac = std::cos(DEG2RAD * wangle);
2159 "bonds that rotated more than %g degrees:\n"
2160 " atom 1 atom 2 angle previous, current, constraint length\n",
2163 for (int b = 0; b < ncons; b++)
2165 const int i = atoms[b].index1;
2166 const int j = atoms[b].index2;
2171 pbc_dx_aiuc(pbc, x[i], x[j], v0);
2172 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2176 rvec_sub(x[i], x[j], v0);
2177 rvec_sub(xprime[i], xprime[j], v1);
2181 real cosine = ::iprod(v0, v1) / (d0 * d1);
2184 fprintf(stderr, " %6d %6d %5.1f %8.4f %8.4f %8.4f\n", ddglatnr(dd, i),
2185 ddglatnr(dd, j), RAD2DEG * std::acos(cosine), d0, d1, bllen[b]);
2186 if (!std::isfinite(d1))
2188 gmx_fatal(FARGS, "Bond length not finite.");
2194 if (*warncount > maxwarn)
2196 too_many_constraint_warnings(econtLINCS, *warncount);
2200 //! Status information about how well LINCS satisified the constraints in this domain
2201 struct LincsDeviations
2203 //! The maximum over all bonds in this domain of the relative deviation in bond lengths
2204 real maxDeviation = 0;
2205 //! Sum over all bonds in this domain of the squared relative deviation
2206 real sumSquaredDeviation = 0;
2207 //! Index of bond with max deviation
2208 int indexOfMaxDeviation = -1;
2209 //! Number of bonds constrained in this domain
2210 int numConstraints = 0;
2213 //! Determine how well the constraints have been satisfied.
2214 static LincsDeviations makeLincsDeviations(const Lincs& lincsd, const rvec* x, const t_pbc* pbc)
2216 LincsDeviations result;
2217 const ArrayRef<const AtomPair> atoms = lincsd.atoms;
2218 const ArrayRef<const real> bllen = lincsd.bllen;
2219 const ArrayRef<const int> nlocat = lincsd.nlocat;
2221 for (int task = 0; task < lincsd.ntask; task++)
2223 for (int b = lincsd.task[task].b0; b < lincsd.task[task].b1; b++)
2228 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
2232 rvec_sub(x[atoms[b].index1], x[atoms[b].index2], dx);
2234 real r2 = ::norm2(dx);
2235 real len = r2 * gmx::invsqrt(r2);
2236 real d = std::abs(len / bllen[b] - 1.0_real);
2237 if (d > result.maxDeviation && (nlocat.empty() || nlocat[b]))
2239 result.maxDeviation = d;
2240 result.indexOfMaxDeviation = b;
2244 result.sumSquaredDeviation += d * d;
2245 result.numConstraints++;
2249 result.sumSquaredDeviation += nlocat[b] * d * d;
2250 result.numConstraints += nlocat[b];
2255 if (!nlocat.empty())
2257 result.numConstraints /= 2;
2258 result.sumSquaredDeviation *= 0.5;
2263 bool constrain_lincs(bool computeRmsd,
2264 const t_inputrec& ir,
2267 const t_mdatoms& md,
2268 const t_commrec* cr,
2269 const gmx_multisim_t* ms,
2281 ConstraintVariable econq,
2288 /* This boolean should be set by a flag passed to this routine.
2289 * We can also easily check if any constraint length is changed,
2290 * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2292 bool bCalcDHDL = (ir.efep != efepNO && dvdlambda != nullptr);
2294 if (lincsd->nc == 0 && cr->dd == nullptr)
2298 lincsd->rmsdData = { { 0 } };
2304 if (econq == ConstraintVariable::Positions)
2306 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2307 * also with efep!=fepNO.
2309 if (ir.efep != efepNO)
2311 if (md.nMassPerturbed && lincsd->matlam != md.lambda)
2313 set_lincs_matrix(lincsd, md.invmass, md.lambda);
2316 for (int i = 0; i < lincsd->nc; i++)
2318 lincsd->bllen[i] = lincsd->bllen0[i] + lambda * lincsd->ddist[i];
2322 if (lincsd->ncg_flex)
2324 /* Set the flexible constraint lengths to the old lengths */
2327 for (int i = 0; i < lincsd->nc; i++)
2329 if (lincsd->bllen[i] == 0)
2332 pbc_dx_aiuc(pbc, x[lincsd->atoms[i].index1], x[lincsd->atoms[i].index2], dx);
2333 lincsd->bllen[i] = norm(dx);
2339 for (int i = 0; i < lincsd->nc; i++)
2341 if (lincsd->bllen[i] == 0)
2343 lincsd->bllen[i] = std::sqrt(
2344 distance2(x[lincsd->atoms[i].index1], x[lincsd->atoms[i].index2]));
2350 const bool printDebugOutput = ((debug != nullptr) && lincsd->nc > 0);
2351 if (printDebugOutput)
2353 LincsDeviations deviations = makeLincsDeviations(*lincsd, xprime, pbc);
2354 fprintf(debug, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2355 fprintf(debug, " Before LINCS %.6f %.6f %6d %6d\n",
2356 std::sqrt(deviations.sumSquaredDeviation / deviations.numConstraints),
2357 deviations.maxDeviation,
2358 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index1),
2359 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index2));
2362 /* This bWarn var can be updated by multiple threads
2363 * at the same time. But as we only need to detect
2364 * if a warning occurred or not, this is not an issue.
2368 /* The OpenMP parallel region of constrain_lincs for coords */
2369 #pragma omp parallel num_threads(lincsd->ntask)
2373 int th = gmx_omp_get_thread_num();
2375 clear_mat(lincsd->task[th].vir_r_m_dr);
2377 do_lincs(x, xprime, box, pbc, lincsd, th, md.invmass, cr, bCalcDHDL, ir.LincsWarnAngle,
2378 &bWarn, invdt, v, bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2380 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2383 if (computeRmsd || printDebugOutput || bWarn)
2385 LincsDeviations deviations = makeLincsDeviations(*lincsd, xprime, pbc);
2389 // This is reduced across domains in compute_globals and
2390 // reported to the log file.
2391 lincsd->rmsdData[0] = deviations.numConstraints;
2392 lincsd->rmsdData[1] = deviations.sumSquaredDeviation;
2396 // This is never read
2397 lincsd->rmsdData = { { 0 } };
2399 if (printDebugOutput)
2401 fprintf(debug, " After LINCS %.6f %.6f %6d %6d\n\n",
2402 std::sqrt(deviations.sumSquaredDeviation / deviations.numConstraints),
2403 deviations.maxDeviation,
2404 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index1),
2405 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index2));
2410 if (maxwarn < INT_MAX)
2412 std::string simMesg;
2415 simMesg += gmx::formatString(" in simulation %d", ms->sim);
2419 ", time %g (ps) LINCS WARNING%s\n"
2420 "relative constraint deviation after LINCS:\n"
2421 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2422 step, ir.init_t + step * ir.delta_t, simMesg.c_str(),
2423 std::sqrt(deviations.sumSquaredDeviation / deviations.numConstraints),
2424 deviations.maxDeviation,
2425 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index1),
2426 ddglatnr(cr->dd, lincsd->atoms[deviations.indexOfMaxDeviation].index2));
2428 lincs_warning(cr->dd, x, xprime, pbc, lincsd->nc, lincsd->atoms, lincsd->bllen,
2429 ir.LincsWarnAngle, maxwarn, warncount);
2431 bOK = (deviations.maxDeviation < 0.5);
2435 if (lincsd->ncg_flex)
2437 for (int i = 0; (i < lincsd->nc); i++)
2439 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2441 lincsd->bllen[i] = 0;
2448 /* The OpenMP parallel region of constrain_lincs for derivatives */
2449 #pragma omp parallel num_threads(lincsd->ntask)
2453 int th = gmx_omp_get_thread_num();
2455 do_lincsp(x, xprime, min_proj, pbc, lincsd, th, md.invmass, econq, bCalcDHDL,
2456 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2458 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2464 /* Reduce the dH/dlambda contributions over the threads */
2469 for (th = 0; th < lincsd->ntask; th++)
2471 dhdlambda += lincsd->task[th].dhdlambda;
2473 if (econq == ConstraintVariable::Positions)
2475 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2476 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2477 dhdlambda /= ir.delta_t * ir.delta_t;
2479 *dvdlambda += dhdlambda;
2482 if (bCalcVir && lincsd->ntask > 1)
2484 for (int i = 1; i < lincsd->ntask; i++)
2486 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2490 /* count assuming nit=1 */
2491 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2492 inc_nrnb(nrnb, eNR_LINCSMAT, (2 + lincsd->nOrder) * lincsd->ncc);
2493 if (lincsd->ntriangle > 0)
2495 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder * lincsd->ncc_triangle);
2499 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real * 2);
2503 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);