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 //! Unit of work within LINCS.
96 //! First constraint for this task.
98 //! b1-1 is the last constraint for this task.
100 //! The number of constraints in triangles.
102 //! The list of triangle constraints.
103 std::vector<int> triangle;
104 //! The bits tell if the matrix element should be used.
105 std::vector<int> tri_bits;
106 //! Constraint index for updating atom data.
107 std::vector<int> ind;
108 //! Constraint index for updating atom data.
109 std::vector<int> ind_r;
110 //! Temporary variable for virial calculation.
111 tensor vir_r_m_dr = {{0}};
112 //! Temporary variable for lambda derivative.
116 /*! \brief Data for LINCS algorithm.
121 //! The global number of constraints.
123 //! The global number of flexible constraints.
125 //! The global number of constraints in triangles.
126 int ncg_triangle = 0;
127 //! The number of iterations.
129 //! The order of the matrix expansion.
131 //! The maximum number of constraints connected to a single atom.
134 //! The number of real constraints.
136 //! The number of constraints including padding for SIMD.
138 //! The number of constraint connections.
140 //! The FE lambda value used for filling blc and blmf.
142 //! mapping from topology to LINCS constraints.
143 std::vector<int> con_index;
144 //! The reference distance in topology A.
145 std::vector < real, AlignedAllocator < real>> bllen0;
146 //! The reference distance in top B - the r.d. in top A.
147 std::vector < real, AlignedAllocator < real>> ddist;
148 //! The atom pairs involved in the constraints.
149 std::vector<int> bla;
150 //! 1/sqrt(invmass1 invmass2).
151 std::vector < real, AlignedAllocator < real>> blc;
152 //! As blc, but with all masses 1.
153 std::vector < real, AlignedAllocator < real>> blc1;
154 //! Index into blbnb and blmf.
155 std::vector<int> blnr;
156 //! List of constraint connections.
157 std::vector<int> blbnb;
158 //! The local number of constraints in triangles.
160 //! The number of constraint connections in triangles.
161 int ncc_triangle = 0;
162 //! Communicate before each LINCS interation.
163 bool bCommIter = false;
164 //! Matrix of mass factors for constraint connections.
165 std::vector<real> blmf;
166 //! As blmf, but with all masses 1.
167 std::vector<real> blmf1;
168 //! The reference bond length.
169 std::vector < real, AlignedAllocator < real>> bllen;
170 //! The local atom count per constraint, can be NULL.
171 std::vector<int> nlocat;
173 /*! \brief The number of tasks used for LINCS work.
175 * \todo This is mostly used to loop over \c task, which would
176 * be nicer to do with range-based for loops, but the thread
177 * index is used for constructing bit masks and organizing the
178 * virial output buffer, so other things need to change,
181 /*! \brief LINCS thread division */
182 std::vector<Task> task;
183 //! Atom flags for thread parallelization.
184 std::vector<gmx_bitmask_t> atf;
185 //! Are the LINCS tasks interdependent?
186 bool bTaskDep = false;
187 //! Are there triangle constraints that cross task borders?
188 bool bTaskDepTri = false;
189 //! Arrays for temporary storage in the LINCS algorithm.
191 PaddedVector<gmx::RVec> tmpv;
192 std::vector<real> tmpncc;
193 std::vector < real, AlignedAllocator < real>> tmp1;
194 std::vector < real, AlignedAllocator < real>> tmp2;
195 std::vector < real, AlignedAllocator < real>> tmp3;
196 std::vector < real, AlignedAllocator < real>> tmp4;
198 //! The Lagrange multipliers times -1.
199 std::vector < real, AlignedAllocator < real>> mlambda;
200 //! Storage for the constraint RMS relative deviation output.
201 std::array<real, 2> rmsdData = {{0}};
204 /*! \brief Define simd_width for memory allocation used for SIMD code */
205 #if GMX_SIMD_HAVE_REAL
206 static const int simd_width = GMX_SIMD_REAL_WIDTH;
208 static const int simd_width = 1;
211 ArrayRef<real> lincs_rmsdData(Lincs *lincsd)
213 return lincsd->rmsdData;
216 real lincs_rmsd(const Lincs *lincsd)
218 if (lincsd->rmsdData[0] > 0)
220 return std::sqrt(lincsd->rmsdData[1]/lincsd->rmsdData[0]);
228 /*! \brief Do a set of nrec LINCS matrix multiplications.
230 * This function will return with up to date thread-local
231 * constraint data, without an OpenMP barrier.
233 static void lincs_matrix_expand(const Lincs *lincsd,
236 real *rhs1, real *rhs2, real *sol)
238 gmx::ArrayRef<const int> blnr = lincsd->blnr;
239 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
241 const int b0 = li_task->b0;
242 const int b1 = li_task->b1;
243 const int nrec = lincsd->nOrder;
245 for (int rec = 0; rec < nrec; rec++)
247 if (lincsd->bTaskDep)
251 for (int b = b0; b < b1; b++)
257 for (n = blnr[b]; n < blnr[b+1]; n++)
259 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
262 sol[b] = sol[b] + mvb;
270 } /* nrec*(ncons+2*nrtot) flops */
272 if (lincsd->ntriangle > 0)
274 /* Perform an extra nrec recursions for only the constraints
275 * involved in rigid triangles.
276 * In this way their accuracy should come close to those of the other
277 * constraints, since traingles of constraints can produce eigenvalues
278 * around 0.7, while the effective eigenvalue for bond constraints
279 * is around 0.4 (and 0.7*0.7=0.5).
282 if (lincsd->bTaskDep)
284 /* We need a barrier here, since other threads might still be
285 * reading the contents of rhs1 and/o rhs2.
286 * We could avoid this barrier by introducing two extra rhs
287 * arrays for the triangle constraints only.
292 /* Constraints involved in a triangle are ensured to be in the same
293 * LINCS task. This means no barriers are required during the extra
294 * iterations for the triangle constraints.
296 gmx::ArrayRef<const int> triangle = li_task->triangle;
297 gmx::ArrayRef<const int> tri_bits = li_task->tri_bits;
299 for (int rec = 0; rec < nrec; rec++)
301 for (int tb = 0; tb < li_task->ntriangle; tb++)
303 int b, bits, nr0, nr1, n;
311 for (n = nr0; n < nr1; n++)
313 if (bits & (1 << (n - nr0)))
315 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
319 sol[b] = sol[b] + mvb;
327 } /* nrec*(ntriangle + ncc_triangle*2) flops */
329 if (lincsd->bTaskDepTri)
331 /* The constraints triangles are decoupled from each other,
332 * but constraints in one triangle cross thread task borders.
333 * We could probably avoid this with more advanced setup code.
340 //! Update atomic coordinates when an index is not required.
341 static void lincs_update_atoms_noind(int ncons,
342 gmx::ArrayRef<const int> bla,
344 const real *fac, rvec *r,
349 real mvb, im1, im2, tmp0, tmp1, tmp2;
351 if (invmass != nullptr)
353 for (b = 0; b < ncons; b++)
369 } /* 16 ncons flops */
373 for (b = 0; b < ncons; b++)
391 //! Update atomic coordinates when an index is required.
392 static void lincs_update_atoms_ind(gmx::ArrayRef<const int> ind,
393 gmx::ArrayRef<const int> bla,
395 const real *fac, rvec *r,
400 real mvb, im1, im2, tmp0, tmp1, tmp2;
402 if (invmass != nullptr)
420 } /* 16 ncons flops */
438 } /* 16 ncons flops */
442 //! Update coordinates for atoms.
443 static void lincs_update_atoms(Lincs *li, int th,
445 const real *fac, rvec *r,
451 /* Single thread, we simply update for all constraints */
452 lincs_update_atoms_noind(li->nc_real,
453 li->bla, prefac, fac, r, invmass, x);
457 /* Update the atom vector components for our thread local
458 * constraints that only access our local atom range.
459 * This can be done without a barrier.
461 lincs_update_atoms_ind(li->task[th].ind,
462 li->bla, prefac, fac, r, invmass, x);
464 if (!li->task[li->ntask].ind.empty())
466 /* Update the constraints that operate on atoms
467 * in multiple thread atom blocks on the master thread.
472 lincs_update_atoms_ind(li->task[li->ntask].ind,
473 li->bla, prefac, fac, r, invmass, x);
479 #if GMX_SIMD_HAVE_REAL
480 //! Helper function so that we can run TSAN with SIMD support (where implemented).
482 static inline void gmx_simdcall
483 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
503 calc_dr_x_f_simd(int b0,
506 const rvec * gmx_restrict x,
507 const rvec * gmx_restrict f,
508 const real * gmx_restrict blc,
509 const real * pbc_simd,
510 rvec * gmx_restrict r,
511 real * gmx_restrict rhs,
512 real * gmx_restrict sol)
514 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
516 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
518 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
523 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
525 SimdReal x0_S, y0_S, z0_S;
526 SimdReal x1_S, y1_S, z1_S;
527 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
528 SimdReal fx_S, fy_S, fz_S, ip_S, rhs_S;
529 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
530 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
532 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
534 offset0[i] = bla[bs*2 + i*2];
535 offset1[i] = bla[bs*2 + i*2 + 1];
538 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
539 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
544 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
546 n2_S = norm2(rx_S, ry_S, rz_S);
547 il_S = invsqrt(n2_S);
553 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
555 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(f), offset0, &x0_S, &y0_S, &z0_S);
556 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(f), offset1, &x1_S, &y1_S, &z1_S);
561 ip_S = iprod(rx_S, ry_S, rz_S, fx_S, fy_S, fz_S);
563 rhs_S = load<SimdReal>(blc + bs) * ip_S;
565 store(rhs + bs, rhs_S);
566 store(sol + bs, rhs_S);
569 #endif // GMX_SIMD_HAVE_REAL
571 /*! \brief LINCS projection, works on derivatives of the coordinates. */
572 static void do_lincsp(const rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
573 Lincs *lincsd, int th,
575 ConstraintVariable econq, bool bCalcDHDL,
576 bool bCalcVir, tensor rmdf)
579 int *bla, *blnr, *blbnb;
581 real *blc, *blmf, *blcc, *rhs1, *rhs2, *sol;
583 b0 = lincsd->task[th].b0;
584 b1 = lincsd->task[th].b1;
586 bla = lincsd->bla.data();
587 r = as_rvec_array(lincsd->tmpv.data());
588 blnr = lincsd->blnr.data();
589 blbnb = lincsd->blbnb.data();
590 if (econq != ConstraintVariable::Force)
592 /* Use mass-weighted parameters */
593 blc = lincsd->blc.data();
594 blmf = lincsd->blmf.data();
598 /* Use non mass-weighted parameters */
599 blc = lincsd->blc1.data();
600 blmf = lincsd->blmf1.data();
602 blcc = lincsd->tmpncc.data();
603 rhs1 = lincsd->tmp1.data();
604 rhs2 = lincsd->tmp2.data();
605 sol = lincsd->tmp3.data();
607 #if GMX_SIMD_HAVE_REAL
608 /* This SIMD code does the same as the plain-C code after the #else.
609 * The only difference is that we always call pbc code, as with SIMD
610 * the overhead of pbc computation (when not needed) is small.
612 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
614 /* Convert the pbc struct for SIMD */
615 set_pbc_simd(pbc, pbc_simd);
617 /* Compute normalized x i-j vectors, store in r.
618 * Compute the inner product of r and xp i-j and store in rhs1.
620 calc_dr_x_f_simd(b0, b1, bla, x, f, blc,
624 #else // GMX_SIMD_HAVE_REAL
626 /* Compute normalized i-j vectors */
629 for (b = b0; b < b1; b++)
633 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
639 for (b = b0; b < b1; b++)
643 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
645 } /* 16 ncons flops */
648 for (b = b0; b < b1; b++)
655 mvb = blc[b]*(r[b][0]*(f[i][0] - f[j][0]) +
656 r[b][1]*(f[i][1] - f[j][1]) +
657 r[b][2]*(f[i][2] - f[j][2]));
663 #endif // GMX_SIMD_HAVE_REAL
665 if (lincsd->bTaskDep)
667 /* We need a barrier, since the matrix construction below
668 * can access entries in r of other threads.
673 /* Construct the (sparse) LINCS matrix */
674 for (b = b0; b < b1; b++)
678 for (n = blnr[b]; n < blnr[b+1]; n++)
680 blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
683 /* Together: 23*ncons + 6*nrtot flops */
685 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
686 /* nrec*(ncons+2*nrtot) flops */
688 if (econq == ConstraintVariable::Deriv_FlexCon)
690 /* We only want to constraint the flexible constraints,
691 * so we mask out the normal ones by setting sol to 0.
693 for (b = b0; b < b1; b++)
695 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
702 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
703 for (b = b0; b < b1; b++)
708 /* When constraining forces, we should not use mass weighting,
709 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
711 lincs_update_atoms(lincsd, th, 1.0, sol, r,
712 (econq != ConstraintVariable::Force) ? invmass : nullptr, fp);
719 for (b = b0; b < b1; b++)
721 dhdlambda -= sol[b]*lincsd->ddist[b];
724 lincsd->task[th].dhdlambda = dhdlambda;
729 /* Constraint virial,
730 * determines sum r_bond x delta f,
731 * where delta f is the constraint correction
732 * of the quantity that is being constrained.
734 for (b = b0; b < b1; b++)
739 mvb = lincsd->bllen[b]*sol[b];
740 for (i = 0; i < DIM; i++)
743 for (j = 0; j < DIM; j++)
745 rmdf[i][j] += tmp1*r[b][j];
748 } /* 23 ncons flops */
752 #if GMX_SIMD_HAVE_REAL
754 /*! \brief Calculate the constraint distance vectors r to project on from x.
756 * Determine the right-hand side of the matrix equation using coordinates xp. */
757 static void gmx_simdcall
758 calc_dr_x_xp_simd(int b0,
761 const rvec * gmx_restrict x,
762 const rvec * gmx_restrict xp,
763 const real * gmx_restrict bllen,
764 const real * gmx_restrict blc,
765 const real * pbc_simd,
766 rvec * gmx_restrict r,
767 real * gmx_restrict rhs,
768 real * gmx_restrict sol)
770 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
771 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
773 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
778 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
780 SimdReal x0_S, y0_S, z0_S;
781 SimdReal x1_S, y1_S, z1_S;
782 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
783 SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
784 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
785 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
787 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
789 offset0[i] = bla[bs*2 + i*2];
790 offset1[i] = bla[bs*2 + i*2 + 1];
793 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
794 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
799 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
801 n2_S = norm2(rx_S, ry_S, rz_S);
802 il_S = invsqrt(n2_S);
808 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
810 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(xp), offset0, &x0_S, &y0_S, &z0_S);
811 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(xp), offset1, &x1_S, &y1_S, &z1_S);
817 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
819 ip_S = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
821 rhs_S = load<SimdReal>(blc + bs) * (ip_S - load<SimdReal>(bllen + bs));
823 store(rhs + bs, rhs_S);
824 store(sol + bs, rhs_S);
827 #endif // GMX_SIMD_HAVE_REAL
829 /*! \brief Determine the distances and right-hand side for the next iteration. */
830 gmx_unused static void calc_dist_iter(
834 const rvec * gmx_restrict xp,
835 const real * gmx_restrict bllen,
836 const real * gmx_restrict blc,
839 real * gmx_restrict rhs,
840 real * gmx_restrict sol,
845 for (b = b0; b < b1; b++)
847 real len, len2, dlen2, mvb;
853 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
857 rvec_sub(xp[bla[2*b]], xp[bla[2*b+1]], dx);
860 dlen2 = 2*len2 - ::norm2(dx);
861 if (dlen2 < wfac*len2)
863 /* not race free - see detailed comment in caller */
868 mvb = blc[b]*(len - dlen2*gmx::invsqrt(dlen2));
876 } /* 20*ncons flops */
879 #if GMX_SIMD_HAVE_REAL
880 /*! \brief As calc_dist_iter(), but using SIMD intrinsics. */
881 static void gmx_simdcall
882 calc_dist_iter_simd(int b0,
885 const rvec * gmx_restrict x,
886 const real * gmx_restrict bllen,
887 const real * gmx_restrict blc,
888 const real * pbc_simd,
890 real * gmx_restrict rhs,
891 real * gmx_restrict sol,
894 SimdReal min_S(GMX_REAL_MIN);
896 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 (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] = bla[bs*2 + i*2];
918 offset1[i] = bla[bs*2 + i*2 + 1];
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 int b0, b1, b, i, j, n, iter;
973 int *bla, *blnr, *blbnb;
975 real *blc, *blmf, *bllen, *blcc, *rhs1, *rhs2, *sol, *blc_sol, *mlambda;
978 b0 = lincsd->task[th].b0;
979 b1 = lincsd->task[th].b1;
981 bla = lincsd->bla.data();
982 r = as_rvec_array(lincsd->tmpv.data());
983 blnr = lincsd->blnr.data();
984 blbnb = lincsd->blbnb.data();
985 blc = lincsd->blc.data();
986 blmf = lincsd->blmf.data();
987 bllen = lincsd->bllen.data();
988 blcc = lincsd->tmpncc.data();
989 rhs1 = lincsd->tmp1.data();
990 rhs2 = lincsd->tmp2.data();
991 sol = lincsd->tmp3.data();
992 blc_sol = lincsd->tmp4.data();
993 mlambda = lincsd->mlambda.data();
994 nlocat = lincsd->nlocat.data();
996 #if GMX_SIMD_HAVE_REAL
998 /* This SIMD code does the same as the plain-C code after the #else.
999 * The only difference is that we always call pbc code, as with SIMD
1000 * the overhead of pbc computation (when not needed) is small.
1002 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
1004 /* Convert the pbc struct for SIMD */
1005 set_pbc_simd(pbc, pbc_simd);
1007 /* Compute normalized x i-j vectors, store in r.
1008 * Compute the inner product of r and xp i-j and store in rhs1.
1010 calc_dr_x_xp_simd(b0, b1, bla, x, xp, bllen, blc,
1014 #else // GMX_SIMD_HAVE_REAL
1018 /* Compute normalized i-j vectors */
1019 for (b = b0; b < b1; b++)
1024 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
1027 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
1028 mvb = blc[b]*(::iprod(r[b], dx) - bllen[b]);
1035 /* Compute normalized i-j vectors */
1036 for (b = b0; b < b1; b++)
1038 real tmp0, tmp1, tmp2, rlen, mvb;
1042 tmp0 = x[i][0] - x[j][0];
1043 tmp1 = x[i][1] - x[j][1];
1044 tmp2 = x[i][2] - x[j][2];
1045 rlen = gmx::invsqrt(tmp0*tmp0 + tmp1*tmp1 + tmp2*tmp2);
1046 r[b][0] = rlen*tmp0;
1047 r[b][1] = rlen*tmp1;
1048 r[b][2] = rlen*tmp2;
1049 /* 16 ncons flops */
1053 mvb = blc[b]*(r[b][0]*(xp[i][0] - xp[j][0]) +
1054 r[b][1]*(xp[i][1] - xp[j][1]) +
1055 r[b][2]*(xp[i][2] - xp[j][2]) - bllen[b]);
1060 /* Together: 26*ncons + 6*nrtot flops */
1063 #endif // GMX_SIMD_HAVE_REAL
1065 if (lincsd->bTaskDep)
1067 /* We need a barrier, since the matrix construction below
1068 * can access entries in r of other threads.
1073 /* Construct the (sparse) LINCS matrix */
1074 for (b = b0; b < b1; b++)
1076 for (n = blnr[b]; n < blnr[b+1]; n++)
1078 blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
1081 /* Together: 26*ncons + 6*nrtot flops */
1083 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1084 /* nrec*(ncons+2*nrtot) flops */
1086 #if GMX_SIMD_HAVE_REAL
1087 for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1089 SimdReal t1 = load<SimdReal>(blc + b);
1090 SimdReal t2 = load<SimdReal>(sol + b);
1091 store(mlambda + b, t1 * t2);
1094 for (b = b0; b < b1; b++)
1096 mlambda[b] = blc[b]*sol[b];
1098 #endif // GMX_SIMD_HAVE_REAL
1100 /* Update the coordinates */
1101 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1104 ******** Correction for centripetal effects ********
1109 wfac = std::cos(DEG2RAD*wangle);
1112 for (iter = 0; iter < lincsd->nIter; iter++)
1114 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
1119 /* Communicate the corrected non-local coordinates */
1120 if (DOMAINDECOMP(cr))
1122 dd_move_x_constraints(cr->dd, box, xp, nullptr, FALSE);
1127 else if (lincsd->bTaskDep)
1132 #if GMX_SIMD_HAVE_REAL
1133 calc_dist_iter_simd(b0, b1, bla, xp, bllen, blc, pbc_simd, wfac,
1136 calc_dist_iter(b0, b1, bla, xp, bllen, blc, pbc, wfac,
1138 /* 20*ncons flops */
1139 #endif // GMX_SIMD_HAVE_REAL
1141 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1142 /* nrec*(ncons+2*nrtot) flops */
1144 #if GMX_SIMD_HAVE_REAL
1145 for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1147 SimdReal t1 = load<SimdReal>(blc + b);
1148 SimdReal t2 = load<SimdReal>(sol + b);
1149 SimdReal mvb = t1 * t2;
1150 store(blc_sol + b, mvb);
1151 store(mlambda + b, load<SimdReal>(mlambda + b) + mvb);
1154 for (b = b0; b < b1; b++)
1158 mvb = blc[b]*sol[b];
1162 #endif // GMX_SIMD_HAVE_REAL
1164 /* Update the coordinates */
1165 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1167 /* nit*ncons*(37+9*nrec) flops */
1171 /* Update the velocities */
1172 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1173 /* 16 ncons flops */
1176 if (nlocat != nullptr && (bCalcDHDL || bCalcVir))
1178 if (lincsd->bTaskDep)
1180 /* In lincs_update_atoms threads might cross-read mlambda */
1184 /* Only account for local atoms */
1185 for (b = b0; b < b1; b++)
1187 mlambda[b] *= 0.5*nlocat[b];
1196 for (b = b0; b < b1; b++)
1198 /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1199 * later after the contributions are reduced over the threads.
1201 dhdl -= lincsd->mlambda[b]*lincsd->ddist[b];
1203 lincsd->task[th].dhdlambda = dhdl;
1208 /* Constraint virial */
1209 for (b = b0; b < b1; b++)
1213 tmp0 = -bllen[b]*mlambda[b];
1214 for (i = 0; i < DIM; i++)
1216 tmp1 = tmp0*r[b][i];
1217 for (j = 0; j < DIM; j++)
1219 vir_r_m_dr[i][j] -= tmp1*r[b][j];
1222 } /* 22 ncons flops */
1226 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1227 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1229 * (26+nrec)*ncons + (6+2*nrec)*nrtot
1230 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1232 * (63+nrec)*ncons + (6+4*nrec)*nrtot
1236 /*! \brief Sets the elements in the LINCS matrix for task task. */
1237 static void set_lincs_matrix_task(Lincs *li,
1239 const real *invmass,
1241 int *nCrossTaskTriangles)
1245 /* Construct the coupling coefficient matrix blmf */
1246 li_task->ntriangle = 0;
1248 *nCrossTaskTriangles = 0;
1249 for (i = li_task->b0; i < li_task->b1; i++)
1254 a2 = li->bla[2*i+1];
1255 for (n = li->blnr[i]; (n < li->blnr[i+1]); n++)
1257 int k, sign, center, end;
1261 /* If we are using multiple, independent tasks for LINCS,
1262 * the calls to check_assign_connected should have
1263 * put all connected constraints in our task.
1265 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1267 if (a1 == li->bla[2*k] || a2 == li->bla[2*k+1])
1275 if (a1 == li->bla[2*k] || a1 == li->bla[2*k+1])
1285 li->blmf[n] = sign*invmass[center]*li->blc[i]*li->blc[k];
1286 li->blmf1[n] = sign*0.5;
1287 if (li->ncg_triangle > 0)
1291 /* Look for constraint triangles */
1292 for (nk = li->blnr[k]; (nk < li->blnr[k+1]); nk++)
1295 if (kk != i && kk != k &&
1296 (li->bla[2*kk] == end || li->bla[2*kk+1] == end))
1298 /* Check if the constraints in this triangle actually
1299 * belong to a different task. We still assign them
1300 * here, since it's convenient for the triangle
1301 * iterations, but we then need an extra barrier.
1303 if (k < li_task->b0 || k >= li_task->b1 ||
1304 kk < li_task->b0 || kk >= li_task->b1)
1306 (*nCrossTaskTriangles)++;
1309 if (li_task->ntriangle == 0 ||
1310 li_task->triangle[li_task->ntriangle - 1] < i)
1312 /* Add this constraint to the triangle list */
1313 li_task->triangle[li_task->ntriangle] = i;
1314 li_task->tri_bits[li_task->ntriangle] = 0;
1315 li_task->ntriangle++;
1316 if (li->blnr[i+1] - li->blnr[i] > static_cast<int>(sizeof(li_task->tri_bits[0])*8 - 1))
1318 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %zu allowed for constraints participating in triangles",
1319 li->blnr[i+1] - li->blnr[i],
1320 sizeof(li_task->tri_bits[0])*8-1);
1323 li_task->tri_bits[li_task->ntriangle-1] |= (1 << (n - li->blnr[i]));
1332 /*! \brief Sets the elements in the LINCS matrix. */
1333 static void set_lincs_matrix(Lincs *li, real *invmass, real lambda)
1336 const real invsqrt2 = 0.7071067811865475244;
1338 for (i = 0; (i < li->nc); i++)
1343 a2 = li->bla[2*i+1];
1344 li->blc[i] = gmx::invsqrt(invmass[a1] + invmass[a2]);
1345 li->blc1[i] = invsqrt2;
1348 /* Construct the coupling coefficient matrix blmf */
1349 int th, ntriangle = 0, ncc_triangle = 0, nCrossTaskTriangles = 0;
1350 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle, nCrossTaskTriangles) num_threads(li->ntask) schedule(static)
1351 for (th = 0; th < li->ntask; th++)
1355 set_lincs_matrix_task(li, &li->task[th], invmass,
1356 &ncc_triangle, &nCrossTaskTriangles);
1357 ntriangle += li->task[th].ntriangle;
1359 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1361 li->ntriangle = ntriangle;
1362 li->ncc_triangle = ncc_triangle;
1363 li->bTaskDepTri = (nCrossTaskTriangles > 0);
1367 fprintf(debug, "The %d constraints participate in %d triangles\n",
1368 li->nc, li->ntriangle);
1369 fprintf(debug, "There are %d constraint couplings, of which %d in triangles\n",
1370 li->ncc, li->ncc_triangle);
1371 if (li->ntriangle > 0 && li->ntask > 1)
1373 fprintf(debug, "%d constraint triangles contain constraints assigned to different tasks\n",
1374 nCrossTaskTriangles);
1379 * so we know with which lambda value the masses have been set.
1381 li->matlam = lambda;
1384 //! Finds all triangles of atoms that share constraints to a central atom.
1385 static int count_triangle_constraints(const InteractionLists &ilist,
1386 const t_blocka &at2con)
1388 int ncon1, ncon_tot;
1389 int c0, n1, c1, ac1, n2, c2;
1392 ncon1 = ilist[F_CONSTR].size()/3;
1393 ncon_tot = ncon1 + ilist[F_CONSTRNC].size()/3;
1395 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1396 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1399 for (c0 = 0; c0 < ncon_tot; c0++)
1401 bool bTriangle = FALSE;
1402 const int *iap = constr_iatomptr(ia1, ia2, c0);
1403 const int a00 = iap[1];
1404 const int a01 = iap[2];
1405 for (n1 = at2con.index[a01]; n1 < at2con.index[a01+1]; n1++)
1410 const int *iap = constr_iatomptr(ia1, ia2, c1);
1411 const int a10 = iap[1];
1412 const int a11 = iap[2];
1421 for (n2 = at2con.index[ac1]; n2 < at2con.index[ac1+1]; n2++)
1424 if (c2 != c0 && c2 != c1)
1426 const int *iap = constr_iatomptr(ia1, ia2, c2);
1427 const int a20 = iap[1];
1428 const int a21 = iap[2];
1429 if (a20 == a00 || a21 == a00)
1443 return ncon_triangle;
1446 //! Finds sequences of sequential constraints.
1447 static bool more_than_two_sequential_constraints(const InteractionLists &ilist,
1448 const t_blocka &at2con)
1450 int ncon1, ncon_tot, c;
1451 bool bMoreThanTwoSequentialConstraints;
1453 ncon1 = ilist[F_CONSTR].size()/3;
1454 ncon_tot = ncon1 + ilist[F_CONSTRNC].size()/3;
1456 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1457 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1459 bMoreThanTwoSequentialConstraints = FALSE;
1460 for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
1462 const int *iap = constr_iatomptr(ia1, ia2, c);
1463 const int a1 = iap[1];
1464 const int a2 = iap[2];
1465 /* Check if this constraint has constraints connected at both atoms */
1466 if (at2con.index[a1+1] - at2con.index[a1] > 1 &&
1467 at2con.index[a2+1] - at2con.index[a2] > 1)
1469 bMoreThanTwoSequentialConstraints = TRUE;
1473 return bMoreThanTwoSequentialConstraints;
1476 Lincs *init_lincs(FILE *fplog, const gmx_mtop_t &mtop,
1477 int nflexcon_global, ArrayRef<const t_blocka> at2con,
1478 bool bPLINCS, int nIter, int nProjOrder)
1480 // TODO this should become a unique_ptr
1482 bool bMoreThanTwoSeq;
1486 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
1487 bPLINCS ? " Parallel" : "");
1493 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1494 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1495 li->ncg_flex = nflexcon_global;
1498 li->nOrder = nProjOrder;
1500 li->max_connect = 0;
1501 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1503 for (int a = 0; a < mtop.moltype[mt].atoms.nr; a++)
1505 li->max_connect = std::max(li->max_connect,
1506 at2con[mt].index[a + 1] - at2con[mt].index[a]);
1510 li->ncg_triangle = 0;
1511 bMoreThanTwoSeq = FALSE;
1512 for (const gmx_molblock_t &molb : mtop.molblock)
1514 const gmx_moltype_t &molt = mtop.moltype[molb.type];
1518 count_triangle_constraints(molt.ilist, at2con[molb.type]);
1520 if (!bMoreThanTwoSeq &&
1521 more_than_two_sequential_constraints(molt.ilist, at2con[molb.type]))
1523 bMoreThanTwoSeq = TRUE;
1527 /* Check if we need to communicate not only before LINCS,
1528 * but also before each iteration.
1529 * The check for only two sequential constraints is only
1530 * useful for the common case of H-bond only constraints.
1531 * With more effort we could also make it useful for small
1532 * molecules with nr. sequential constraints <= nOrder-1.
1534 li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1536 if (debug && bPLINCS)
1538 fprintf(debug, "PLINCS communication before each iteration: %d\n",
1539 static_cast<int>(li->bCommIter));
1542 /* LINCS can run on any number of threads.
1543 * Currently the number is fixed for the whole simulation,
1544 * but it could be set in set_lincs().
1545 * The current constraint to task assignment code can create independent
1546 * tasks only when not more than two constraints are connected sequentially.
1548 li->ntask = gmx_omp_nthreads_get(emntLINCS);
1549 li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1552 fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n",
1553 li->ntask, li->bTaskDep ? "" : "in");
1561 /* Allocate an extra elements for "task-overlap" constraints */
1562 li->task.resize(li->ntask + 1);
1565 if (bPLINCS || li->ncg_triangle > 0)
1567 please_cite(fplog, "Hess2008a");
1571 please_cite(fplog, "Hess97a");
1576 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1579 fprintf(fplog, "There are constraints between atoms in different decomposition domains,\n"
1580 "will communicate selected coordinates each lincs iteration\n");
1582 if (li->ncg_triangle > 0)
1585 "%d constraints are involved in constraint triangles,\n"
1586 "will apply an additional matrix expansion of order %d for couplings\n"
1587 "between constraints inside triangles\n",
1588 li->ncg_triangle, li->nOrder);
1595 void done_lincs(Lincs *li)
1600 /*! \brief Sets up the work division over the threads. */
1601 static void lincs_thread_setup(Lincs *li, int natoms)
1603 li->atf.resize(natoms);
1605 gmx::ArrayRef<gmx_bitmask_t> atf = li->atf;
1607 /* Clear the atom flags */
1608 for (gmx_bitmask_t &mask : atf)
1610 bitmask_clear(&mask);
1613 if (li->ntask > BITMASK_SIZE)
1615 gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1618 for (int th = 0; th < li->ntask; th++)
1620 const Task *li_task = &li->task[th];
1622 /* For each atom set a flag for constraints from each */
1623 for (int b = li_task->b0; b < li_task->b1; b++)
1625 bitmask_set_bit(&atf[li->bla[b*2 ]], th);
1626 bitmask_set_bit(&atf[li->bla[b*2 + 1]], th);
1630 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1631 for (int th = 0; th < li->ntask; th++)
1639 li_task = &li->task[th];
1641 bitmask_init_low_bits(&mask, th);
1643 li_task->ind.clear();
1644 li_task->ind_r.clear();
1645 for (b = li_task->b0; b < li_task->b1; b++)
1647 /* We let the constraint with the lowest thread index
1648 * operate on atoms with constraints from multiple threads.
1650 if (bitmask_is_disjoint(atf[li->bla[b*2]], mask) &&
1651 bitmask_is_disjoint(atf[li->bla[b*2+1]], mask))
1653 /* Add the constraint to the local atom update index */
1654 li_task->ind.push_back(b);
1658 /* Add the constraint to the rest block */
1659 li_task->ind_r.push_back(b);
1663 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1666 /* We need to copy all constraints which have not be assigned
1667 * to a thread to a separate list which will be handled by one thread.
1669 Task *li_m = &li->task[li->ntask];
1672 for (int th = 0; th < li->ntask; th++)
1674 const Task &li_task = li->task[th];
1676 for (int ind_r : li_task.ind_r)
1678 li_m->ind.push_back(ind_r);
1683 fprintf(debug, "LINCS thread %d: %zu constraints\n",
1684 th, li_task.ind.size());
1690 fprintf(debug, "LINCS thread r: %zu constraints\n",
1695 //! Assign a constraint.
1696 static void assign_constraint(Lincs *li,
1697 int constraint_index,
1699 real lenA, real lenB,
1700 const t_blocka *at2con)
1706 /* Make an mapping of local topology constraint index to LINCS index */
1707 li->con_index[constraint_index] = con;
1709 li->bllen0[con] = lenA;
1710 li->ddist[con] = lenB - lenA;
1711 /* Set the length to the topology A length */
1712 li->bllen[con] = lenA;
1713 li->bla[2*con] = a1;
1714 li->bla[2*con+1] = a2;
1716 /* Make space in the constraint connection matrix for constraints
1717 * connected to both end of the current constraint.
1720 at2con->index[a1 + 1] - at2con->index[a1] - 1 +
1721 at2con->index[a2 + 1] - at2con->index[a2] - 1;
1723 li->blnr[con + 1] = li->ncc;
1725 /* Increase the constraint count */
1729 /*! \brief Check if constraint with topology index constraint_index is connected
1730 * to other constraints, and if so add those connected constraints to our task. */
1731 static void check_assign_connected(Lincs *li,
1732 const t_iatom *iatom,
1736 const t_blocka *at2con)
1738 /* Currently this function only supports constraint groups
1739 * in which all constraints share at least one atom
1740 * (e.g. H-bond constraints).
1741 * Check both ends of the current constraint for
1742 * connected constraints. We need to assign those
1747 for (end = 0; end < 2; end++)
1751 a = (end == 0 ? a1 : a2);
1753 for (k = at2con->index[a]; k < at2con->index[a + 1]; k++)
1758 /* Check if constraint cc has not yet been assigned */
1759 if (li->con_index[cc] == -1)
1765 lenA = idef.iparams[type].constr.dA;
1766 lenB = idef.iparams[type].constr.dB;
1768 if (bDynamics || lenA != 0 || lenB != 0)
1770 assign_constraint(li, cc, iatom[3*cc + 1], iatom[3*cc + 2], lenA, lenB, at2con);
1777 /*! \brief Check if constraint with topology index constraint_index is involved
1778 * in a constraint triangle, and if so add the other two constraints
1779 * in the triangle to our task. */
1780 static void check_assign_triangle(Lincs *li,
1781 const t_iatom *iatom,
1784 int constraint_index,
1786 const t_blocka *at2con)
1788 int nca, cc[32], ca[32], k;
1789 int c_triangle[2] = { -1, -1 };
1792 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1797 if (c != constraint_index)
1801 aa1 = iatom[c*3 + 1];
1802 aa2 = iatom[c*3 + 2];
1818 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1823 if (c != constraint_index)
1827 aa1 = iatom[c*3 + 1];
1828 aa2 = iatom[c*3 + 2];
1831 for (i = 0; i < nca; i++)
1835 c_triangle[0] = cc[i];
1842 for (i = 0; i < nca; i++)
1846 c_triangle[0] = cc[i];
1854 if (c_triangle[0] >= 0)
1858 for (end = 0; end < 2; end++)
1860 /* Check if constraint c_triangle[end] has not yet been assigned */
1861 if (li->con_index[c_triangle[end]] == -1)
1866 i = c_triangle[end]*3;
1868 lenA = idef.iparams[type].constr.dA;
1869 lenB = idef.iparams[type].constr.dB;
1871 if (bDynamics || lenA != 0 || lenB != 0)
1873 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1880 //! Sets matrix indices.
1881 static void set_matrix_indices(Lincs *li,
1882 const Task *li_task,
1883 const t_blocka *at2con,
1888 for (b = li_task->b0; b < li_task->b1; b++)
1893 a2 = li->bla[b*2 + 1];
1896 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1900 concon = li->con_index[at2con->a[k]];
1903 li->blbnb[i++] = concon;
1906 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1910 concon = li->con_index[at2con->a[k]];
1913 li->blbnb[i++] = concon;
1919 /* Order the blbnb matrix to optimize memory access */
1920 std::sort(&(li->blbnb[li->blnr[b]]), &(li->blbnb[li->blnr[b+1]]));
1925 void set_lincs(const t_idef &idef,
1926 const t_mdatoms &md,
1928 const t_commrec *cr,
1938 /* Zero the thread index ranges.
1939 * Otherwise without local constraints we could return with old ranges.
1941 for (int i = 0; i < li->ntask; i++)
1945 li->task[i].ind.clear();
1949 li->task[li->ntask].ind.clear();
1952 /* This is the local topology, so there are only F_CONSTR constraints */
1953 if (idef.il[F_CONSTR].nr == 0)
1955 /* There are no constraints,
1956 * we do not need to fill any data structures.
1963 fprintf(debug, "Building the LINCS connectivity\n");
1966 if (DOMAINDECOMP(cr))
1968 if (cr->dd->constraints)
1972 dd_get_constraint_range(cr->dd, &start, &natoms);
1976 natoms = dd_numHomeAtoms(*cr->dd);
1984 at2con = make_at2con(natoms, idef.il, idef.iparams,
1985 flexibleConstraintTreatment(bDynamics));
1987 const int ncon_tot = idef.il[F_CONSTR].nr/3;
1989 /* Ensure we have enough padding for aligned loads for each thread */
1990 const int numEntries = ncon_tot + li->ntask*simd_width;
1991 li->con_index.resize(numEntries);
1992 li->bllen0.resize(numEntries);
1993 li->ddist.resize(numEntries);
1994 li->bla.resize(numEntries*2);
1995 li->blc.resize(numEntries);
1996 li->blc1.resize(numEntries);
1997 li->blnr.resize(numEntries);
1998 li->bllen.resize(numEntries);
1999 li->tmpv.resizeWithPadding(numEntries);
2000 if (DOMAINDECOMP(cr))
2002 li->nlocat.resize(numEntries);
2004 li->tmp1.resize(numEntries);
2005 li->tmp2.resize(numEntries);
2006 li->tmp3.resize(numEntries);
2007 li->tmp4.resize(numEntries);
2008 li->mlambda.resize(numEntries);
2010 iatom = idef.il[F_CONSTR].iatoms;
2012 li->blnr[0] = li->ncc;
2014 /* Assign the constraints for li->ntask LINCS tasks.
2015 * We target a uniform distribution of constraints over the tasks.
2016 * Note that when flexible constraints are present, but are removed here
2017 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
2018 * happen during normal MD, that's ok.
2020 int ncon_assign, ncon_target, con, th;
2022 /* Determine the number of constraints we need to assign here */
2023 ncon_assign = ncon_tot;
2026 /* With energy minimization, flexible constraints are ignored
2027 * (and thus minimized, as they should be).
2029 ncon_assign -= countFlexibleConstraints(idef.il, idef.iparams);
2032 /* Set the target constraint count per task to exactly uniform,
2033 * this might be overridden below.
2035 ncon_target = (ncon_assign + li->ntask - 1)/li->ntask;
2037 /* Mark all constraints as unassigned by setting their index to -1 */
2038 for (con = 0; con < ncon_tot; con++)
2040 li->con_index[con] = -1;
2044 for (th = 0; th < li->ntask; th++)
2048 li_task = &li->task[th];
2050 #if GMX_SIMD_HAVE_REAL
2051 /* With indepedent tasks we likely have H-bond constraints or constraint
2052 * pairs. The connected constraints will be pulled into the task, so the
2053 * constraints per task will often exceed ncon_target.
2054 * Triangle constraints can also increase the count, but there are
2055 * relatively few of those, so we usually expect to get ncon_target.
2059 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
2060 * since otherwise a lot of operations can be wasted.
2061 * There are several ways to round here, we choose the one
2062 * that alternates block sizes, which helps with Intel HT.
2064 ncon_target = ((ncon_assign*(th + 1))/li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1) & ~(GMX_SIMD_REAL_WIDTH - 1);
2066 #endif // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
2068 /* Continue filling the arrays where we left off with the previous task,
2069 * including padding for SIMD.
2071 li_task->b0 = li->nc;
2073 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2075 if (li->con_index[con] == -1)
2080 type = iatom[3*con];
2081 a1 = iatom[3*con + 1];
2082 a2 = iatom[3*con + 2];
2083 lenA = idef.iparams[type].constr.dA;
2084 lenB = idef.iparams[type].constr.dB;
2085 /* Skip the flexible constraints when not doing dynamics */
2086 if (bDynamics || lenA != 0 || lenB != 0)
2088 assign_constraint(li, con, a1, a2, lenA, lenB, &at2con);
2090 if (li->ntask > 1 && !li->bTaskDep)
2092 /* We can generate independent tasks. Check if we
2093 * need to assign connected constraints to our task.
2095 check_assign_connected(li, iatom, idef, bDynamics,
2098 if (li->ntask > 1 && li->ncg_triangle > 0)
2100 /* Ensure constraints in one triangle are assigned
2103 check_assign_triangle(li, iatom, idef, bDynamics,
2104 con, a1, a2, &at2con);
2112 li_task->b1 = li->nc;
2116 /* Copy the last atom pair indices and lengths for constraints
2117 * up to a multiple of simd_width, such that we can do all
2118 * SIMD operations without having to worry about end effects.
2122 li->nc = ((li_task->b1 + simd_width - 1)/simd_width)*simd_width;
2123 last = li_task->b1 - 1;
2124 for (i = li_task->b1; i < li->nc; i++)
2126 li->bla[i*2 ] = li->bla[last*2 ];
2127 li->bla[i*2 + 1] = li->bla[last*2 + 1];
2128 li->bllen0[i] = li->bllen0[last];
2129 li->ddist[i] = li->ddist[last];
2130 li->bllen[i] = li->bllen[last];
2131 li->blnr[i + 1] = li->blnr[last + 1];
2135 /* Keep track of how many constraints we assigned */
2136 li->nc_real += li_task->b1 - li_task->b0;
2140 fprintf(debug, "LINCS task %d constraints %d - %d\n",
2141 th, li_task->b0, li_task->b1);
2145 assert(li->nc_real == ncon_assign);
2149 /* Without DD we order the blbnb matrix to optimize memory access.
2150 * With DD the overhead of sorting is more than the gain during access.
2152 bSortMatrix = !DOMAINDECOMP(cr);
2154 li->blbnb.resize(li->ncc);
2156 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2157 for (th = 0; th < li->ntask; th++)
2163 li_task = &li->task[th];
2165 if (li->ncg_triangle > 0)
2167 /* This is allocating too much, but it is difficult to improve */
2168 li_task->triangle.resize(li_task->b1 - li_task->b0);
2169 li_task->tri_bits.resize(li_task->b1 - li_task->b0);
2172 set_matrix_indices(li, li_task, &at2con, bSortMatrix);
2174 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2177 done_blocka(&at2con);
2179 if (cr->dd == nullptr)
2181 /* Since the matrix is static, we should free some memory */
2182 li->blbnb.resize(li->ncc);
2185 li->blmf.resize(li->ncc);
2186 li->blmf1.resize(li->ncc);
2187 li->tmpncc.resize(li->ncc);
2189 gmx::ArrayRef<const int> nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2190 if (!nlocat_dd.empty())
2192 /* Convert nlocat from local topology to LINCS constraint indexing */
2193 for (con = 0; con < ncon_tot; con++)
2195 li->nlocat[li->con_index[con]] = nlocat_dd[con];
2205 fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n",
2206 li->nc_real, li->nc, li->ncc);
2211 lincs_thread_setup(li, md.nr);
2214 set_lincs_matrix(li, md.invmass, md.lambda);
2217 //! Issues a warning when LINCS constraints cannot be satisfied.
2218 static void lincs_warning(gmx_domdec_t *dd, const rvec *x, rvec *xprime, t_pbc *pbc,
2219 int ncons, gmx::ArrayRef<const int> bla,
2220 gmx::ArrayRef<const real> bllen, real wangle,
2221 int maxwarn, int *warncount)
2225 real wfac, d0, d1, cosine;
2227 wfac = std::cos(DEG2RAD*wangle);
2230 "bonds that rotated more than %g degrees:\n"
2231 " atom 1 atom 2 angle previous, current, constraint length\n",
2234 for (b = 0; b < ncons; b++)
2240 pbc_dx_aiuc(pbc, x[i], x[j], v0);
2241 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2245 rvec_sub(x[i], x[j], v0);
2246 rvec_sub(xprime[i], xprime[j], v1);
2250 cosine = ::iprod(v0, v1)/(d0*d1);
2254 " %6d %6d %5.1f %8.4f %8.4f %8.4f\n",
2255 ddglatnr(dd, i), ddglatnr(dd, j),
2256 RAD2DEG*std::acos(cosine), d0, d1, bllen[b]);
2257 if (!std::isfinite(d1))
2259 gmx_fatal(FARGS, "Bond length not finite.");
2265 if (*warncount > maxwarn)
2267 too_many_constraint_warnings(econtLINCS, *warncount);
2271 //! Determine how well the constraints have been satisfied.
2272 static void cconerr(const Lincs *lincsd,
2273 rvec *x, t_pbc *pbc,
2274 real *ncons_loc, real *ssd, real *max, int *imax)
2276 gmx::ArrayRef<const int> bla = lincsd->bla;
2277 gmx::ArrayRef<const real> bllen = lincsd->bllen;
2278 gmx::ArrayRef<const int> nlocat = lincsd->nlocat;
2284 for (int task = 0; task < lincsd->ntask; task++)
2288 for (b = lincsd->task[task].b0; b < lincsd->task[task].b1; b++)
2295 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
2299 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
2302 len = r2*gmx::invsqrt(r2);
2303 d = std::abs(len/bllen[b]-1);
2304 if (d > ma && (nlocat.empty() || nlocat[b]))
2316 ssd2 += nlocat[b]*d*d;
2322 *ncons_loc = (nlocat.empty() ? 1 : 0.5)*count;
2323 *ssd = (nlocat.empty() ? 1 : 0.5)*ssd2;
2328 bool constrain_lincs(bool computeRmsd,
2329 const t_inputrec &ir,
2331 Lincs *lincsd, const t_mdatoms &md,
2332 const t_commrec *cr,
2333 const gmx_multisim_t *ms,
2334 const rvec *x, rvec *xprime, rvec *min_proj,
2335 matrix box, t_pbc *pbc,
2336 real lambda, real *dvdlambda,
2337 real invdt, rvec *v,
2338 bool bCalcVir, tensor vir_r_m_dr,
2339 ConstraintVariable econq,
2341 int maxwarn, int *warncount)
2344 char buf2[22], buf3[STRLEN];
2346 real ncons_loc, p_ssd, p_max = 0;
2352 /* This boolean should be set by a flag passed to this routine.
2353 * We can also easily check if any constraint length is changed,
2354 * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2356 bCalcDHDL = (ir.efep != efepNO && dvdlambda != nullptr);
2358 if (lincsd->nc == 0 && cr->dd == nullptr)
2362 lincsd->rmsdData = {{0}};
2368 if (econq == ConstraintVariable::Positions)
2370 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2371 * also with efep!=fepNO.
2373 if (ir.efep != efepNO)
2375 if (md.nMassPerturbed && lincsd->matlam != md.lambda)
2377 set_lincs_matrix(lincsd, md.invmass, md.lambda);
2380 for (i = 0; i < lincsd->nc; i++)
2382 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
2386 if (lincsd->ncg_flex)
2388 /* Set the flexible constraint lengths to the old lengths */
2391 for (i = 0; i < lincsd->nc; i++)
2393 if (lincsd->bllen[i] == 0)
2395 pbc_dx_aiuc(pbc, x[lincsd->bla[2*i]], x[lincsd->bla[2*i+1]], dx);
2396 lincsd->bllen[i] = norm(dx);
2402 for (i = 0; i < lincsd->nc; i++)
2404 if (lincsd->bllen[i] == 0)
2407 std::sqrt(distance2(x[lincsd->bla[2*i]],
2408 x[lincsd->bla[2*i+1]]));
2416 cconerr(lincsd, xprime, pbc,
2417 &ncons_loc, &p_ssd, &p_max, &p_imax);
2420 /* This bWarn var can be updated by multiple threads
2421 * at the same time. But as we only need to detect
2422 * if a warning occurred or not, this is not an issue.
2426 /* The OpenMP parallel region of constrain_lincs for coords */
2427 #pragma omp parallel num_threads(lincsd->ntask)
2431 int th = gmx_omp_get_thread_num();
2433 clear_mat(lincsd->task[th].vir_r_m_dr);
2435 do_lincs(x, xprime, box, pbc, lincsd, th,
2438 ir.LincsWarnAngle, &bWarn,
2440 th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2442 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2445 if (debug && lincsd->nc > 0)
2447 fprintf(debug, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2448 fprintf(debug, " Before LINCS %.6f %.6f %6d %6d\n",
2449 std::sqrt(p_ssd/ncons_loc), p_max,
2450 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2451 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2453 if (computeRmsd || debug)
2455 cconerr(lincsd, xprime, pbc,
2456 &ncons_loc, &p_ssd, &p_max, &p_imax);
2457 lincsd->rmsdData[0] = ncons_loc;
2458 lincsd->rmsdData[1] = p_ssd;
2462 lincsd->rmsdData = {{0}};
2464 if (debug && lincsd->nc > 0)
2467 " After LINCS %.6f %.6f %6d %6d\n\n",
2468 std::sqrt(p_ssd/ncons_loc), p_max,
2469 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2470 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2475 if (maxwarn < INT_MAX)
2477 cconerr(lincsd, xprime, pbc,
2478 &ncons_loc, &p_ssd, &p_max, &p_imax);
2481 sprintf(buf3, " in simulation %d", ms->sim);
2488 "\nStep %s, time %g (ps) LINCS WARNING%s\n"
2489 "relative constraint deviation after LINCS:\n"
2490 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2491 gmx_step_str(step, buf2), ir.init_t+step*ir.delta_t,
2493 std::sqrt(p_ssd/ncons_loc), p_max,
2494 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2495 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2497 lincs_warning(cr->dd, x, xprime, pbc,
2498 lincsd->nc, lincsd->bla, lincsd->bllen,
2499 ir.LincsWarnAngle, maxwarn, warncount);
2501 bOK = (p_max < 0.5);
2504 if (lincsd->ncg_flex)
2506 for (i = 0; (i < lincsd->nc); i++)
2508 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2510 lincsd->bllen[i] = 0;
2517 /* The OpenMP parallel region of constrain_lincs for derivatives */
2518 #pragma omp parallel num_threads(lincsd->ntask)
2522 int th = gmx_omp_get_thread_num();
2524 do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
2525 md.invmass, econq, bCalcDHDL,
2526 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2528 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2534 /* Reduce the dH/dlambda contributions over the threads */
2539 for (th = 0; th < lincsd->ntask; th++)
2541 dhdlambda += lincsd->task[th].dhdlambda;
2543 if (econq == ConstraintVariable::Positions)
2545 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2546 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2547 dhdlambda /= ir.delta_t*ir.delta_t;
2549 *dvdlambda += dhdlambda;
2552 if (bCalcVir && lincsd->ntask > 1)
2554 for (i = 1; i < lincsd->ntask; i++)
2556 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2560 /* count assuming nit=1 */
2561 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2562 inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
2563 if (lincsd->ntriangle > 0)
2565 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
2569 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real*2);
2573 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);