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,
235 gmx::ArrayRef<const real> blcc,
236 gmx::ArrayRef<real> rhs1,
237 gmx::ArrayRef<real> rhs2,
238 gmx::ArrayRef<real> sol)
240 gmx::ArrayRef<const int> blnr = lincsd.blnr;
241 gmx::ArrayRef<const int> blbnb = lincsd.blbnb;
243 const int b0 = li_task.b0;
244 const int b1 = li_task.b1;
245 const int nrec = lincsd.nOrder;
247 for (int rec = 0; rec < nrec; rec++)
253 for (int b = b0; b < b1; b++)
259 for (n = blnr[b]; n < blnr[b+1]; n++)
261 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
264 sol[b] = sol[b] + mvb;
267 gmx::ArrayRef<real> swap;
272 } /* nrec*(ncons+2*nrtot) flops */
274 if (lincsd.ntriangle > 0)
276 /* Perform an extra nrec recursions for only the constraints
277 * involved in rigid triangles.
278 * In this way their accuracy should come close to those of the other
279 * constraints, since traingles of constraints can produce eigenvalues
280 * around 0.7, while the effective eigenvalue for bond constraints
281 * is around 0.4 (and 0.7*0.7=0.5).
286 /* We need a barrier here, since other threads might still be
287 * reading the contents of rhs1 and/o rhs2.
288 * We could avoid this barrier by introducing two extra rhs
289 * arrays for the triangle constraints only.
294 /* Constraints involved in a triangle are ensured to be in the same
295 * LINCS task. This means no barriers are required during the extra
296 * iterations for the triangle constraints.
298 gmx::ArrayRef<const int> triangle = li_task.triangle;
299 gmx::ArrayRef<const int> tri_bits = li_task.tri_bits;
301 for (int rec = 0; rec < nrec; rec++)
303 for (int tb = 0; tb < li_task.ntriangle; tb++)
305 int b, bits, nr0, nr1, n;
313 for (n = nr0; n < nr1; n++)
315 if (bits & (1 << (n - nr0)))
317 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
321 sol[b] = sol[b] + mvb;
324 gmx::ArrayRef<real> swap;
329 } /* nrec*(ntriangle + ncc_triangle*2) flops */
331 if (lincsd.bTaskDepTri)
333 /* The constraints triangles are decoupled from each other,
334 * but constraints in one triangle cross thread task borders.
335 * We could probably avoid this with more advanced setup code.
342 //! Update atomic coordinates when an index is not required.
344 lincs_update_atoms_noind(int ncons,
345 gmx::ArrayRef<const int> bla,
347 gmx::ArrayRef<const real> fac,
348 gmx::ArrayRef<const gmx::RVec> r,
352 if (invmass != nullptr)
354 for (int b = 0; b < ncons; b++)
358 real mvb = preFactor*fac[b];
359 real im1 = invmass[i];
360 real im2 = invmass[j];
361 real tmp0 = r[b][0]*mvb;
362 real tmp1 = r[b][1]*mvb;
363 real tmp2 = r[b][2]*mvb;
370 } /* 16 ncons flops */
374 for (int b = 0; b < ncons; b++)
378 real mvb = preFactor*fac[b];
379 real tmp0 = r[b][0]*mvb;
380 real tmp1 = r[b][1]*mvb;
381 real tmp2 = r[b][2]*mvb;
392 //! Update atomic coordinates when an index is required.
394 lincs_update_atoms_ind(gmx::ArrayRef<const int> ind,
395 gmx::ArrayRef<const int> bla,
397 gmx::ArrayRef<const real> fac,
398 gmx::ArrayRef<const gmx::RVec> r,
402 if (invmass != nullptr)
408 real mvb = preFactor*fac[b];
409 real im1 = invmass[i];
410 real im2 = invmass[j];
411 real tmp0 = r[b][0]*mvb;
412 real tmp1 = r[b][1]*mvb;
413 real tmp2 = r[b][2]*mvb;
420 } /* 16 ncons flops */
428 real mvb = preFactor*fac[b];
429 real tmp0 = r[b][0]*mvb;
430 real tmp1 = r[b][1]*mvb;
431 real tmp2 = r[b][2]*mvb;
438 } /* 16 ncons flops */
442 //! Update coordinates for atoms.
444 lincs_update_atoms(Lincs *li,
447 gmx::ArrayRef<const real> fac,
448 gmx::ArrayRef<const gmx::RVec> r,
454 /* Single thread, we simply update for all constraints */
455 lincs_update_atoms_noind(li->nc_real,
456 li->bla, preFactor, fac, r, invmass, x);
460 /* Update the atom vector components for our thread local
461 * constraints that only access our local atom range.
462 * This can be done without a barrier.
464 lincs_update_atoms_ind(li->task[th].ind,
465 li->bla, preFactor, fac, r, invmass, x);
467 if (!li->task[li->ntask].ind.empty())
469 /* Update the constraints that operate on atoms
470 * in multiple thread atom blocks on the master thread.
475 lincs_update_atoms_ind(li->task[li->ntask].ind,
476 li->bla, preFactor, fac, r, invmass, x);
482 #if GMX_SIMD_HAVE_REAL
483 //! Helper function so that we can run TSAN with SIMD support (where implemented).
485 static inline void gmx_simdcall
486 gatherLoadUTransposeTSANSafe(const real *base,
487 const std::int32_t *offset,
492 #if (CMAKE_BUILD_TYPE == CMAKE_BUILD_TYPE_TSAN) && GMX_SIMD_X86_AVX2_256
493 // This function is only implemented in this case
494 gatherLoadUTransposeSafe<align>(base, offset, v0, v1, v2);
496 gatherLoadUTranspose<align>(base, offset, v0, v1, v2);
500 /*! \brief Calculate the constraint distance vectors r to project on from x.
502 * Determine the right-hand side of the matrix equation using quantity f.
503 * This function only differs from calc_dr_x_xp_simd below in that
504 * no constraint length is subtracted and no PBC is used for f. */
505 static void gmx_simdcall
506 calc_dr_x_f_simd(int b0,
508 gmx::ArrayRef<const int> bla,
509 const rvec * gmx_restrict x,
510 const rvec * gmx_restrict f,
511 const real * gmx_restrict blc,
512 const real * pbc_simd,
513 rvec * gmx_restrict r,
514 real * gmx_restrict rhs,
515 real * gmx_restrict sol)
517 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
519 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
521 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
526 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
528 SimdReal x0_S, y0_S, z0_S;
529 SimdReal x1_S, y1_S, z1_S;
530 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
531 SimdReal fx_S, fy_S, fz_S, ip_S, rhs_S;
532 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
533 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
535 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
537 offset0[i] = bla[bs*2 + i*2];
538 offset1[i] = bla[bs*2 + i*2 + 1];
541 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
542 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
547 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
549 n2_S = norm2(rx_S, ry_S, rz_S);
550 il_S = invsqrt(n2_S);
556 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
558 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(f), offset0, &x0_S, &y0_S, &z0_S);
559 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(f), offset1, &x1_S, &y1_S, &z1_S);
564 ip_S = iprod(rx_S, ry_S, rz_S, fx_S, fy_S, fz_S);
566 rhs_S = load<SimdReal>(blc + bs) * ip_S;
568 store(rhs + bs, rhs_S);
569 store(sol + bs, rhs_S);
572 #endif // GMX_SIMD_HAVE_REAL
574 /*! \brief LINCS projection, works on derivatives of the coordinates. */
575 static void do_lincsp(const rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
576 Lincs *lincsd, int th,
578 ConstraintVariable econq, bool bCalcDHDL,
579 bool bCalcVir, tensor rmdf)
581 const int b0 = lincsd->task[th].b0;
582 const int b1 = lincsd->task[th].b1;
584 gmx::ArrayRef<const int> bla = lincsd->bla;
585 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
586 gmx::ArrayRef<const int> blnr = lincsd->blnr;
587 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
589 gmx::ArrayRef<const real> blc;
590 gmx::ArrayRef<const real> blmf;
591 if (econq != ConstraintVariable::Force)
593 /* Use mass-weighted parameters */
599 /* Use non mass-weighted parameters */
601 blmf = lincsd->blmf1;
603 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
604 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
605 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
606 gmx::ArrayRef<real> sol = lincsd->tmp3;
608 #if GMX_SIMD_HAVE_REAL
609 /* This SIMD code does the same as the plain-C code after the #else.
610 * The only difference is that we always call pbc code, as with SIMD
611 * the overhead of pbc computation (when not needed) is small.
613 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
615 /* Convert the pbc struct for SIMD */
616 set_pbc_simd(pbc, pbc_simd);
618 /* Compute normalized x i-j vectors, store in r.
619 * Compute the inner product of r and xp i-j and store in rhs1.
621 calc_dr_x_f_simd(b0, b1, bla, x, f, blc.data(),
623 as_rvec_array(r.data()), rhs1.data(), sol.data());
625 #else // GMX_SIMD_HAVE_REAL
627 /* Compute normalized i-j vectors */
630 for (int b = b0; b < b1; b++)
634 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
640 for (int b = b0; b < b1; b++)
644 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
646 } /* 16 ncons flops */
649 for (int b = b0; b < b1; b++)
653 real mvb = blc[b]*(r[b][0]*(f[i][0] - f[j][0]) +
654 r[b][1]*(f[i][1] - f[j][1]) +
655 r[b][2]*(f[i][2] - f[j][2]));
661 #endif // GMX_SIMD_HAVE_REAL
663 if (lincsd->bTaskDep)
665 /* We need a barrier, since the matrix construction below
666 * can access entries in r of other threads.
671 /* Construct the (sparse) LINCS matrix */
672 for (int b = b0; b < b1; b++)
674 for (int n = blnr[b]; n < blnr[b+1]; n++)
676 blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
679 /* Together: 23*ncons + 6*nrtot flops */
681 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
682 /* nrec*(ncons+2*nrtot) flops */
684 if (econq == ConstraintVariable::Deriv_FlexCon)
686 /* We only want to constraint the flexible constraints,
687 * so we mask out the normal ones by setting sol to 0.
689 for (int b = b0; b < b1; b++)
691 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
698 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
699 for (int b = b0; b < b1; b++)
704 /* When constraining forces, we should not use mass weighting,
705 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
707 lincs_update_atoms(lincsd, th, 1.0, sol, r,
708 (econq != ConstraintVariable::Force) ? invmass : nullptr, fp);
713 for (int b = b0; b < b1; b++)
715 dhdlambda -= sol[b]*lincsd->ddist[b];
718 lincsd->task[th].dhdlambda = dhdlambda;
723 /* Constraint virial,
724 * determines sum r_bond x delta f,
725 * where delta f is the constraint correction
726 * of the quantity that is being constrained.
728 for (int b = b0; b < b1; b++)
730 const real mvb = lincsd->bllen[b]*sol[b];
731 for (int i = 0; i < DIM; i++)
733 const real tmp1 = mvb*r[b][i];
734 for (int j = 0; j < DIM; j++)
736 rmdf[i][j] += tmp1*r[b][j];
739 } /* 23 ncons flops */
743 #if GMX_SIMD_HAVE_REAL
745 /*! \brief Calculate the constraint distance vectors r to project on from x.
747 * Determine the right-hand side of the matrix equation using coordinates xp. */
748 static void gmx_simdcall
749 calc_dr_x_xp_simd(int b0,
751 gmx::ArrayRef<const int> bla,
752 const rvec * gmx_restrict x,
753 const rvec * gmx_restrict xp,
754 const real * gmx_restrict bllen,
755 const real * gmx_restrict blc,
756 const real * pbc_simd,
757 rvec * gmx_restrict r,
758 real * gmx_restrict rhs,
759 real * gmx_restrict sol)
761 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
762 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
764 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
769 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
771 SimdReal x0_S, y0_S, z0_S;
772 SimdReal x1_S, y1_S, z1_S;
773 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
774 SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
775 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
776 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
778 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
780 offset0[i] = bla[bs*2 + i*2];
781 offset1[i] = bla[bs*2 + i*2 + 1];
784 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
785 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
790 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
792 n2_S = norm2(rx_S, ry_S, rz_S);
793 il_S = invsqrt(n2_S);
799 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
801 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(xp), offset0, &x0_S, &y0_S, &z0_S);
802 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(xp), offset1, &x1_S, &y1_S, &z1_S);
808 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
810 ip_S = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
812 rhs_S = load<SimdReal>(blc + bs) * (ip_S - load<SimdReal>(bllen + bs));
814 store(rhs + bs, rhs_S);
815 store(sol + bs, rhs_S);
818 #endif // GMX_SIMD_HAVE_REAL
820 /*! \brief Determine the distances and right-hand side for the next iteration. */
821 gmx_unused static void calc_dist_iter(
824 gmx::ArrayRef<const int> bla,
825 const rvec * gmx_restrict xp,
826 const real * gmx_restrict bllen,
827 const real * gmx_restrict blc,
830 real * gmx_restrict rhs,
831 real * gmx_restrict sol,
836 for (b = b0; b < b1; b++)
838 real len, len2, dlen2, mvb;
844 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
848 rvec_sub(xp[bla[2*b]], xp[bla[2*b+1]], dx);
851 dlen2 = 2*len2 - ::norm2(dx);
852 if (dlen2 < wfac*len2)
854 /* not race free - see detailed comment in caller */
859 mvb = blc[b]*(len - dlen2*gmx::invsqrt(dlen2));
867 } /* 20*ncons flops */
870 #if GMX_SIMD_HAVE_REAL
871 /*! \brief As calc_dist_iter(), but using SIMD intrinsics. */
872 static void gmx_simdcall
873 calc_dist_iter_simd(int b0,
875 gmx::ArrayRef<const int> bla,
876 const rvec * gmx_restrict x,
877 const real * gmx_restrict bllen,
878 const real * gmx_restrict blc,
879 const real * pbc_simd,
881 real * gmx_restrict rhs,
882 real * gmx_restrict sol,
885 SimdReal min_S(GMX_REAL_MIN);
887 SimdReal wfac_S(wfac);
892 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
894 /* Initialize all to FALSE */
895 warn_S = (two_S < setZero());
897 for (bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
899 SimdReal x0_S, y0_S, z0_S;
900 SimdReal x1_S, y1_S, z1_S;
901 SimdReal rx_S, ry_S, rz_S, n2_S;
902 SimdReal len_S, len2_S, dlen2_S, lc_S, blc_S;
903 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
904 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
906 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
908 offset0[i] = bla[bs*2 + i*2];
909 offset1[i] = bla[bs*2 + i*2 + 1];
912 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
913 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
919 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
921 n2_S = norm2(rx_S, ry_S, rz_S);
923 len_S = load<SimdReal>(bllen + bs);
924 len2_S = len_S * len_S;
926 dlen2_S = fms(two_S, len2_S, n2_S);
928 warn_S = warn_S || (dlen2_S < (wfac_S * len2_S));
930 /* Avoid 1/0 by taking the max with REAL_MIN.
931 * Note: when dlen2 is close to zero (90 degree constraint rotation),
932 * the accuracy of the algorithm is no longer relevant.
934 dlen2_S = max(dlen2_S, min_S);
936 lc_S = fnma(dlen2_S, invsqrt(dlen2_S), len_S);
938 blc_S = load<SimdReal>(blc + bs);
942 store(rhs + bs, lc_S);
943 store(sol + bs, lc_S);
951 #endif // GMX_SIMD_HAVE_REAL
953 //! Implements LINCS constraining.
954 static void do_lincs(const rvec *x, rvec *xp, matrix box, t_pbc *pbc,
955 Lincs *lincsd, int th,
959 real wangle, bool *bWarn,
960 real invdt, rvec * gmx_restrict v,
961 bool bCalcVir, tensor vir_r_m_dr)
963 const int b0 = lincsd->task[th].b0;
964 const int b1 = lincsd->task[th].b1;
966 gmx::ArrayRef<const int> bla = lincsd->bla;
967 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
968 gmx::ArrayRef<const int> blnr = lincsd->blnr;
969 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
970 gmx::ArrayRef<const real> blc = lincsd->blc;
971 gmx::ArrayRef<const real> blmf = lincsd->blmf;
972 gmx::ArrayRef<const real> bllen = lincsd->bllen;
973 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
974 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
975 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
976 gmx::ArrayRef<real> sol = lincsd->tmp3;
977 gmx::ArrayRef<real> blc_sol = lincsd->tmp4;
978 gmx::ArrayRef<real> mlambda = lincsd->mlambda;
979 gmx::ArrayRef<const int> nlocat = lincsd->nlocat;
981 #if GMX_SIMD_HAVE_REAL
983 /* This SIMD code does the same as the plain-C code after the #else.
984 * The only difference is that we always call pbc code, as with SIMD
985 * the overhead of pbc computation (when not needed) is small.
987 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
989 /* Convert the pbc struct for SIMD */
990 set_pbc_simd(pbc, pbc_simd);
992 /* Compute normalized x i-j vectors, store in r.
993 * Compute the inner product of r and xp i-j and store in rhs1.
995 calc_dr_x_xp_simd(b0, b1, bla, x, xp, bllen.data(), blc.data(),
997 as_rvec_array(r.data()), rhs1.data(), sol.data());
999 #else // GMX_SIMD_HAVE_REAL
1003 /* Compute normalized i-j vectors */
1004 for (int b = b0; b < b1; b++)
1007 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
1010 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
1011 real mvb = blc[b]*(::iprod(r[b], dx) - bllen[b]);
1018 /* Compute normalized i-j vectors */
1019 for (int b = b0; b < b1; b++)
1023 real tmp0 = x[i][0] - x[j][0];
1024 real tmp1 = x[i][1] - x[j][1];
1025 real tmp2 = x[i][2] - x[j][2];
1026 real rlen = gmx::invsqrt(tmp0*tmp0 + tmp1*tmp1 + tmp2*tmp2);
1027 r[b][0] = rlen*tmp0;
1028 r[b][1] = rlen*tmp1;
1029 r[b][2] = rlen*tmp2;
1030 /* 16 ncons flops */
1032 real mvb = blc[b]*(r[b][0]*(xp[i][0] - xp[j][0]) +
1033 r[b][1]*(xp[i][1] - xp[j][1]) +
1034 r[b][2]*(xp[i][2] - xp[j][2]) - bllen[b]);
1039 /* Together: 26*ncons + 6*nrtot flops */
1042 #endif // GMX_SIMD_HAVE_REAL
1044 if (lincsd->bTaskDep)
1046 /* We need a barrier, since the matrix construction below
1047 * can access entries in r of other threads.
1052 /* Construct the (sparse) LINCS matrix */
1053 for (int b = b0; b < b1; b++)
1055 for (int n = blnr[b]; n < blnr[b+1]; n++)
1057 blcc[n] = blmf[n]*gmx::dot(r[b], r[blbnb[n]]);
1060 /* Together: 26*ncons + 6*nrtot flops */
1062 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
1063 /* nrec*(ncons+2*nrtot) flops */
1065 #if GMX_SIMD_HAVE_REAL
1066 for (int b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1068 SimdReal t1 = load<SimdReal>(blc.data() + b);
1069 SimdReal t2 = load<SimdReal>(sol.data() + b);
1070 store(mlambda.data() + b, t1 * t2);
1073 for (int b = b0; b < b1; b++)
1075 mlambda[b] = blc[b]*sol[b];
1077 #endif // GMX_SIMD_HAVE_REAL
1079 /* Update the coordinates */
1080 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1083 ******** Correction for centripetal effects ********
1088 wfac = std::cos(DEG2RAD*wangle);
1091 for (int iter = 0; iter < lincsd->nIter; iter++)
1093 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
1098 /* Communicate the corrected non-local coordinates */
1099 if (DOMAINDECOMP(cr))
1101 dd_move_x_constraints(cr->dd, box, xp, nullptr, FALSE);
1106 else if (lincsd->bTaskDep)
1111 #if GMX_SIMD_HAVE_REAL
1112 calc_dist_iter_simd(b0, b1, bla,
1113 xp, bllen.data(), blc.data(), pbc_simd, wfac,
1114 rhs1.data(), sol.data(), bWarn);
1116 calc_dist_iter(b0, b1, bla, xp,
1117 bllen.data(), blc.data(), pbc, wfac,
1118 rhs1.data(), sol.data(), bWarn);
1119 /* 20*ncons flops */
1120 #endif // GMX_SIMD_HAVE_REAL
1122 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
1123 /* nrec*(ncons+2*nrtot) flops */
1125 #if GMX_SIMD_HAVE_REAL
1126 for (int b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1128 SimdReal t1 = load<SimdReal>(blc.data() + b);
1129 SimdReal t2 = load<SimdReal>(sol.data() + b);
1130 SimdReal mvb = t1 * t2;
1131 store(blc_sol.data() + b, mvb);
1132 store(mlambda.data() + b, load<SimdReal>(mlambda.data() + b) + mvb);
1135 for (int b = b0; b < b1; b++)
1137 real mvb = blc[b]*sol[b];
1141 #endif // GMX_SIMD_HAVE_REAL
1143 /* Update the coordinates */
1144 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1146 /* nit*ncons*(37+9*nrec) flops */
1150 /* Update the velocities */
1151 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1152 /* 16 ncons flops */
1155 if (!nlocat.empty() && (bCalcDHDL || bCalcVir))
1157 if (lincsd->bTaskDep)
1159 /* In lincs_update_atoms threads might cross-read mlambda */
1163 /* Only account for local atoms */
1164 for (int b = b0; b < b1; b++)
1166 mlambda[b] *= 0.5*nlocat[b];
1173 for (int b = b0; b < b1; b++)
1175 /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1176 * later after the contributions are reduced over the threads.
1178 dhdl -= lincsd->mlambda[b]*lincsd->ddist[b];
1180 lincsd->task[th].dhdlambda = dhdl;
1185 /* Constraint virial */
1186 for (int b = b0; b < b1; b++)
1188 real tmp0 = -bllen[b]*mlambda[b];
1189 for (int i = 0; i < DIM; i++)
1191 real tmp1 = tmp0*r[b][i];
1192 for (int j = 0; j < DIM; j++)
1194 vir_r_m_dr[i][j] -= tmp1*r[b][j];
1197 } /* 22 ncons flops */
1201 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1202 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1204 * (26+nrec)*ncons + (6+2*nrec)*nrtot
1205 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1207 * (63+nrec)*ncons + (6+4*nrec)*nrtot
1211 /*! \brief Sets the elements in the LINCS matrix for task task. */
1212 static void set_lincs_matrix_task(Lincs *li,
1214 const real *invmass,
1216 int *nCrossTaskTriangles)
1220 /* Construct the coupling coefficient matrix blmf */
1221 li_task->ntriangle = 0;
1223 *nCrossTaskTriangles = 0;
1224 for (i = li_task->b0; i < li_task->b1; i++)
1229 a2 = li->bla[2*i+1];
1230 for (n = li->blnr[i]; (n < li->blnr[i+1]); n++)
1232 int k, sign, center, end;
1236 /* If we are using multiple, independent tasks for LINCS,
1237 * the calls to check_assign_connected should have
1238 * put all connected constraints in our task.
1240 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1242 if (a1 == li->bla[2*k] || a2 == li->bla[2*k+1])
1250 if (a1 == li->bla[2*k] || a1 == li->bla[2*k+1])
1260 li->blmf[n] = sign*invmass[center]*li->blc[i]*li->blc[k];
1261 li->blmf1[n] = sign*0.5;
1262 if (li->ncg_triangle > 0)
1266 /* Look for constraint triangles */
1267 for (nk = li->blnr[k]; (nk < li->blnr[k+1]); nk++)
1270 if (kk != i && kk != k &&
1271 (li->bla[2*kk] == end || li->bla[2*kk+1] == end))
1273 /* Check if the constraints in this triangle actually
1274 * belong to a different task. We still assign them
1275 * here, since it's convenient for the triangle
1276 * iterations, but we then need an extra barrier.
1278 if (k < li_task->b0 || k >= li_task->b1 ||
1279 kk < li_task->b0 || kk >= li_task->b1)
1281 (*nCrossTaskTriangles)++;
1284 if (li_task->ntriangle == 0 ||
1285 li_task->triangle[li_task->ntriangle - 1] < i)
1287 /* Add this constraint to the triangle list */
1288 li_task->triangle[li_task->ntriangle] = i;
1289 li_task->tri_bits[li_task->ntriangle] = 0;
1290 li_task->ntriangle++;
1291 if (li->blnr[i+1] - li->blnr[i] > static_cast<int>(sizeof(li_task->tri_bits[0])*8 - 1))
1293 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %zu allowed for constraints participating in triangles",
1294 li->blnr[i+1] - li->blnr[i],
1295 sizeof(li_task->tri_bits[0])*8-1);
1298 li_task->tri_bits[li_task->ntriangle-1] |= (1 << (n - li->blnr[i]));
1307 /*! \brief Sets the elements in the LINCS matrix. */
1308 static void set_lincs_matrix(Lincs *li, real *invmass, real lambda)
1311 const real invsqrt2 = 0.7071067811865475244;
1313 for (i = 0; (i < li->nc); i++)
1318 a2 = li->bla[2*i+1];
1319 li->blc[i] = gmx::invsqrt(invmass[a1] + invmass[a2]);
1320 li->blc1[i] = invsqrt2;
1323 /* Construct the coupling coefficient matrix blmf */
1324 int th, ntriangle = 0, ncc_triangle = 0, nCrossTaskTriangles = 0;
1325 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle, nCrossTaskTriangles) num_threads(li->ntask) schedule(static)
1326 for (th = 0; th < li->ntask; th++)
1330 set_lincs_matrix_task(li, &li->task[th], invmass,
1331 &ncc_triangle, &nCrossTaskTriangles);
1332 ntriangle += li->task[th].ntriangle;
1334 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1336 li->ntriangle = ntriangle;
1337 li->ncc_triangle = ncc_triangle;
1338 li->bTaskDepTri = (nCrossTaskTriangles > 0);
1342 fprintf(debug, "The %d constraints participate in %d triangles\n",
1343 li->nc, li->ntriangle);
1344 fprintf(debug, "There are %d constraint couplings, of which %d in triangles\n",
1345 li->ncc, li->ncc_triangle);
1346 if (li->ntriangle > 0 && li->ntask > 1)
1348 fprintf(debug, "%d constraint triangles contain constraints assigned to different tasks\n",
1349 nCrossTaskTriangles);
1354 * so we know with which lambda value the masses have been set.
1356 li->matlam = lambda;
1359 //! Finds all triangles of atoms that share constraints to a central atom.
1360 static int count_triangle_constraints(const InteractionLists &ilist,
1361 const t_blocka &at2con)
1363 int ncon1, ncon_tot;
1364 int c0, n1, c1, ac1, n2, c2;
1367 ncon1 = ilist[F_CONSTR].size()/3;
1368 ncon_tot = ncon1 + ilist[F_CONSTRNC].size()/3;
1370 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1371 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1374 for (c0 = 0; c0 < ncon_tot; c0++)
1376 bool bTriangle = FALSE;
1377 const int *iap = constr_iatomptr(ia1, ia2, c0);
1378 const int a00 = iap[1];
1379 const int a01 = iap[2];
1380 for (n1 = at2con.index[a01]; n1 < at2con.index[a01+1]; n1++)
1385 const int *iap = constr_iatomptr(ia1, ia2, c1);
1386 const int a10 = iap[1];
1387 const int a11 = iap[2];
1396 for (n2 = at2con.index[ac1]; n2 < at2con.index[ac1+1]; n2++)
1399 if (c2 != c0 && c2 != c1)
1401 const int *iap = constr_iatomptr(ia1, ia2, c2);
1402 const int a20 = iap[1];
1403 const int a21 = iap[2];
1404 if (a20 == a00 || a21 == a00)
1418 return ncon_triangle;
1421 //! Finds sequences of sequential constraints.
1422 static bool more_than_two_sequential_constraints(const InteractionLists &ilist,
1423 const t_blocka &at2con)
1425 int ncon1, ncon_tot, c;
1426 bool bMoreThanTwoSequentialConstraints;
1428 ncon1 = ilist[F_CONSTR].size()/3;
1429 ncon_tot = ncon1 + ilist[F_CONSTRNC].size()/3;
1431 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1432 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1434 bMoreThanTwoSequentialConstraints = FALSE;
1435 for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
1437 const int *iap = constr_iatomptr(ia1, ia2, c);
1438 const int a1 = iap[1];
1439 const int a2 = iap[2];
1440 /* Check if this constraint has constraints connected at both atoms */
1441 if (at2con.index[a1+1] - at2con.index[a1] > 1 &&
1442 at2con.index[a2+1] - at2con.index[a2] > 1)
1444 bMoreThanTwoSequentialConstraints = TRUE;
1448 return bMoreThanTwoSequentialConstraints;
1451 Lincs *init_lincs(FILE *fplog, const gmx_mtop_t &mtop,
1452 int nflexcon_global, ArrayRef<const t_blocka> at2con,
1453 bool bPLINCS, int nIter, int nProjOrder)
1455 // TODO this should become a unique_ptr
1457 bool bMoreThanTwoSeq;
1461 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
1462 bPLINCS ? " Parallel" : "");
1468 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1469 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1470 li->ncg_flex = nflexcon_global;
1473 li->nOrder = nProjOrder;
1475 li->max_connect = 0;
1476 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1478 for (int a = 0; a < mtop.moltype[mt].atoms.nr; a++)
1480 li->max_connect = std::max(li->max_connect,
1481 at2con[mt].index[a + 1] - at2con[mt].index[a]);
1485 li->ncg_triangle = 0;
1486 bMoreThanTwoSeq = FALSE;
1487 for (const gmx_molblock_t &molb : mtop.molblock)
1489 const gmx_moltype_t &molt = mtop.moltype[molb.type];
1493 count_triangle_constraints(molt.ilist, at2con[molb.type]);
1495 if (!bMoreThanTwoSeq &&
1496 more_than_two_sequential_constraints(molt.ilist, at2con[molb.type]))
1498 bMoreThanTwoSeq = TRUE;
1502 /* Check if we need to communicate not only before LINCS,
1503 * but also before each iteration.
1504 * The check for only two sequential constraints is only
1505 * useful for the common case of H-bond only constraints.
1506 * With more effort we could also make it useful for small
1507 * molecules with nr. sequential constraints <= nOrder-1.
1509 li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1511 if (debug && bPLINCS)
1513 fprintf(debug, "PLINCS communication before each iteration: %d\n",
1514 static_cast<int>(li->bCommIter));
1517 /* LINCS can run on any number of threads.
1518 * Currently the number is fixed for the whole simulation,
1519 * but it could be set in set_lincs().
1520 * The current constraint to task assignment code can create independent
1521 * tasks only when not more than two constraints are connected sequentially.
1523 li->ntask = gmx_omp_nthreads_get(emntLINCS);
1524 li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1527 fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n",
1528 li->ntask, li->bTaskDep ? "" : "in");
1536 /* Allocate an extra elements for "task-overlap" constraints */
1537 li->task.resize(li->ntask + 1);
1540 if (bPLINCS || li->ncg_triangle > 0)
1542 please_cite(fplog, "Hess2008a");
1546 please_cite(fplog, "Hess97a");
1551 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1554 fprintf(fplog, "There are constraints between atoms in different decomposition domains,\n"
1555 "will communicate selected coordinates each lincs iteration\n");
1557 if (li->ncg_triangle > 0)
1560 "%d constraints are involved in constraint triangles,\n"
1561 "will apply an additional matrix expansion of order %d for couplings\n"
1562 "between constraints inside triangles\n",
1563 li->ncg_triangle, li->nOrder);
1570 void done_lincs(Lincs *li)
1575 /*! \brief Sets up the work division over the threads. */
1576 static void lincs_thread_setup(Lincs *li, int natoms)
1578 li->atf.resize(natoms);
1580 gmx::ArrayRef<gmx_bitmask_t> atf = li->atf;
1582 /* Clear the atom flags */
1583 for (gmx_bitmask_t &mask : atf)
1585 bitmask_clear(&mask);
1588 if (li->ntask > BITMASK_SIZE)
1590 gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1593 for (int th = 0; th < li->ntask; th++)
1595 const Task *li_task = &li->task[th];
1597 /* For each atom set a flag for constraints from each */
1598 for (int b = li_task->b0; b < li_task->b1; b++)
1600 bitmask_set_bit(&atf[li->bla[b*2 ]], th);
1601 bitmask_set_bit(&atf[li->bla[b*2 + 1]], th);
1605 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1606 for (int th = 0; th < li->ntask; th++)
1614 li_task = &li->task[th];
1616 bitmask_init_low_bits(&mask, th);
1618 li_task->ind.clear();
1619 li_task->ind_r.clear();
1620 for (b = li_task->b0; b < li_task->b1; b++)
1622 /* We let the constraint with the lowest thread index
1623 * operate on atoms with constraints from multiple threads.
1625 if (bitmask_is_disjoint(atf[li->bla[b*2]], mask) &&
1626 bitmask_is_disjoint(atf[li->bla[b*2+1]], mask))
1628 /* Add the constraint to the local atom update index */
1629 li_task->ind.push_back(b);
1633 /* Add the constraint to the rest block */
1634 li_task->ind_r.push_back(b);
1638 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1641 /* We need to copy all constraints which have not be assigned
1642 * to a thread to a separate list which will be handled by one thread.
1644 Task *li_m = &li->task[li->ntask];
1647 for (int th = 0; th < li->ntask; th++)
1649 const Task &li_task = li->task[th];
1651 for (int ind_r : li_task.ind_r)
1653 li_m->ind.push_back(ind_r);
1658 fprintf(debug, "LINCS thread %d: %zu constraints\n",
1659 th, li_task.ind.size());
1665 fprintf(debug, "LINCS thread r: %zu constraints\n",
1670 //! Assign a constraint.
1671 static void assign_constraint(Lincs *li,
1672 int constraint_index,
1674 real lenA, real lenB,
1675 const t_blocka *at2con)
1681 /* Make an mapping of local topology constraint index to LINCS index */
1682 li->con_index[constraint_index] = con;
1684 li->bllen0[con] = lenA;
1685 li->ddist[con] = lenB - lenA;
1686 /* Set the length to the topology A length */
1687 li->bllen[con] = lenA;
1688 li->bla[2*con] = a1;
1689 li->bla[2*con+1] = a2;
1691 /* Make space in the constraint connection matrix for constraints
1692 * connected to both end of the current constraint.
1695 at2con->index[a1 + 1] - at2con->index[a1] - 1 +
1696 at2con->index[a2 + 1] - at2con->index[a2] - 1;
1698 li->blnr[con + 1] = li->ncc;
1700 /* Increase the constraint count */
1704 /*! \brief Check if constraint with topology index constraint_index is connected
1705 * to other constraints, and if so add those connected constraints to our task. */
1706 static void check_assign_connected(Lincs *li,
1707 const t_iatom *iatom,
1711 const t_blocka *at2con)
1713 /* Currently this function only supports constraint groups
1714 * in which all constraints share at least one atom
1715 * (e.g. H-bond constraints).
1716 * Check both ends of the current constraint for
1717 * connected constraints. We need to assign those
1722 for (end = 0; end < 2; end++)
1726 a = (end == 0 ? a1 : a2);
1728 for (k = at2con->index[a]; k < at2con->index[a + 1]; k++)
1733 /* Check if constraint cc has not yet been assigned */
1734 if (li->con_index[cc] == -1)
1740 lenA = idef.iparams[type].constr.dA;
1741 lenB = idef.iparams[type].constr.dB;
1743 if (bDynamics || lenA != 0 || lenB != 0)
1745 assign_constraint(li, cc, iatom[3*cc + 1], iatom[3*cc + 2], lenA, lenB, at2con);
1752 /*! \brief Check if constraint with topology index constraint_index is involved
1753 * in a constraint triangle, and if so add the other two constraints
1754 * in the triangle to our task. */
1755 static void check_assign_triangle(Lincs *li,
1756 const t_iatom *iatom,
1759 int constraint_index,
1761 const t_blocka *at2con)
1763 int nca, cc[32], ca[32], k;
1764 int c_triangle[2] = { -1, -1 };
1767 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1772 if (c != constraint_index)
1776 aa1 = iatom[c*3 + 1];
1777 aa2 = iatom[c*3 + 2];
1793 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1798 if (c != constraint_index)
1802 aa1 = iatom[c*3 + 1];
1803 aa2 = iatom[c*3 + 2];
1806 for (i = 0; i < nca; i++)
1810 c_triangle[0] = cc[i];
1817 for (i = 0; i < nca; i++)
1821 c_triangle[0] = cc[i];
1829 if (c_triangle[0] >= 0)
1833 for (end = 0; end < 2; end++)
1835 /* Check if constraint c_triangle[end] has not yet been assigned */
1836 if (li->con_index[c_triangle[end]] == -1)
1841 i = c_triangle[end]*3;
1843 lenA = idef.iparams[type].constr.dA;
1844 lenB = idef.iparams[type].constr.dB;
1846 if (bDynamics || lenA != 0 || lenB != 0)
1848 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1855 //! Sets matrix indices.
1856 static void set_matrix_indices(Lincs *li,
1857 const Task *li_task,
1858 const t_blocka *at2con,
1863 for (b = li_task->b0; b < li_task->b1; b++)
1868 a2 = li->bla[b*2 + 1];
1871 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1875 concon = li->con_index[at2con->a[k]];
1878 li->blbnb[i++] = concon;
1881 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1885 concon = li->con_index[at2con->a[k]];
1888 li->blbnb[i++] = concon;
1894 /* Order the blbnb matrix to optimize memory access */
1895 std::sort(&(li->blbnb[li->blnr[b]]), &(li->blbnb[li->blnr[b+1]]));
1900 void set_lincs(const t_idef &idef,
1901 const t_mdatoms &md,
1903 const t_commrec *cr,
1913 /* Zero the thread index ranges.
1914 * Otherwise without local constraints we could return with old ranges.
1916 for (int i = 0; i < li->ntask; i++)
1920 li->task[i].ind.clear();
1924 li->task[li->ntask].ind.clear();
1927 /* This is the local topology, so there are only F_CONSTR constraints */
1928 if (idef.il[F_CONSTR].nr == 0)
1930 /* There are no constraints,
1931 * we do not need to fill any data structures.
1938 fprintf(debug, "Building the LINCS connectivity\n");
1941 if (DOMAINDECOMP(cr))
1943 if (cr->dd->constraints)
1947 dd_get_constraint_range(cr->dd, &start, &natoms);
1951 natoms = dd_numHomeAtoms(*cr->dd);
1959 at2con = make_at2con(natoms, idef.il, idef.iparams,
1960 flexibleConstraintTreatment(bDynamics));
1962 const int ncon_tot = idef.il[F_CONSTR].nr/3;
1964 /* Ensure we have enough padding for aligned loads for each thread */
1965 const int numEntries = ncon_tot + li->ntask*simd_width;
1966 li->con_index.resize(numEntries);
1967 li->bllen0.resize(numEntries);
1968 li->ddist.resize(numEntries);
1969 li->bla.resize(numEntries*2);
1970 li->blc.resize(numEntries);
1971 li->blc1.resize(numEntries);
1972 li->blnr.resize(numEntries);
1973 li->bllen.resize(numEntries);
1974 li->tmpv.resizeWithPadding(numEntries);
1975 if (DOMAINDECOMP(cr))
1977 li->nlocat.resize(numEntries);
1979 li->tmp1.resize(numEntries);
1980 li->tmp2.resize(numEntries);
1981 li->tmp3.resize(numEntries);
1982 li->tmp4.resize(numEntries);
1983 li->mlambda.resize(numEntries);
1985 iatom = idef.il[F_CONSTR].iatoms;
1987 li->blnr[0] = li->ncc;
1989 /* Assign the constraints for li->ntask LINCS tasks.
1990 * We target a uniform distribution of constraints over the tasks.
1991 * Note that when flexible constraints are present, but are removed here
1992 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
1993 * happen during normal MD, that's ok.
1995 int ncon_assign, ncon_target, con, th;
1997 /* Determine the number of constraints we need to assign here */
1998 ncon_assign = ncon_tot;
2001 /* With energy minimization, flexible constraints are ignored
2002 * (and thus minimized, as they should be).
2004 ncon_assign -= countFlexibleConstraints(idef.il, idef.iparams);
2007 /* Set the target constraint count per task to exactly uniform,
2008 * this might be overridden below.
2010 ncon_target = (ncon_assign + li->ntask - 1)/li->ntask;
2012 /* Mark all constraints as unassigned by setting their index to -1 */
2013 for (con = 0; con < ncon_tot; con++)
2015 li->con_index[con] = -1;
2019 for (th = 0; th < li->ntask; th++)
2023 li_task = &li->task[th];
2025 #if GMX_SIMD_HAVE_REAL
2026 /* With indepedent tasks we likely have H-bond constraints or constraint
2027 * pairs. The connected constraints will be pulled into the task, so the
2028 * constraints per task will often exceed ncon_target.
2029 * Triangle constraints can also increase the count, but there are
2030 * relatively few of those, so we usually expect to get ncon_target.
2034 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
2035 * since otherwise a lot of operations can be wasted.
2036 * There are several ways to round here, we choose the one
2037 * that alternates block sizes, which helps with Intel HT.
2039 ncon_target = ((ncon_assign*(th + 1))/li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1) & ~(GMX_SIMD_REAL_WIDTH - 1);
2041 #endif // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
2043 /* Continue filling the arrays where we left off with the previous task,
2044 * including padding for SIMD.
2046 li_task->b0 = li->nc;
2048 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2050 if (li->con_index[con] == -1)
2055 type = iatom[3*con];
2056 a1 = iatom[3*con + 1];
2057 a2 = iatom[3*con + 2];
2058 lenA = idef.iparams[type].constr.dA;
2059 lenB = idef.iparams[type].constr.dB;
2060 /* Skip the flexible constraints when not doing dynamics */
2061 if (bDynamics || lenA != 0 || lenB != 0)
2063 assign_constraint(li, con, a1, a2, lenA, lenB, &at2con);
2065 if (li->ntask > 1 && !li->bTaskDep)
2067 /* We can generate independent tasks. Check if we
2068 * need to assign connected constraints to our task.
2070 check_assign_connected(li, iatom, idef, bDynamics,
2073 if (li->ntask > 1 && li->ncg_triangle > 0)
2075 /* Ensure constraints in one triangle are assigned
2078 check_assign_triangle(li, iatom, idef, bDynamics,
2079 con, a1, a2, &at2con);
2087 li_task->b1 = li->nc;
2091 /* Copy the last atom pair indices and lengths for constraints
2092 * up to a multiple of simd_width, such that we can do all
2093 * SIMD operations without having to worry about end effects.
2097 li->nc = ((li_task->b1 + simd_width - 1)/simd_width)*simd_width;
2098 last = li_task->b1 - 1;
2099 for (i = li_task->b1; i < li->nc; i++)
2101 li->bla[i*2 ] = li->bla[last*2 ];
2102 li->bla[i*2 + 1] = li->bla[last*2 + 1];
2103 li->bllen0[i] = li->bllen0[last];
2104 li->ddist[i] = li->ddist[last];
2105 li->bllen[i] = li->bllen[last];
2106 li->blnr[i + 1] = li->blnr[last + 1];
2110 /* Keep track of how many constraints we assigned */
2111 li->nc_real += li_task->b1 - li_task->b0;
2115 fprintf(debug, "LINCS task %d constraints %d - %d\n",
2116 th, li_task->b0, li_task->b1);
2120 assert(li->nc_real == ncon_assign);
2124 /* Without DD we order the blbnb matrix to optimize memory access.
2125 * With DD the overhead of sorting is more than the gain during access.
2127 bSortMatrix = !DOMAINDECOMP(cr);
2129 li->blbnb.resize(li->ncc);
2131 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2132 for (th = 0; th < li->ntask; th++)
2138 li_task = &li->task[th];
2140 if (li->ncg_triangle > 0)
2142 /* This is allocating too much, but it is difficult to improve */
2143 li_task->triangle.resize(li_task->b1 - li_task->b0);
2144 li_task->tri_bits.resize(li_task->b1 - li_task->b0);
2147 set_matrix_indices(li, li_task, &at2con, bSortMatrix);
2149 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2152 done_blocka(&at2con);
2154 if (cr->dd == nullptr)
2156 /* Since the matrix is static, we should free some memory */
2157 li->blbnb.resize(li->ncc);
2160 li->blmf.resize(li->ncc);
2161 li->blmf1.resize(li->ncc);
2162 li->tmpncc.resize(li->ncc);
2164 gmx::ArrayRef<const int> nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2165 if (!nlocat_dd.empty())
2167 /* Convert nlocat from local topology to LINCS constraint indexing */
2168 for (con = 0; con < ncon_tot; con++)
2170 li->nlocat[li->con_index[con]] = nlocat_dd[con];
2180 fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n",
2181 li->nc_real, li->nc, li->ncc);
2186 lincs_thread_setup(li, md.nr);
2189 set_lincs_matrix(li, md.invmass, md.lambda);
2192 //! Issues a warning when LINCS constraints cannot be satisfied.
2193 static void lincs_warning(gmx_domdec_t *dd, const rvec *x, rvec *xprime, t_pbc *pbc,
2194 int ncons, gmx::ArrayRef<const int> bla,
2195 gmx::ArrayRef<const real> bllen, real wangle,
2196 int maxwarn, int *warncount)
2200 real wfac, d0, d1, cosine;
2202 wfac = std::cos(DEG2RAD*wangle);
2205 "bonds that rotated more than %g degrees:\n"
2206 " atom 1 atom 2 angle previous, current, constraint length\n",
2209 for (b = 0; b < ncons; b++)
2215 pbc_dx_aiuc(pbc, x[i], x[j], v0);
2216 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2220 rvec_sub(x[i], x[j], v0);
2221 rvec_sub(xprime[i], xprime[j], v1);
2225 cosine = ::iprod(v0, v1)/(d0*d1);
2229 " %6d %6d %5.1f %8.4f %8.4f %8.4f\n",
2230 ddglatnr(dd, i), ddglatnr(dd, j),
2231 RAD2DEG*std::acos(cosine), d0, d1, bllen[b]);
2232 if (!std::isfinite(d1))
2234 gmx_fatal(FARGS, "Bond length not finite.");
2240 if (*warncount > maxwarn)
2242 too_many_constraint_warnings(econtLINCS, *warncount);
2246 //! Determine how well the constraints have been satisfied.
2247 static void cconerr(const Lincs &lincsd,
2248 const rvec *x, const t_pbc *pbc,
2249 real *ncons_loc, real *ssd, real *max, int *imax)
2251 gmx::ArrayRef<const int> bla = lincsd.bla;
2252 gmx::ArrayRef<const real> bllen = lincsd.bllen;
2253 gmx::ArrayRef<const int> nlocat = lincsd.nlocat;
2259 for (int task = 0; task < lincsd.ntask; task++)
2261 for (int b = lincsd.task[task].b0; b < lincsd.task[task].b1; b++)
2266 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
2270 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
2272 real r2 = ::norm2(dx);
2273 real len = r2*gmx::invsqrt(r2);
2274 real d = std::abs(len/bllen[b]-1);
2275 if (d > ma && (nlocat.empty() || nlocat[b]))
2287 ssd2 += nlocat[b]*d*d;
2293 *ncons_loc = (nlocat.empty() ? 1 : 0.5)*count;
2294 *ssd = (nlocat.empty() ? 1 : 0.5)*ssd2;
2299 bool constrain_lincs(bool computeRmsd,
2300 const t_inputrec &ir,
2302 Lincs *lincsd, const t_mdatoms &md,
2303 const t_commrec *cr,
2304 const gmx_multisim_t *ms,
2305 const rvec *x, rvec *xprime, rvec *min_proj,
2306 matrix box, t_pbc *pbc,
2307 real lambda, real *dvdlambda,
2308 real invdt, rvec *v,
2309 bool bCalcVir, tensor vir_r_m_dr,
2310 ConstraintVariable econq,
2312 int maxwarn, int *warncount)
2315 char buf2[22], buf3[STRLEN];
2317 real ncons_loc, p_ssd, p_max = 0;
2323 /* This boolean should be set by a flag passed to this routine.
2324 * We can also easily check if any constraint length is changed,
2325 * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2327 bCalcDHDL = (ir.efep != efepNO && dvdlambda != nullptr);
2329 if (lincsd->nc == 0 && cr->dd == nullptr)
2333 lincsd->rmsdData = {{0}};
2339 if (econq == ConstraintVariable::Positions)
2341 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2342 * also with efep!=fepNO.
2344 if (ir.efep != efepNO)
2346 if (md.nMassPerturbed && lincsd->matlam != md.lambda)
2348 set_lincs_matrix(lincsd, md.invmass, md.lambda);
2351 for (i = 0; i < lincsd->nc; i++)
2353 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
2357 if (lincsd->ncg_flex)
2359 /* Set the flexible constraint lengths to the old lengths */
2362 for (i = 0; i < lincsd->nc; i++)
2364 if (lincsd->bllen[i] == 0)
2366 pbc_dx_aiuc(pbc, x[lincsd->bla[2*i]], x[lincsd->bla[2*i+1]], dx);
2367 lincsd->bllen[i] = norm(dx);
2373 for (i = 0; i < lincsd->nc; i++)
2375 if (lincsd->bllen[i] == 0)
2378 std::sqrt(distance2(x[lincsd->bla[2*i]],
2379 x[lincsd->bla[2*i+1]]));
2387 cconerr(*lincsd, xprime, pbc,
2388 &ncons_loc, &p_ssd, &p_max, &p_imax);
2391 /* This bWarn var can be updated by multiple threads
2392 * at the same time. But as we only need to detect
2393 * if a warning occurred or not, this is not an issue.
2397 /* The OpenMP parallel region of constrain_lincs for coords */
2398 #pragma omp parallel num_threads(lincsd->ntask)
2402 int th = gmx_omp_get_thread_num();
2404 clear_mat(lincsd->task[th].vir_r_m_dr);
2406 do_lincs(x, xprime, box, pbc, lincsd, th,
2409 ir.LincsWarnAngle, &bWarn,
2411 th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2413 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2416 if (debug && lincsd->nc > 0)
2418 fprintf(debug, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2419 fprintf(debug, " Before LINCS %.6f %.6f %6d %6d\n",
2420 std::sqrt(p_ssd/ncons_loc), p_max,
2421 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2422 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2424 if (computeRmsd || debug)
2426 cconerr(*lincsd, xprime, pbc,
2427 &ncons_loc, &p_ssd, &p_max, &p_imax);
2428 lincsd->rmsdData[0] = ncons_loc;
2429 lincsd->rmsdData[1] = p_ssd;
2433 lincsd->rmsdData = {{0}};
2435 if (debug && lincsd->nc > 0)
2438 " After LINCS %.6f %.6f %6d %6d\n\n",
2439 std::sqrt(p_ssd/ncons_loc), p_max,
2440 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2441 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2446 if (maxwarn < INT_MAX)
2448 cconerr(*lincsd, xprime, pbc,
2449 &ncons_loc, &p_ssd, &p_max, &p_imax);
2452 sprintf(buf3, " in simulation %d", ms->sim);
2459 "\nStep %s, time %g (ps) LINCS WARNING%s\n"
2460 "relative constraint deviation after LINCS:\n"
2461 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2462 gmx_step_str(step, buf2), ir.init_t+step*ir.delta_t,
2464 std::sqrt(p_ssd/ncons_loc), p_max,
2465 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2466 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2468 lincs_warning(cr->dd, x, xprime, pbc,
2469 lincsd->nc, lincsd->bla, lincsd->bllen,
2470 ir.LincsWarnAngle, maxwarn, warncount);
2472 bOK = (p_max < 0.5);
2475 if (lincsd->ncg_flex)
2477 for (i = 0; (i < lincsd->nc); i++)
2479 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2481 lincsd->bllen[i] = 0;
2488 /* The OpenMP parallel region of constrain_lincs for derivatives */
2489 #pragma omp parallel num_threads(lincsd->ntask)
2493 int th = gmx_omp_get_thread_num();
2495 do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
2496 md.invmass, econq, bCalcDHDL,
2497 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2499 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2505 /* Reduce the dH/dlambda contributions over the threads */
2510 for (th = 0; th < lincsd->ntask; th++)
2512 dhdlambda += lincsd->task[th].dhdlambda;
2514 if (econq == ConstraintVariable::Positions)
2516 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2517 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2518 dhdlambda /= ir.delta_t*ir.delta_t;
2520 *dvdlambda += dhdlambda;
2523 if (bCalcVir && lincsd->ntask > 1)
2525 for (i = 1; i < lincsd->ntask; i++)
2527 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2531 /* count assuming nit=1 */
2532 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2533 inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
2534 if (lincsd->ntriangle > 0)
2536 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
2540 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real*2);
2544 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);