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, 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/mdlib/mdrun.h"
67 #include "gromacs/mdtypes/commrec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/mdatom.h"
71 #include "gromacs/pbcutil/pbc.h"
72 #include "gromacs/pbcutil/pbc-simd.h"
73 #include "gromacs/simd/simd.h"
74 #include "gromacs/simd/simd_math.h"
75 #include "gromacs/simd/vector_operations.h"
76 #include "gromacs/topology/block.h"
77 #include "gromacs/topology/mtop_util.h"
78 #include "gromacs/utility/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
107 //! Unit of work within LINCS.
110 //! First constraint for this task.
112 //! b1-1 is the last constraint for this task.
114 //! The number of constraints in triangles.
116 //! The list of triangle constraints.
117 std::vector<int> triangle;
118 //! The bits tell if the matrix element should be used.
119 std::vector<int> tri_bits;
120 //! Constraint index for updating atom data.
121 std::vector<int> ind;
122 //! Constraint index for updating atom data.
123 std::vector<int> ind_r;
124 //! Temporary variable for virial calculation.
125 tensor vir_r_m_dr = {{0}};
126 //! Temporary variable for lambda derivative.
130 /*! \brief Data for LINCS algorithm.
135 //! The global number of constraints.
137 //! The global number of flexible constraints.
139 //! The global number of constraints in triangles.
140 int ncg_triangle = 0;
141 //! The number of iterations.
143 //! The order of the matrix expansion.
145 //! The maximum number of constraints connected to a single atom.
148 //! The number of real constraints.
150 //! The number of constraints including padding for SIMD.
152 //! The number of constraint connections.
154 //! The FE lambda value used for filling blc and blmf.
156 //! mapping from topology to LINCS constraints.
157 std::vector<int> con_index;
158 //! The reference distance in topology A.
159 std::vector < real, AlignedAllocator < real>> bllen0;
160 //! The reference distance in top B - the r.d. in top A.
161 std::vector < real, AlignedAllocator < real>> ddist;
162 //! The atom pairs involved in the constraints.
163 std::vector<AtomPair> atoms;
164 //! 1/sqrt(invmass1 invmass2).
165 std::vector < real, AlignedAllocator < real>> blc;
166 //! As blc, but with all masses 1.
167 std::vector < real, AlignedAllocator < real>> blc1;
168 //! Index into blbnb and blmf.
169 std::vector<int> blnr;
170 //! List of constraint connections.
171 std::vector<int> blbnb;
172 //! The local number of constraints in triangles.
174 //! The number of constraint connections in triangles.
175 int ncc_triangle = 0;
176 //! Communicate before each LINCS interation.
177 bool bCommIter = false;
178 //! Matrix of mass factors for constraint connections.
179 std::vector<real> blmf;
180 //! As blmf, but with all masses 1.
181 std::vector<real> blmf1;
182 //! The reference bond length.
183 std::vector < real, AlignedAllocator < real>> bllen;
184 //! The local atom count per constraint, can be NULL.
185 std::vector<int> nlocat;
187 /*! \brief The number of tasks used for LINCS work.
189 * \todo This is mostly used to loop over \c task, which would
190 * be nicer to do with range-based for loops, but the thread
191 * index is used for constructing bit masks and organizing the
192 * virial output buffer, so other things need to change,
195 /*! \brief LINCS thread division */
196 std::vector<Task> task;
197 //! Atom flags for thread parallelization.
198 std::vector<gmx_bitmask_t> atf;
199 //! Are the LINCS tasks interdependent?
200 bool bTaskDep = false;
201 //! Are there triangle constraints that cross task borders?
202 bool bTaskDepTri = false;
203 //! Arrays for temporary storage in the LINCS algorithm.
205 PaddedVector<gmx::RVec> tmpv;
206 std::vector<real> tmpncc;
207 std::vector < real, AlignedAllocator < real>> tmp1;
208 std::vector < real, AlignedAllocator < real>> tmp2;
209 std::vector < real, AlignedAllocator < real>> tmp3;
210 std::vector < real, AlignedAllocator < real>> tmp4;
212 //! The Lagrange multipliers times -1.
213 std::vector < real, AlignedAllocator < real>> mlambda;
214 //! Storage for the constraint RMS relative deviation output.
215 std::array<real, 2> rmsdData = {{0}};
218 /*! \brief Define simd_width for memory allocation used for SIMD code */
219 #if GMX_SIMD_HAVE_REAL
220 static const int simd_width = GMX_SIMD_REAL_WIDTH;
222 static const int simd_width = 1;
225 ArrayRef<real> lincs_rmsdData(Lincs *lincsd)
227 return lincsd->rmsdData;
230 real lincs_rmsd(const Lincs *lincsd)
232 if (lincsd->rmsdData[0] > 0)
234 return std::sqrt(lincsd->rmsdData[1]/lincsd->rmsdData[0]);
242 /*! \brief Do a set of nrec LINCS matrix multiplications.
244 * This function will return with up to date thread-local
245 * constraint data, without an OpenMP barrier.
247 static void lincs_matrix_expand(const Lincs &lincsd,
249 gmx::ArrayRef<const real> blcc,
250 gmx::ArrayRef<real> rhs1,
251 gmx::ArrayRef<real> rhs2,
252 gmx::ArrayRef<real> sol)
254 gmx::ArrayRef<const int> blnr = lincsd.blnr;
255 gmx::ArrayRef<const int> blbnb = lincsd.blbnb;
257 const int b0 = li_task.b0;
258 const int b1 = li_task.b1;
259 const int nrec = lincsd.nOrder;
261 for (int rec = 0; rec < nrec; rec++)
267 for (int b = b0; b < b1; b++)
273 for (n = blnr[b]; n < blnr[b+1]; n++)
275 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
278 sol[b] = sol[b] + mvb;
281 gmx::ArrayRef<real> swap;
286 } /* nrec*(ncons+2*nrtot) flops */
288 if (lincsd.ntriangle > 0)
290 /* Perform an extra nrec recursions for only the constraints
291 * involved in rigid triangles.
292 * In this way their accuracy should come close to those of the other
293 * constraints, since traingles of constraints can produce eigenvalues
294 * around 0.7, while the effective eigenvalue for bond constraints
295 * is around 0.4 (and 0.7*0.7=0.5).
300 /* We need a barrier here, since other threads might still be
301 * reading the contents of rhs1 and/o rhs2.
302 * We could avoid this barrier by introducing two extra rhs
303 * arrays for the triangle constraints only.
308 /* Constraints involved in a triangle are ensured to be in the same
309 * LINCS task. This means no barriers are required during the extra
310 * iterations for the triangle constraints.
312 gmx::ArrayRef<const int> triangle = li_task.triangle;
313 gmx::ArrayRef<const int> tri_bits = li_task.tri_bits;
315 for (int rec = 0; rec < nrec; rec++)
317 for (int tb = 0; tb < li_task.ntriangle; tb++)
319 int b, bits, nr0, nr1, n;
327 for (n = nr0; n < nr1; n++)
329 if (bits & (1 << (n - nr0)))
331 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
335 sol[b] = sol[b] + mvb;
338 gmx::ArrayRef<real> swap;
343 } /* nrec*(ntriangle + ncc_triangle*2) flops */
345 if (lincsd.bTaskDepTri)
347 /* The constraints triangles are decoupled from each other,
348 * but constraints in one triangle cross thread task borders.
349 * We could probably avoid this with more advanced setup code.
356 //! Update atomic coordinates when an index is not required.
358 lincs_update_atoms_noind(int ncons,
359 gmx::ArrayRef<const AtomPair> atoms,
361 gmx::ArrayRef<const real> fac,
362 gmx::ArrayRef<const gmx::RVec> r,
366 if (invmass != nullptr)
368 for (int b = 0; b < ncons; b++)
370 int i = atoms[b].index1;
371 int j = atoms[b].index2;
372 real mvb = preFactor*fac[b];
373 real im1 = invmass[i];
374 real im2 = invmass[j];
375 real tmp0 = r[b][0]*mvb;
376 real tmp1 = r[b][1]*mvb;
377 real tmp2 = r[b][2]*mvb;
384 } /* 16 ncons flops */
388 for (int b = 0; b < ncons; b++)
390 int i = atoms[b].index1;
391 int j = atoms[b].index2;
392 real mvb = preFactor*fac[b];
393 real tmp0 = r[b][0]*mvb;
394 real tmp1 = r[b][1]*mvb;
395 real tmp2 = r[b][2]*mvb;
406 //! Update atomic coordinates when an index is required.
408 lincs_update_atoms_ind(gmx::ArrayRef<const int> ind,
409 gmx::ArrayRef<const AtomPair> atoms,
411 gmx::ArrayRef<const real> fac,
412 gmx::ArrayRef<const gmx::RVec> r,
416 if (invmass != nullptr)
420 int i = atoms[b].index1;
421 int j = atoms[b].index2;
422 real mvb = preFactor*fac[b];
423 real im1 = invmass[i];
424 real im2 = invmass[j];
425 real tmp0 = r[b][0]*mvb;
426 real tmp1 = r[b][1]*mvb;
427 real tmp2 = r[b][2]*mvb;
434 } /* 16 ncons flops */
440 int i = atoms[b].index1;
441 int j = atoms[b].index2;
442 real mvb = preFactor*fac[b];
443 real tmp0 = r[b][0]*mvb;
444 real tmp1 = r[b][1]*mvb;
445 real tmp2 = r[b][2]*mvb;
452 } /* 16 ncons flops */
456 //! Update coordinates for atoms.
458 lincs_update_atoms(Lincs *li,
461 gmx::ArrayRef<const real> fac,
462 gmx::ArrayRef<const gmx::RVec> r,
468 /* Single thread, we simply update for all constraints */
469 lincs_update_atoms_noind(li->nc_real,
470 li->atoms, preFactor, fac, r, invmass, x);
474 /* Update the atom vector components for our thread local
475 * constraints that only access our local atom range.
476 * This can be done without a barrier.
478 lincs_update_atoms_ind(li->task[th].ind,
479 li->atoms, preFactor, fac, r, invmass, x);
481 if (!li->task[li->ntask].ind.empty())
483 /* Update the constraints that operate on atoms
484 * in multiple thread atom blocks on the master thread.
489 lincs_update_atoms_ind(li->task[li->ntask].ind,
490 li->atoms, preFactor, fac, r, invmass, x);
496 #if GMX_SIMD_HAVE_REAL
497 //! Helper function so that we can run TSAN with SIMD support (where implemented).
499 static inline void gmx_simdcall
500 gatherLoadUTransposeTSANSafe(const real *base,
501 const std::int32_t *offset,
506 #if (CMAKE_BUILD_TYPE == CMAKE_BUILD_TYPE_TSAN) && GMX_SIMD_X86_AVX2_256
507 // This function is only implemented in this case
508 gatherLoadUTransposeSafe<align>(base, offset, v0, v1, v2);
510 gatherLoadUTranspose<align>(base, offset, v0, v1, v2);
514 /*! \brief Calculate the constraint distance vectors r to project on from x.
516 * Determine the right-hand side of the matrix equation using quantity f.
517 * This function only differs from calc_dr_x_xp_simd below in that
518 * no constraint length is subtracted and no PBC is used for f. */
519 static void gmx_simdcall
520 calc_dr_x_f_simd(int b0,
522 gmx::ArrayRef<const AtomPair> atoms,
523 const rvec * gmx_restrict x,
524 const rvec * gmx_restrict f,
525 const real * gmx_restrict blc,
526 const real * pbc_simd,
527 rvec * gmx_restrict r,
528 real * gmx_restrict rhs,
529 real * gmx_restrict sol)
531 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
533 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
535 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
540 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
542 SimdReal x0_S, y0_S, z0_S;
543 SimdReal x1_S, y1_S, z1_S;
544 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
545 SimdReal fx_S, fy_S, fz_S, ip_S, rhs_S;
546 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
547 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
549 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
551 offset0[i] = atoms[bs + i].index1;
552 offset1[i] = atoms[bs + i].index2;
555 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
556 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
561 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
563 n2_S = norm2(rx_S, ry_S, rz_S);
564 il_S = invsqrt(n2_S);
570 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
572 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(f), offset0, &x0_S, &y0_S, &z0_S);
573 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(f), offset1, &x1_S, &y1_S, &z1_S);
578 ip_S = iprod(rx_S, ry_S, rz_S, fx_S, fy_S, fz_S);
580 rhs_S = load<SimdReal>(blc + bs) * ip_S;
582 store(rhs + bs, rhs_S);
583 store(sol + bs, rhs_S);
586 #endif // GMX_SIMD_HAVE_REAL
588 /*! \brief LINCS projection, works on derivatives of the coordinates. */
589 static void do_lincsp(const rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
590 Lincs *lincsd, int th,
592 ConstraintVariable econq, bool bCalcDHDL,
593 bool bCalcVir, tensor rmdf)
595 const int b0 = lincsd->task[th].b0;
596 const int b1 = lincsd->task[th].b1;
598 gmx::ArrayRef<const AtomPair> atoms = lincsd->atoms;
599 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
600 gmx::ArrayRef<const int> blnr = lincsd->blnr;
601 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
603 gmx::ArrayRef<const real> blc;
604 gmx::ArrayRef<const real> blmf;
605 if (econq != ConstraintVariable::Force)
607 /* Use mass-weighted parameters */
613 /* Use non mass-weighted parameters */
615 blmf = lincsd->blmf1;
617 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
618 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
619 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
620 gmx::ArrayRef<real> sol = lincsd->tmp3;
622 #if GMX_SIMD_HAVE_REAL
623 /* This SIMD code does the same as the plain-C code after the #else.
624 * The only difference is that we always call pbc code, as with SIMD
625 * the overhead of pbc computation (when not needed) is small.
627 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
629 /* Convert the pbc struct for SIMD */
630 set_pbc_simd(pbc, pbc_simd);
632 /* Compute normalized x i-j vectors, store in r.
633 * Compute the inner product of r and xp i-j and store in rhs1.
635 calc_dr_x_f_simd(b0, b1, atoms, x, f, blc.data(),
637 as_rvec_array(r.data()), rhs1.data(), sol.data());
639 #else // GMX_SIMD_HAVE_REAL
641 /* Compute normalized i-j vectors */
644 for (int b = b0; b < b1; b++)
648 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
654 for (int b = b0; b < b1; b++)
658 rvec_sub(x[atoms[b].index1], x[atoms[b].index2], dx);
660 } /* 16 ncons flops */
663 for (int b = b0; b < b1; b++)
665 int i = atoms[b].index1;
666 int j = atoms[b].index2;
667 real mvb = blc[b]*(r[b][0]*(f[i][0] - f[j][0]) +
668 r[b][1]*(f[i][1] - f[j][1]) +
669 r[b][2]*(f[i][2] - f[j][2]));
675 #endif // GMX_SIMD_HAVE_REAL
677 if (lincsd->bTaskDep)
679 /* We need a barrier, since the matrix construction below
680 * can access entries in r of other threads.
685 /* Construct the (sparse) LINCS matrix */
686 for (int b = b0; b < b1; b++)
688 for (int n = blnr[b]; n < blnr[b+1]; n++)
690 blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
693 /* Together: 23*ncons + 6*nrtot flops */
695 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
696 /* nrec*(ncons+2*nrtot) flops */
698 if (econq == ConstraintVariable::Deriv_FlexCon)
700 /* We only want to constraint the flexible constraints,
701 * so we mask out the normal ones by setting sol to 0.
703 for (int b = b0; b < b1; b++)
705 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
712 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
713 for (int b = b0; b < b1; b++)
718 /* When constraining forces, we should not use mass weighting,
719 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
721 lincs_update_atoms(lincsd, th, 1.0, sol, r,
722 (econq != ConstraintVariable::Force) ? invmass : nullptr, fp);
727 for (int b = b0; b < b1; b++)
729 dhdlambda -= sol[b]*lincsd->ddist[b];
732 lincsd->task[th].dhdlambda = dhdlambda;
737 /* Constraint virial,
738 * determines sum r_bond x delta f,
739 * where delta f is the constraint correction
740 * of the quantity that is being constrained.
742 for (int b = b0; b < b1; b++)
744 const real mvb = lincsd->bllen[b]*sol[b];
745 for (int i = 0; i < DIM; i++)
747 const real tmp1 = mvb*r[b][i];
748 for (int j = 0; j < DIM; j++)
750 rmdf[i][j] += tmp1*r[b][j];
753 } /* 23 ncons flops */
757 #if GMX_SIMD_HAVE_REAL
759 /*! \brief Calculate the constraint distance vectors r to project on from x.
761 * Determine the right-hand side of the matrix equation using coordinates xp. */
762 static void gmx_simdcall
763 calc_dr_x_xp_simd(int b0,
765 gmx::ArrayRef<const AtomPair> atoms,
766 const rvec * gmx_restrict x,
767 const rvec * gmx_restrict xp,
768 const real * gmx_restrict bllen,
769 const real * gmx_restrict blc,
770 const real * pbc_simd,
771 rvec * gmx_restrict r,
772 real * gmx_restrict rhs,
773 real * gmx_restrict sol)
775 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
776 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
778 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
783 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
785 SimdReal x0_S, y0_S, z0_S;
786 SimdReal x1_S, y1_S, z1_S;
787 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
788 SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
789 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
790 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
792 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
794 offset0[i] = atoms[bs + i].index1;
795 offset1[i] = atoms[bs + i].index2;
798 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
799 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
804 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
806 n2_S = norm2(rx_S, ry_S, rz_S);
807 il_S = invsqrt(n2_S);
813 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
815 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(xp), offset0, &x0_S, &y0_S, &z0_S);
816 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(xp), offset1, &x1_S, &y1_S, &z1_S);
822 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
824 ip_S = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
826 rhs_S = load<SimdReal>(blc + bs) * (ip_S - load<SimdReal>(bllen + bs));
828 store(rhs + bs, rhs_S);
829 store(sol + bs, rhs_S);
832 #endif // GMX_SIMD_HAVE_REAL
834 /*! \brief Determine the distances and right-hand side for the next iteration. */
835 gmx_unused static void calc_dist_iter(
838 gmx::ArrayRef<const AtomPair> atoms,
839 const rvec * gmx_restrict xp,
840 const real * gmx_restrict bllen,
841 const real * gmx_restrict blc,
844 real * gmx_restrict rhs,
845 real * gmx_restrict sol,
848 for (int b = b0; b < b1; b++)
854 pbc_dx_aiuc(pbc, xp[atoms[b].index1], xp[atoms[b].index2], dx);
858 rvec_sub(xp[atoms[b].index1], xp[atoms[b].index2], dx);
861 real dlen2 = 2*len2 - ::norm2(dx);
862 if (dlen2 < wfac*len2)
864 /* not race free - see detailed comment in caller */
870 mvb = blc[b]*(len - dlen2*gmx::invsqrt(dlen2));
878 } /* 20*ncons flops */
881 #if GMX_SIMD_HAVE_REAL
882 /*! \brief As calc_dist_iter(), but using SIMD intrinsics. */
883 static void gmx_simdcall
884 calc_dist_iter_simd(int b0,
886 gmx::ArrayRef<const AtomPair> atoms,
887 const rvec * gmx_restrict x,
888 const real * gmx_restrict bllen,
889 const real * gmx_restrict blc,
890 const real * pbc_simd,
892 real * gmx_restrict rhs,
893 real * gmx_restrict sol,
896 SimdReal min_S(GMX_REAL_MIN);
898 SimdReal wfac_S(wfac);
901 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
903 /* Initialize all to FALSE */
904 warn_S = (two_S < setZero());
906 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
908 SimdReal x0_S, y0_S, z0_S;
909 SimdReal x1_S, y1_S, z1_S;
910 SimdReal rx_S, ry_S, rz_S, n2_S;
911 SimdReal len_S, len2_S, dlen2_S, lc_S, blc_S;
912 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
913 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
915 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
917 offset0[i] = atoms[bs + i].index1;
918 offset1[i] = atoms[bs + i].index2;
921 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
922 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
928 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
930 n2_S = norm2(rx_S, ry_S, rz_S);
932 len_S = load<SimdReal>(bllen + bs);
933 len2_S = len_S * len_S;
935 dlen2_S = fms(two_S, len2_S, n2_S);
937 warn_S = warn_S || (dlen2_S < (wfac_S * len2_S));
939 /* Avoid 1/0 by taking the max with REAL_MIN.
940 * Note: when dlen2 is close to zero (90 degree constraint rotation),
941 * the accuracy of the algorithm is no longer relevant.
943 dlen2_S = max(dlen2_S, min_S);
945 lc_S = fnma(dlen2_S, invsqrt(dlen2_S), len_S);
947 blc_S = load<SimdReal>(blc + bs);
951 store(rhs + bs, lc_S);
952 store(sol + bs, lc_S);
960 #endif // GMX_SIMD_HAVE_REAL
962 //! Implements LINCS constraining.
963 static void do_lincs(const rvec *x, rvec *xp, matrix box, t_pbc *pbc,
964 Lincs *lincsd, int th,
968 real wangle, bool *bWarn,
969 real invdt, rvec * gmx_restrict v,
970 bool bCalcVir, tensor vir_r_m_dr)
972 const int b0 = lincsd->task[th].b0;
973 const int b1 = lincsd->task[th].b1;
975 gmx::ArrayRef<const AtomPair> atoms = lincsd->atoms;
976 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
977 gmx::ArrayRef<const int> blnr = lincsd->blnr;
978 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
979 gmx::ArrayRef<const real> blc = lincsd->blc;
980 gmx::ArrayRef<const real> blmf = lincsd->blmf;
981 gmx::ArrayRef<const real> bllen = lincsd->bllen;
982 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
983 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
984 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
985 gmx::ArrayRef<real> sol = lincsd->tmp3;
986 gmx::ArrayRef<real> blc_sol = lincsd->tmp4;
987 gmx::ArrayRef<real> mlambda = lincsd->mlambda;
988 gmx::ArrayRef<const int> nlocat = lincsd->nlocat;
990 #if GMX_SIMD_HAVE_REAL
992 /* This SIMD code does the same as the plain-C code after the #else.
993 * The only difference is that we always call pbc code, as with SIMD
994 * the overhead of pbc computation (when not needed) is small.
996 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
998 /* Convert the pbc struct for SIMD */
999 set_pbc_simd(pbc, pbc_simd);
1001 /* Compute normalized x i-j vectors, store in r.
1002 * Compute the inner product of r and xp i-j and store in rhs1.
1004 calc_dr_x_xp_simd(b0, b1, atoms, x, xp, bllen.data(), blc.data(),
1006 as_rvec_array(r.data()), rhs1.data(), sol.data());
1008 #else // GMX_SIMD_HAVE_REAL
1012 /* Compute normalized i-j vectors */
1013 for (int b = b0; b < b1; b++)
1016 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
1019 pbc_dx_aiuc(pbc, xp[atoms[b].index1], xp[atoms[b].index2], dx);
1020 real mvb = blc[b]*(::iprod(r[b], dx) - bllen[b]);
1027 /* Compute normalized i-j vectors */
1028 for (int b = b0; b < b1; b++)
1030 int i = atoms[b].index1;
1031 int j = atoms[b].index2;
1032 real tmp0 = x[i][0] - x[j][0];
1033 real tmp1 = x[i][1] - x[j][1];
1034 real tmp2 = x[i][2] - x[j][2];
1035 real rlen = gmx::invsqrt(tmp0*tmp0 + tmp1*tmp1 + tmp2*tmp2);
1036 r[b][0] = rlen*tmp0;
1037 r[b][1] = rlen*tmp1;
1038 r[b][2] = rlen*tmp2;
1039 /* 16 ncons flops */
1041 real mvb = blc[b]*(r[b][0]*(xp[i][0] - xp[j][0]) +
1042 r[b][1]*(xp[i][1] - xp[j][1]) +
1043 r[b][2]*(xp[i][2] - xp[j][2]) - bllen[b]);
1048 /* Together: 26*ncons + 6*nrtot flops */
1051 #endif // GMX_SIMD_HAVE_REAL
1053 if (lincsd->bTaskDep)
1055 /* We need a barrier, since the matrix construction below
1056 * can access entries in r of other threads.
1061 /* Construct the (sparse) LINCS matrix */
1062 for (int b = b0; b < b1; b++)
1064 for (int n = blnr[b]; n < blnr[b+1]; n++)
1066 blcc[n] = blmf[n]*gmx::dot(r[b], r[blbnb[n]]);
1069 /* Together: 26*ncons + 6*nrtot flops */
1071 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
1072 /* nrec*(ncons+2*nrtot) flops */
1074 #if GMX_SIMD_HAVE_REAL
1075 for (int b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1077 SimdReal t1 = load<SimdReal>(blc.data() + b);
1078 SimdReal t2 = load<SimdReal>(sol.data() + b);
1079 store(mlambda.data() + b, t1 * t2);
1082 for (int b = b0; b < b1; b++)
1084 mlambda[b] = blc[b]*sol[b];
1086 #endif // GMX_SIMD_HAVE_REAL
1088 /* Update the coordinates */
1089 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1092 ******** Correction for centripetal effects ********
1097 wfac = std::cos(DEG2RAD*wangle);
1100 for (int iter = 0; iter < lincsd->nIter; iter++)
1102 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
1107 /* Communicate the corrected non-local coordinates */
1108 if (DOMAINDECOMP(cr))
1110 dd_move_x_constraints(cr->dd, box, xp, nullptr, FALSE);
1115 else if (lincsd->bTaskDep)
1120 #if GMX_SIMD_HAVE_REAL
1121 calc_dist_iter_simd(b0, b1, atoms,
1122 xp, bllen.data(), blc.data(), pbc_simd, wfac,
1123 rhs1.data(), sol.data(), bWarn);
1125 calc_dist_iter(b0, b1, atoms, xp,
1126 bllen.data(), blc.data(), pbc, wfac,
1127 rhs1.data(), sol.data(), bWarn);
1128 /* 20*ncons flops */
1129 #endif // GMX_SIMD_HAVE_REAL
1131 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
1132 /* nrec*(ncons+2*nrtot) flops */
1134 #if GMX_SIMD_HAVE_REAL
1135 for (int b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1137 SimdReal t1 = load<SimdReal>(blc.data() + b);
1138 SimdReal t2 = load<SimdReal>(sol.data() + b);
1139 SimdReal mvb = t1 * t2;
1140 store(blc_sol.data() + b, mvb);
1141 store(mlambda.data() + b, load<SimdReal>(mlambda.data() + b) + mvb);
1144 for (int b = b0; b < b1; b++)
1146 real mvb = blc[b]*sol[b];
1150 #endif // GMX_SIMD_HAVE_REAL
1152 /* Update the coordinates */
1153 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1155 /* nit*ncons*(37+9*nrec) flops */
1159 /* Update the velocities */
1160 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1161 /* 16 ncons flops */
1164 if (!nlocat.empty() && (bCalcDHDL || bCalcVir))
1166 if (lincsd->bTaskDep)
1168 /* In lincs_update_atoms threads might cross-read mlambda */
1172 /* Only account for local atoms */
1173 for (int b = b0; b < b1; b++)
1175 mlambda[b] *= 0.5*nlocat[b];
1182 for (int b = b0; b < b1; b++)
1184 /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1185 * later after the contributions are reduced over the threads.
1187 dhdl -= lincsd->mlambda[b]*lincsd->ddist[b];
1189 lincsd->task[th].dhdlambda = dhdl;
1194 /* Constraint virial */
1195 for (int b = b0; b < b1; b++)
1197 real tmp0 = -bllen[b]*mlambda[b];
1198 for (int i = 0; i < DIM; i++)
1200 real tmp1 = tmp0*r[b][i];
1201 for (int j = 0; j < DIM; j++)
1203 vir_r_m_dr[i][j] -= tmp1*r[b][j];
1206 } /* 22 ncons flops */
1210 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1211 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1213 * (26+nrec)*ncons + (6+2*nrec)*nrtot
1214 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1216 * (63+nrec)*ncons + (6+4*nrec)*nrtot
1220 /*! \brief Sets the elements in the LINCS matrix for task task. */
1221 static void set_lincs_matrix_task(Lincs *li,
1223 const real *invmass,
1225 int *nCrossTaskTriangles)
1227 /* Construct the coupling coefficient matrix blmf */
1228 li_task->ntriangle = 0;
1230 *nCrossTaskTriangles = 0;
1231 for (int i = li_task->b0; i < li_task->b1; i++)
1233 const int a1 = li->atoms[i].index1;
1234 const int a2 = li->atoms[i].index2;
1235 for (int n = li->blnr[i]; n < li->blnr[i + 1]; n++)
1237 const int k = li->blbnb[n];
1239 /* If we are using multiple, independent tasks for LINCS,
1240 * the calls to check_assign_connected should have
1241 * put all connected constraints in our task.
1243 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1246 if (a1 == li->atoms[k].index1 || a2 == li->atoms[k].index2)
1256 if (a1 == li->atoms[k].index1 || a1 == li->atoms[k].index2)
1266 li->blmf[n] = sign*invmass[center]*li->blc[i]*li->blc[k];
1267 li->blmf1[n] = sign*0.5;
1268 if (li->ncg_triangle > 0)
1270 /* Look for constraint triangles */
1271 for (int nk = li->blnr[k]; nk < li->blnr[k + 1]; nk++)
1273 const int kk = li->blbnb[nk];
1274 if (kk != i && kk != k &&
1275 (li->atoms[kk].index1 == end ||
1276 li->atoms[kk].index2 == end))
1278 /* Check if the constraints in this triangle actually
1279 * belong to a different task. We still assign them
1280 * here, since it's convenient for the triangle
1281 * iterations, but we then need an extra barrier.
1283 if (k < li_task->b0 || k >= li_task->b1 ||
1284 kk < li_task->b0 || kk >= li_task->b1)
1286 (*nCrossTaskTriangles)++;
1289 if (li_task->ntriangle == 0 ||
1290 li_task->triangle[li_task->ntriangle - 1] < i)
1292 /* Add this constraint to the triangle list */
1293 li_task->triangle[li_task->ntriangle] = i;
1294 li_task->tri_bits[li_task->ntriangle] = 0;
1295 li_task->ntriangle++;
1296 if (li->blnr[i+1] - li->blnr[i] > static_cast<int>(sizeof(li_task->tri_bits[0])*8 - 1))
1298 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %zu allowed for constraints participating in triangles",
1299 li->blnr[i+1] - li->blnr[i],
1300 sizeof(li_task->tri_bits[0])*8-1);
1303 li_task->tri_bits[li_task->ntriangle-1] |= (1 << (n - li->blnr[i]));
1312 /*! \brief Sets the elements in the LINCS matrix. */
1313 static void set_lincs_matrix(Lincs *li, real *invmass, real lambda)
1315 const real invsqrt2 = 0.7071067811865475244;
1317 for (int i = 0; (i < li->nc); i++)
1319 const int a1 = li->atoms[i].index1;
1320 const int a2 = li->atoms[i].index2;
1321 li->blc[i] = gmx::invsqrt(invmass[a1] + invmass[a2]);
1322 li->blc1[i] = invsqrt2;
1325 /* Construct the coupling coefficient matrix blmf */
1326 int ntriangle = 0, ncc_triangle = 0, nCrossTaskTriangles = 0;
1327 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle, nCrossTaskTriangles) num_threads(li->ntask) schedule(static)
1328 for (int th = 0; th < li->ntask; th++)
1332 set_lincs_matrix_task(li, &li->task[th], invmass,
1333 &ncc_triangle, &nCrossTaskTriangles);
1334 ntriangle += li->task[th].ntriangle;
1336 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1338 li->ntriangle = ntriangle;
1339 li->ncc_triangle = ncc_triangle;
1340 li->bTaskDepTri = (nCrossTaskTriangles > 0);
1344 fprintf(debug, "The %d constraints participate in %d triangles\n",
1345 li->nc, li->ntriangle);
1346 fprintf(debug, "There are %d constraint couplings, of which %d in triangles\n",
1347 li->ncc, li->ncc_triangle);
1348 if (li->ntriangle > 0 && li->ntask > 1)
1350 fprintf(debug, "%d constraint triangles contain constraints assigned to different tasks\n",
1351 nCrossTaskTriangles);
1356 * so we know with which lambda value the masses have been set.
1358 li->matlam = lambda;
1361 //! Finds all triangles of atoms that share constraints to a central atom.
1362 static int count_triangle_constraints(const InteractionLists &ilist,
1363 const t_blocka &at2con)
1365 int ncon1, ncon_tot;
1366 int c0, n1, c1, ac1, n2, c2;
1369 ncon1 = ilist[F_CONSTR].size()/3;
1370 ncon_tot = ncon1 + ilist[F_CONSTRNC].size()/3;
1372 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1373 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1376 for (c0 = 0; c0 < ncon_tot; c0++)
1378 bool bTriangle = FALSE;
1379 const int *iap = constr_iatomptr(ia1, ia2, c0);
1380 const int a00 = iap[1];
1381 const int a01 = iap[2];
1382 for (n1 = at2con.index[a01]; n1 < at2con.index[a01+1]; n1++)
1387 const int *iap = constr_iatomptr(ia1, ia2, c1);
1388 const int a10 = iap[1];
1389 const int a11 = iap[2];
1398 for (n2 = at2con.index[ac1]; n2 < at2con.index[ac1+1]; n2++)
1401 if (c2 != c0 && c2 != c1)
1403 const int *iap = constr_iatomptr(ia1, ia2, c2);
1404 const int a20 = iap[1];
1405 const int a21 = iap[2];
1406 if (a20 == a00 || a21 == a00)
1420 return ncon_triangle;
1423 //! Finds sequences of sequential constraints.
1424 static bool more_than_two_sequential_constraints(const InteractionLists &ilist,
1425 const t_blocka &at2con)
1427 int ncon1, ncon_tot, c;
1428 bool bMoreThanTwoSequentialConstraints;
1430 ncon1 = ilist[F_CONSTR].size()/3;
1431 ncon_tot = ncon1 + ilist[F_CONSTRNC].size()/3;
1433 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1434 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1436 bMoreThanTwoSequentialConstraints = FALSE;
1437 for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
1439 const int *iap = constr_iatomptr(ia1, ia2, c);
1440 const int a1 = iap[1];
1441 const int a2 = iap[2];
1442 /* Check if this constraint has constraints connected at both atoms */
1443 if (at2con.index[a1+1] - at2con.index[a1] > 1 &&
1444 at2con.index[a2+1] - at2con.index[a2] > 1)
1446 bMoreThanTwoSequentialConstraints = TRUE;
1450 return bMoreThanTwoSequentialConstraints;
1453 Lincs *init_lincs(FILE *fplog, const gmx_mtop_t &mtop,
1454 int nflexcon_global, ArrayRef<const t_blocka> at2con,
1455 bool bPLINCS, int nIter, int nProjOrder)
1457 // TODO this should become a unique_ptr
1459 bool bMoreThanTwoSeq;
1463 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
1464 bPLINCS ? " Parallel" : "");
1470 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1471 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1472 li->ncg_flex = nflexcon_global;
1475 li->nOrder = nProjOrder;
1477 li->max_connect = 0;
1478 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1480 for (int a = 0; a < mtop.moltype[mt].atoms.nr; a++)
1482 li->max_connect = std::max(li->max_connect,
1483 at2con[mt].index[a + 1] - at2con[mt].index[a]);
1487 li->ncg_triangle = 0;
1488 bMoreThanTwoSeq = FALSE;
1489 for (const gmx_molblock_t &molb : mtop.molblock)
1491 const gmx_moltype_t &molt = mtop.moltype[molb.type];
1495 count_triangle_constraints(molt.ilist, at2con[molb.type]);
1497 if (!bMoreThanTwoSeq &&
1498 more_than_two_sequential_constraints(molt.ilist, at2con[molb.type]))
1500 bMoreThanTwoSeq = TRUE;
1504 /* Check if we need to communicate not only before LINCS,
1505 * but also before each iteration.
1506 * The check for only two sequential constraints is only
1507 * useful for the common case of H-bond only constraints.
1508 * With more effort we could also make it useful for small
1509 * molecules with nr. sequential constraints <= nOrder-1.
1511 li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1513 if (debug && bPLINCS)
1515 fprintf(debug, "PLINCS communication before each iteration: %d\n",
1516 static_cast<int>(li->bCommIter));
1519 /* LINCS can run on any number of threads.
1520 * Currently the number is fixed for the whole simulation,
1521 * but it could be set in set_lincs().
1522 * The current constraint to task assignment code can create independent
1523 * tasks only when not more than two constraints are connected sequentially.
1525 li->ntask = gmx_omp_nthreads_get(emntLINCS);
1526 li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1529 fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n",
1530 li->ntask, li->bTaskDep ? "" : "in");
1538 /* Allocate an extra elements for "task-overlap" constraints */
1539 li->task.resize(li->ntask + 1);
1542 if (bPLINCS || li->ncg_triangle > 0)
1544 please_cite(fplog, "Hess2008a");
1548 please_cite(fplog, "Hess97a");
1553 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1556 fprintf(fplog, "There are constraints between atoms in different decomposition domains,\n"
1557 "will communicate selected coordinates each lincs iteration\n");
1559 if (li->ncg_triangle > 0)
1562 "%d constraints are involved in constraint triangles,\n"
1563 "will apply an additional matrix expansion of order %d for couplings\n"
1564 "between constraints inside triangles\n",
1565 li->ncg_triangle, li->nOrder);
1572 void done_lincs(Lincs *li)
1577 /*! \brief Sets up the work division over the threads. */
1578 static void lincs_thread_setup(Lincs *li, int natoms)
1580 li->atf.resize(natoms);
1582 gmx::ArrayRef<gmx_bitmask_t> atf = li->atf;
1584 /* Clear the atom flags */
1585 for (gmx_bitmask_t &mask : atf)
1587 bitmask_clear(&mask);
1590 if (li->ntask > BITMASK_SIZE)
1592 gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1595 for (int th = 0; th < li->ntask; th++)
1597 const Task *li_task = &li->task[th];
1599 /* For each atom set a flag for constraints from each */
1600 for (int b = li_task->b0; b < li_task->b1; b++)
1602 bitmask_set_bit(&atf[li->atoms[b].index1], th);
1603 bitmask_set_bit(&atf[li->atoms[b].index2], th);
1607 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1608 for (int th = 0; th < li->ntask; th++)
1616 li_task = &li->task[th];
1618 bitmask_init_low_bits(&mask, th);
1620 li_task->ind.clear();
1621 li_task->ind_r.clear();
1622 for (b = li_task->b0; b < li_task->b1; b++)
1624 /* We let the constraint with the lowest thread index
1625 * operate on atoms with constraints from multiple threads.
1627 if (bitmask_is_disjoint(atf[li->atoms[b].index1], mask) &&
1628 bitmask_is_disjoint(atf[li->atoms[b].index2], mask))
1630 /* Add the constraint to the local atom update index */
1631 li_task->ind.push_back(b);
1635 /* Add the constraint to the rest block */
1636 li_task->ind_r.push_back(b);
1640 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1643 /* We need to copy all constraints which have not be assigned
1644 * to a thread to a separate list which will be handled by one thread.
1646 Task *li_m = &li->task[li->ntask];
1649 for (int th = 0; th < li->ntask; th++)
1651 const Task &li_task = li->task[th];
1653 for (int ind_r : li_task.ind_r)
1655 li_m->ind.push_back(ind_r);
1660 fprintf(debug, "LINCS thread %d: %zu constraints\n",
1661 th, li_task.ind.size());
1667 fprintf(debug, "LINCS thread r: %zu constraints\n",
1672 //! Assign a constraint.
1673 static void assign_constraint(Lincs *li,
1674 int constraint_index,
1676 real lenA, real lenB,
1677 const t_blocka *at2con)
1683 /* Make an mapping of local topology constraint index to LINCS index */
1684 li->con_index[constraint_index] = con;
1686 li->bllen0[con] = lenA;
1687 li->ddist[con] = lenB - lenA;
1688 /* Set the length to the topology A length */
1689 li->bllen[con] = lenA;
1690 li->atoms[con].index1 = a1;
1691 li->atoms[con].index2 = a2;
1693 /* Make space in the constraint connection matrix for constraints
1694 * connected to both end of the current constraint.
1697 at2con->index[a1 + 1] - at2con->index[a1] - 1 +
1698 at2con->index[a2 + 1] - at2con->index[a2] - 1;
1700 li->blnr[con + 1] = li->ncc;
1702 /* Increase the constraint count */
1706 /*! \brief Check if constraint with topology index constraint_index is connected
1707 * to other constraints, and if so add those connected constraints to our task. */
1708 static void check_assign_connected(Lincs *li,
1709 const t_iatom *iatom,
1713 const t_blocka *at2con)
1715 /* Currently this function only supports constraint groups
1716 * in which all constraints share at least one atom
1717 * (e.g. H-bond constraints).
1718 * Check both ends of the current constraint for
1719 * connected constraints. We need to assign those
1724 for (end = 0; end < 2; end++)
1728 a = (end == 0 ? a1 : a2);
1730 for (k = at2con->index[a]; k < at2con->index[a + 1]; k++)
1735 /* Check if constraint cc has not yet been assigned */
1736 if (li->con_index[cc] == -1)
1742 lenA = idef.iparams[type].constr.dA;
1743 lenB = idef.iparams[type].constr.dB;
1745 if (bDynamics || lenA != 0 || lenB != 0)
1747 assign_constraint(li, cc, iatom[3*cc + 1], iatom[3*cc + 2], lenA, lenB, at2con);
1754 /*! \brief Check if constraint with topology index constraint_index is involved
1755 * in a constraint triangle, and if so add the other two constraints
1756 * in the triangle to our task. */
1757 static void check_assign_triangle(Lincs *li,
1758 const t_iatom *iatom,
1761 int constraint_index,
1763 const t_blocka *at2con)
1765 int nca, cc[32], ca[32], k;
1766 int c_triangle[2] = { -1, -1 };
1769 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1774 if (c != constraint_index)
1778 aa1 = iatom[c*3 + 1];
1779 aa2 = iatom[c*3 + 2];
1795 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1800 if (c != constraint_index)
1804 aa1 = iatom[c*3 + 1];
1805 aa2 = iatom[c*3 + 2];
1808 for (i = 0; i < nca; i++)
1812 c_triangle[0] = cc[i];
1819 for (i = 0; i < nca; i++)
1823 c_triangle[0] = cc[i];
1831 if (c_triangle[0] >= 0)
1835 for (end = 0; end < 2; end++)
1837 /* Check if constraint c_triangle[end] has not yet been assigned */
1838 if (li->con_index[c_triangle[end]] == -1)
1843 i = c_triangle[end]*3;
1845 lenA = idef.iparams[type].constr.dA;
1846 lenB = idef.iparams[type].constr.dB;
1848 if (bDynamics || lenA != 0 || lenB != 0)
1850 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1857 //! Sets matrix indices.
1858 static void set_matrix_indices(Lincs *li,
1859 const Task &li_task,
1860 const t_blocka *at2con,
1863 for (int b = li_task.b0; b < li_task.b1; b++)
1865 const int a1 = li->atoms[b].index1;
1866 const int a2 = li->atoms[b].index2;
1868 int i = li->blnr[b];
1869 for (int k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1871 int concon = li->con_index[at2con->a[k]];
1874 li->blbnb[i++] = concon;
1877 for (int k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1879 int concon = li->con_index[at2con->a[k]];
1882 li->blbnb[i++] = concon;
1888 /* Order the blbnb matrix to optimize memory access */
1889 std::sort(&(li->blbnb[li->blnr[b]]), &(li->blbnb[li->blnr[b+1]]));
1894 void set_lincs(const t_idef &idef,
1895 const t_mdatoms &md,
1897 const t_commrec *cr,
1907 /* Zero the thread index ranges.
1908 * Otherwise without local constraints we could return with old ranges.
1910 for (int i = 0; i < li->ntask; i++)
1914 li->task[i].ind.clear();
1918 li->task[li->ntask].ind.clear();
1921 /* This is the local topology, so there are only F_CONSTR constraints */
1922 if (idef.il[F_CONSTR].nr == 0)
1924 /* There are no constraints,
1925 * we do not need to fill any data structures.
1932 fprintf(debug, "Building the LINCS connectivity\n");
1935 if (DOMAINDECOMP(cr))
1937 if (cr->dd->constraints)
1941 dd_get_constraint_range(cr->dd, &start, &natoms);
1945 natoms = dd_numHomeAtoms(*cr->dd);
1953 at2con = make_at2con(natoms, idef.il, idef.iparams,
1954 flexibleConstraintTreatment(bDynamics));
1956 const int ncon_tot = idef.il[F_CONSTR].nr/3;
1958 /* Ensure we have enough padding for aligned loads for each thread */
1959 const int numEntries = ncon_tot + li->ntask*simd_width;
1960 li->con_index.resize(numEntries);
1961 li->bllen0.resize(numEntries);
1962 li->ddist.resize(numEntries);
1963 li->atoms.resize(numEntries);
1964 li->blc.resize(numEntries);
1965 li->blc1.resize(numEntries);
1966 li->blnr.resize(numEntries);
1967 li->bllen.resize(numEntries);
1968 li->tmpv.resizeWithPadding(numEntries);
1969 if (DOMAINDECOMP(cr))
1971 li->nlocat.resize(numEntries);
1973 li->tmp1.resize(numEntries);
1974 li->tmp2.resize(numEntries);
1975 li->tmp3.resize(numEntries);
1976 li->tmp4.resize(numEntries);
1977 li->mlambda.resize(numEntries);
1979 iatom = idef.il[F_CONSTR].iatoms;
1981 li->blnr[0] = li->ncc;
1983 /* Assign the constraints for li->ntask LINCS tasks.
1984 * We target a uniform distribution of constraints over the tasks.
1985 * Note that when flexible constraints are present, but are removed here
1986 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
1987 * happen during normal MD, that's ok.
1990 /* Determine the number of constraints we need to assign here */
1991 int ncon_assign = ncon_tot;
1994 /* With energy minimization, flexible constraints are ignored
1995 * (and thus minimized, as they should be).
1997 ncon_assign -= countFlexibleConstraints(idef.il, idef.iparams);
2000 /* Set the target constraint count per task to exactly uniform,
2001 * this might be overridden below.
2003 int ncon_target = (ncon_assign + li->ntask - 1)/li->ntask;
2005 /* Mark all constraints as unassigned by setting their index to -1 */
2006 for (int con = 0; con < ncon_tot; con++)
2008 li->con_index[con] = -1;
2012 for (int th = 0; th < li->ntask; th++)
2016 li_task = &li->task[th];
2018 #if GMX_SIMD_HAVE_REAL
2019 /* With indepedent tasks we likely have H-bond constraints or constraint
2020 * pairs. The connected constraints will be pulled into the task, so the
2021 * constraints per task will often exceed ncon_target.
2022 * Triangle constraints can also increase the count, but there are
2023 * relatively few of those, so we usually expect to get ncon_target.
2027 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
2028 * since otherwise a lot of operations can be wasted.
2029 * There are several ways to round here, we choose the one
2030 * that alternates block sizes, which helps with Intel HT.
2032 ncon_target = ((ncon_assign*(th + 1))/li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1) & ~(GMX_SIMD_REAL_WIDTH - 1);
2034 #endif // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
2036 /* Continue filling the arrays where we left off with the previous task,
2037 * including padding for SIMD.
2039 li_task->b0 = li->nc;
2041 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2043 if (li->con_index[con] == -1)
2048 type = iatom[3*con];
2049 a1 = iatom[3*con + 1];
2050 a2 = iatom[3*con + 2];
2051 lenA = idef.iparams[type].constr.dA;
2052 lenB = idef.iparams[type].constr.dB;
2053 /* Skip the flexible constraints when not doing dynamics */
2054 if (bDynamics || lenA != 0 || lenB != 0)
2056 assign_constraint(li, con, a1, a2, lenA, lenB, &at2con);
2058 if (li->ntask > 1 && !li->bTaskDep)
2060 /* We can generate independent tasks. Check if we
2061 * need to assign connected constraints to our task.
2063 check_assign_connected(li, iatom, idef, bDynamics,
2066 if (li->ntask > 1 && li->ncg_triangle > 0)
2068 /* Ensure constraints in one triangle are assigned
2071 check_assign_triangle(li, iatom, idef, bDynamics,
2072 con, a1, a2, &at2con);
2080 li_task->b1 = li->nc;
2084 /* Copy the last atom pair indices and lengths for constraints
2085 * up to a multiple of simd_width, such that we can do all
2086 * SIMD operations without having to worry about end effects.
2090 li->nc = ((li_task->b1 + simd_width - 1)/simd_width)*simd_width;
2091 last = li_task->b1 - 1;
2092 for (i = li_task->b1; i < li->nc; i++)
2094 li->atoms[i] = li->atoms[last];
2095 li->bllen0[i] = li->bllen0[last];
2096 li->ddist[i] = li->ddist[last];
2097 li->bllen[i] = li->bllen[last];
2098 li->blnr[i + 1] = li->blnr[last + 1];
2102 /* Keep track of how many constraints we assigned */
2103 li->nc_real += li_task->b1 - li_task->b0;
2107 fprintf(debug, "LINCS task %d constraints %d - %d\n",
2108 th, li_task->b0, li_task->b1);
2112 assert(li->nc_real == ncon_assign);
2116 /* Without DD we order the blbnb matrix to optimize memory access.
2117 * With DD the overhead of sorting is more than the gain during access.
2119 bSortMatrix = !DOMAINDECOMP(cr);
2121 li->blbnb.resize(li->ncc);
2123 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2124 for (int th = 0; th < li->ntask; th++)
2128 Task &li_task = li->task[th];
2130 if (li->ncg_triangle > 0)
2132 /* This is allocating too much, but it is difficult to improve */
2133 li_task.triangle.resize(li_task.b1 - li_task.b0);
2134 li_task.tri_bits.resize(li_task.b1 - li_task.b0);
2137 set_matrix_indices(li, li_task, &at2con, bSortMatrix);
2139 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2142 done_blocka(&at2con);
2144 if (cr->dd == nullptr)
2146 /* Since the matrix is static, we should free some memory */
2147 li->blbnb.resize(li->ncc);
2150 li->blmf.resize(li->ncc);
2151 li->blmf1.resize(li->ncc);
2152 li->tmpncc.resize(li->ncc);
2154 gmx::ArrayRef<const int> nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2155 if (!nlocat_dd.empty())
2157 /* Convert nlocat from local topology to LINCS constraint indexing */
2158 for (con = 0; con < ncon_tot; con++)
2160 li->nlocat[li->con_index[con]] = nlocat_dd[con];
2170 fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n",
2171 li->nc_real, li->nc, li->ncc);
2176 lincs_thread_setup(li, md.nr);
2179 set_lincs_matrix(li, md.invmass, md.lambda);
2182 //! Issues a warning when LINCS constraints cannot be satisfied.
2183 static void lincs_warning(gmx_domdec_t *dd, const rvec *x, rvec *xprime, t_pbc *pbc,
2184 int ncons, gmx::ArrayRef<const AtomPair> atoms,
2185 gmx::ArrayRef<const real> bllen, real wangle,
2186 int maxwarn, int *warncount)
2188 real wfac = std::cos(DEG2RAD*wangle);
2191 "bonds that rotated more than %g degrees:\n"
2192 " atom 1 atom 2 angle previous, current, constraint length\n",
2195 for (int b = 0; b < ncons; b++)
2197 const int i = atoms[b].index1;
2198 const int j = atoms[b].index2;
2203 pbc_dx_aiuc(pbc, x[i], x[j], v0);
2204 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2208 rvec_sub(x[i], x[j], v0);
2209 rvec_sub(xprime[i], xprime[j], v1);
2213 real cosine = ::iprod(v0, v1)/(d0*d1);
2217 " %6d %6d %5.1f %8.4f %8.4f %8.4f\n",
2218 ddglatnr(dd, i), ddglatnr(dd, j),
2219 RAD2DEG*std::acos(cosine), d0, d1, bllen[b]);
2220 if (!std::isfinite(d1))
2222 gmx_fatal(FARGS, "Bond length not finite.");
2228 if (*warncount > maxwarn)
2230 too_many_constraint_warnings(econtLINCS, *warncount);
2234 //! Determine how well the constraints have been satisfied.
2235 static void cconerr(const Lincs &lincsd,
2236 const rvec *x, const t_pbc *pbc,
2237 real *ncons_loc, real *ssd, real *max, int *imax)
2239 gmx::ArrayRef<const AtomPair> atoms = lincsd.atoms;
2240 gmx::ArrayRef<const real> bllen = lincsd.bllen;
2241 gmx::ArrayRef<const int> nlocat = lincsd.nlocat;
2247 for (int task = 0; task < lincsd.ntask; task++)
2249 for (int b = lincsd.task[task].b0; b < lincsd.task[task].b1; b++)
2254 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
2258 rvec_sub(x[atoms[b].index1], x[atoms[b].index2], dx);
2260 real r2 = ::norm2(dx);
2261 real len = r2*gmx::invsqrt(r2);
2262 real d = std::abs(len/bllen[b]-1);
2263 if (d > ma && (nlocat.empty() || nlocat[b]))
2275 ssd2 += nlocat[b]*d*d;
2281 *ncons_loc = (nlocat.empty() ? 1 : 0.5)*count;
2282 *ssd = (nlocat.empty() ? 1 : 0.5)*ssd2;
2287 bool constrain_lincs(bool computeRmsd,
2288 const t_inputrec &ir,
2290 Lincs *lincsd, const t_mdatoms &md,
2291 const t_commrec *cr,
2292 const gmx_multisim_t *ms,
2293 const rvec *x, rvec *xprime, rvec *min_proj,
2294 matrix box, t_pbc *pbc,
2295 real lambda, real *dvdlambda,
2296 real invdt, rvec *v,
2297 bool bCalcVir, tensor vir_r_m_dr,
2298 ConstraintVariable econq,
2300 int maxwarn, int *warncount)
2304 /* This boolean should be set by a flag passed to this routine.
2305 * We can also easily check if any constraint length is changed,
2306 * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2308 bool bCalcDHDL = (ir.efep != efepNO && dvdlambda != nullptr);
2310 if (lincsd->nc == 0 && cr->dd == nullptr)
2314 lincsd->rmsdData = {{0}};
2320 if (econq == ConstraintVariable::Positions)
2322 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2323 * also with efep!=fepNO.
2325 if (ir.efep != efepNO)
2327 if (md.nMassPerturbed && lincsd->matlam != md.lambda)
2329 set_lincs_matrix(lincsd, md.invmass, md.lambda);
2332 for (int i = 0; i < lincsd->nc; i++)
2334 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
2338 if (lincsd->ncg_flex)
2340 /* Set the flexible constraint lengths to the old lengths */
2343 for (int i = 0; i < lincsd->nc; i++)
2345 if (lincsd->bllen[i] == 0)
2348 pbc_dx_aiuc(pbc, x[lincsd->atoms[i].index1], x[lincsd->atoms[i].index2], dx);
2349 lincsd->bllen[i] = norm(dx);
2355 for (int i = 0; i < lincsd->nc; i++)
2357 if (lincsd->bllen[i] == 0)
2360 std::sqrt(distance2(x[lincsd->atoms[i].index1],
2361 x[lincsd->atoms[i].index2]));
2373 cconerr(*lincsd, xprime, pbc,
2374 &ncons_loc, &p_ssd, &p_max, &p_imax);
2377 /* This bWarn var can be updated by multiple threads
2378 * at the same time. But as we only need to detect
2379 * if a warning occurred or not, this is not an issue.
2383 /* The OpenMP parallel region of constrain_lincs for coords */
2384 #pragma omp parallel num_threads(lincsd->ntask)
2388 int th = gmx_omp_get_thread_num();
2390 clear_mat(lincsd->task[th].vir_r_m_dr);
2392 do_lincs(x, xprime, box, pbc, lincsd, th,
2395 ir.LincsWarnAngle, &bWarn,
2397 th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2399 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2402 if (debug && lincsd->nc > 0)
2404 fprintf(debug, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2405 fprintf(debug, " Before LINCS %.6f %.6f %6d %6d\n",
2406 std::sqrt(p_ssd/ncons_loc), p_max,
2407 ddglatnr(cr->dd, lincsd->atoms[p_imax].index1),
2408 ddglatnr(cr->dd, lincsd->atoms[p_imax].index2));
2410 if (computeRmsd || debug)
2412 cconerr(*lincsd, xprime, pbc,
2413 &ncons_loc, &p_ssd, &p_max, &p_imax);
2414 lincsd->rmsdData[0] = ncons_loc;
2415 lincsd->rmsdData[1] = p_ssd;
2419 lincsd->rmsdData = {{0}};
2421 if (debug && lincsd->nc > 0)
2424 " After LINCS %.6f %.6f %6d %6d\n\n",
2425 std::sqrt(p_ssd/ncons_loc), p_max,
2426 ddglatnr(cr->dd, lincsd->atoms[p_imax].index1),
2427 ddglatnr(cr->dd, lincsd->atoms[p_imax].index2));
2432 if (maxwarn < INT_MAX)
2434 cconerr(*lincsd, xprime, pbc,
2435 &ncons_loc, &p_ssd, &p_max, &p_imax);
2436 std::string simMesg;
2439 simMesg += gmx::formatString(" in simulation %d", ms->sim);
2442 "\nStep %" PRId64 ", time %g (ps) LINCS WARNING%s\n"
2443 "relative constraint deviation after LINCS:\n"
2444 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2445 step, ir.init_t+step*ir.delta_t,
2447 std::sqrt(p_ssd/ncons_loc), p_max,
2448 ddglatnr(cr->dd, lincsd->atoms[p_imax].index1),
2449 ddglatnr(cr->dd, lincsd->atoms[p_imax].index2));
2451 lincs_warning(cr->dd, x, xprime, pbc,
2452 lincsd->nc, lincsd->atoms, lincsd->bllen,
2453 ir.LincsWarnAngle, maxwarn, warncount);
2455 bOK = (p_max < 0.5);
2458 if (lincsd->ncg_flex)
2460 for (int i = 0; (i < lincsd->nc); i++)
2462 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2464 lincsd->bllen[i] = 0;
2471 /* The OpenMP parallel region of constrain_lincs for derivatives */
2472 #pragma omp parallel num_threads(lincsd->ntask)
2476 int th = gmx_omp_get_thread_num();
2478 do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
2479 md.invmass, econq, bCalcDHDL,
2480 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2482 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2488 /* Reduce the dH/dlambda contributions over the threads */
2493 for (th = 0; th < lincsd->ntask; th++)
2495 dhdlambda += lincsd->task[th].dhdlambda;
2497 if (econq == ConstraintVariable::Positions)
2499 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2500 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2501 dhdlambda /= ir.delta_t*ir.delta_t;
2503 *dvdlambda += dhdlambda;
2506 if (bCalcVir && lincsd->ntask > 1)
2508 for (int i = 1; i < lincsd->ntask; i++)
2510 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2514 /* count assuming nit=1 */
2515 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2516 inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
2517 if (lincsd->ntriangle > 0)
2519 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
2523 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real*2);
2527 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);