2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief Defines LINCS code.
40 * \author Berk Hess <hess@kth.se>
41 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \ingroup module_mdlib
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/gmxlib/nrnb.h"
60 #include "gromacs/math/functions.h"
61 #include "gromacs/math/paddedvector.h"
62 #include "gromacs/math/units.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/mdlib/constr.h"
65 #include "gromacs/mdlib/gmx_omp_nthreads.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/mdatom.h"
70 #include "gromacs/pbcutil/pbc.h"
71 #include "gromacs/pbcutil/pbc_simd.h"
72 #include "gromacs/simd/simd.h"
73 #include "gromacs/simd/simd_math.h"
74 #include "gromacs/simd/vector_operations.h"
75 #include "gromacs/topology/block.h"
76 #include "gromacs/topology/mtop_util.h"
77 #include "gromacs/utility/alignedallocator.h"
78 #include "gromacs/utility/arrayref.h"
79 #include "gromacs/utility/basedefinitions.h"
80 #include "gromacs/utility/bitmask.h"
81 #include "gromacs/utility/cstringutil.h"
82 #include "gromacs/utility/exceptions.h"
83 #include "gromacs/utility/fatalerror.h"
84 #include "gromacs/utility/gmxomp.h"
85 #include "gromacs/utility/pleasecite.h"
87 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
92 //! Indices of the two atoms involved in a single constraint
95 //! \brief Constructor, does not initialize to catch bugs and faster construction
106 //! Unit of work within LINCS.
109 //! First constraint for this task.
111 //! b1-1 is the last constraint for this task.
113 //! The number of constraints in triangles.
115 //! The list of triangle constraints.
116 std::vector<int> triangle;
117 //! The bits tell if the matrix element should be used.
118 std::vector<int> tri_bits;
119 //! Constraint index for updating atom data.
120 std::vector<int> ind;
121 //! Constraint index for updating atom data.
122 std::vector<int> ind_r;
123 //! Temporary variable for virial calculation.
124 tensor vir_r_m_dr = {{0}};
125 //! Temporary variable for lambda derivative.
129 /*! \brief Data for LINCS algorithm.
134 //! The global number of constraints.
136 //! The global number of flexible constraints.
138 //! The global number of constraints in triangles.
139 int ncg_triangle = 0;
140 //! The number of iterations.
142 //! The order of the matrix expansion.
144 //! The maximum number of constraints connected to a single atom.
147 //! The number of real constraints.
149 //! The number of constraints including padding for SIMD.
151 //! The number of constraint connections.
153 //! The FE lambda value used for filling blc and blmf.
155 //! mapping from topology to LINCS constraints.
156 std::vector<int> con_index;
157 //! The reference distance in topology A.
158 std::vector < real, AlignedAllocator < real>> bllen0;
159 //! The reference distance in top B - the r.d. in top A.
160 std::vector < real, AlignedAllocator < real>> ddist;
161 //! The atom pairs involved in the constraints.
162 std::vector<AtomPair> atoms;
163 //! 1/sqrt(invmass1 invmass2).
164 std::vector < real, AlignedAllocator < real>> blc;
165 //! As blc, but with all masses 1.
166 std::vector < real, AlignedAllocator < real>> blc1;
167 //! Index into blbnb and blmf.
168 std::vector<int> blnr;
169 //! List of constraint connections.
170 std::vector<int> blbnb;
171 //! The local number of constraints in triangles.
173 //! The number of constraint connections in triangles.
174 int ncc_triangle = 0;
175 //! Communicate before each LINCS interation.
176 bool bCommIter = false;
177 //! Matrix of mass factors for constraint connections.
178 std::vector<real> blmf;
179 //! As blmf, but with all masses 1.
180 std::vector<real> blmf1;
181 //! The reference bond length.
182 std::vector < real, AlignedAllocator < real>> bllen;
183 //! The local atom count per constraint, can be NULL.
184 std::vector<int> nlocat;
186 /*! \brief The number of tasks used for LINCS work.
188 * \todo This is mostly used to loop over \c task, which would
189 * be nicer to do with range-based for loops, but the thread
190 * index is used for constructing bit masks and organizing the
191 * virial output buffer, so other things need to change,
194 /*! \brief LINCS thread division */
195 std::vector<Task> task;
196 //! Atom flags for thread parallelization.
197 std::vector<gmx_bitmask_t> atf;
198 //! Are the LINCS tasks interdependent?
199 bool bTaskDep = false;
200 //! Are there triangle constraints that cross task borders?
201 bool bTaskDepTri = false;
202 //! Arrays for temporary storage in the LINCS algorithm.
204 PaddedVector<gmx::RVec> tmpv;
205 std::vector<real> tmpncc;
206 std::vector < real, AlignedAllocator < real>> tmp1;
207 std::vector < real, AlignedAllocator < real>> tmp2;
208 std::vector < real, AlignedAllocator < real>> tmp3;
209 std::vector < real, AlignedAllocator < real>> tmp4;
211 //! The Lagrange multipliers times -1.
212 std::vector < real, AlignedAllocator < real>> mlambda;
213 //! Storage for the constraint RMS relative deviation output.
214 std::array<real, 2> rmsdData = {{0}};
217 /*! \brief Define simd_width for memory allocation used for SIMD code */
218 #if GMX_SIMD_HAVE_REAL
219 static const int simd_width = GMX_SIMD_REAL_WIDTH;
221 static const int simd_width = 1;
224 ArrayRef<real> lincs_rmsdData(Lincs *lincsd)
226 return lincsd->rmsdData;
229 real lincs_rmsd(const Lincs *lincsd)
231 if (lincsd->rmsdData[0] > 0)
233 return std::sqrt(lincsd->rmsdData[1]/lincsd->rmsdData[0]);
241 /*! \brief Do a set of nrec LINCS matrix multiplications.
243 * This function will return with up to date thread-local
244 * constraint data, without an OpenMP barrier.
246 static void lincs_matrix_expand(const Lincs &lincsd,
248 gmx::ArrayRef<const real> blcc,
249 gmx::ArrayRef<real> rhs1,
250 gmx::ArrayRef<real> rhs2,
251 gmx::ArrayRef<real> sol)
253 gmx::ArrayRef<const int> blnr = lincsd.blnr;
254 gmx::ArrayRef<const int> blbnb = lincsd.blbnb;
256 const int b0 = li_task.b0;
257 const int b1 = li_task.b1;
258 const int nrec = lincsd.nOrder;
260 for (int rec = 0; rec < nrec; rec++)
266 for (int b = b0; b < b1; b++)
272 for (n = blnr[b]; n < blnr[b+1]; n++)
274 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
277 sol[b] = sol[b] + mvb;
280 std::swap(rhs1, rhs2);
281 } /* nrec*(ncons+2*nrtot) flops */
283 if (lincsd.ntriangle > 0)
285 /* Perform an extra nrec recursions for only the constraints
286 * involved in rigid triangles.
287 * In this way their accuracy should come close to those of the other
288 * constraints, since traingles of constraints can produce eigenvalues
289 * around 0.7, while the effective eigenvalue for bond constraints
290 * is around 0.4 (and 0.7*0.7=0.5).
295 /* We need a barrier here, since other threads might still be
296 * reading the contents of rhs1 and/o rhs2.
297 * We could avoid this barrier by introducing two extra rhs
298 * arrays for the triangle constraints only.
303 /* Constraints involved in a triangle are ensured to be in the same
304 * LINCS task. This means no barriers are required during the extra
305 * iterations for the triangle constraints.
307 gmx::ArrayRef<const int> triangle = li_task.triangle;
308 gmx::ArrayRef<const int> tri_bits = li_task.tri_bits;
310 for (int rec = 0; rec < nrec; rec++)
312 for (int tb = 0; tb < li_task.ntriangle; tb++)
314 int b, bits, nr0, nr1, n;
322 for (n = nr0; n < nr1; n++)
324 if (bits & (1 << (n - nr0)))
326 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
330 sol[b] = sol[b] + mvb;
333 std::swap(rhs1, rhs2);
334 } /* nrec*(ntriangle + ncc_triangle*2) flops */
336 if (lincsd.bTaskDepTri)
338 /* The constraints triangles are decoupled from each other,
339 * but constraints in one triangle cross thread task borders.
340 * We could probably avoid this with more advanced setup code.
347 //! Update atomic coordinates when an index is not required.
349 lincs_update_atoms_noind(int ncons,
350 gmx::ArrayRef<const AtomPair> atoms,
352 gmx::ArrayRef<const real> fac,
353 gmx::ArrayRef<const gmx::RVec> r,
357 if (invmass != nullptr)
359 for (int b = 0; b < ncons; b++)
361 int i = atoms[b].index1;
362 int j = atoms[b].index2;
363 real mvb = preFactor*fac[b];
364 real im1 = invmass[i];
365 real im2 = invmass[j];
366 real tmp0 = r[b][0]*mvb;
367 real tmp1 = r[b][1]*mvb;
368 real tmp2 = r[b][2]*mvb;
375 } /* 16 ncons flops */
379 for (int b = 0; b < ncons; b++)
381 int i = atoms[b].index1;
382 int j = atoms[b].index2;
383 real mvb = preFactor*fac[b];
384 real tmp0 = r[b][0]*mvb;
385 real tmp1 = r[b][1]*mvb;
386 real tmp2 = r[b][2]*mvb;
397 //! Update atomic coordinates when an index is required.
399 lincs_update_atoms_ind(gmx::ArrayRef<const int> ind,
400 gmx::ArrayRef<const AtomPair> atoms,
402 gmx::ArrayRef<const real> fac,
403 gmx::ArrayRef<const gmx::RVec> r,
407 if (invmass != nullptr)
411 int i = atoms[b].index1;
412 int j = atoms[b].index2;
413 real mvb = preFactor*fac[b];
414 real im1 = invmass[i];
415 real im2 = invmass[j];
416 real tmp0 = r[b][0]*mvb;
417 real tmp1 = r[b][1]*mvb;
418 real tmp2 = r[b][2]*mvb;
425 } /* 16 ncons flops */
431 int i = atoms[b].index1;
432 int j = atoms[b].index2;
433 real mvb = preFactor*fac[b];
434 real tmp0 = r[b][0]*mvb;
435 real tmp1 = r[b][1]*mvb;
436 real tmp2 = r[b][2]*mvb;
443 } /* 16 ncons flops */
447 //! Update coordinates for atoms.
449 lincs_update_atoms(Lincs *li,
452 gmx::ArrayRef<const real> fac,
453 gmx::ArrayRef<const gmx::RVec> r,
459 /* Single thread, we simply update for all constraints */
460 lincs_update_atoms_noind(li->nc_real,
461 li->atoms, preFactor, fac, r, invmass, x);
465 /* Update the atom vector components for our thread local
466 * constraints that only access our local atom range.
467 * This can be done without a barrier.
469 lincs_update_atoms_ind(li->task[th].ind,
470 li->atoms, preFactor, fac, r, invmass, x);
472 if (!li->task[li->ntask].ind.empty())
474 /* Update the constraints that operate on atoms
475 * in multiple thread atom blocks on the master thread.
480 lincs_update_atoms_ind(li->task[li->ntask].ind,
481 li->atoms, preFactor, fac, r, invmass, x);
487 #if GMX_SIMD_HAVE_REAL
488 //! Helper function so that we can run TSAN with SIMD support (where implemented).
490 static inline void gmx_simdcall
491 gatherLoadUTransposeTSANSafe(const real *base,
492 const std::int32_t *offset,
497 #if (CMAKE_BUILD_TYPE == CMAKE_BUILD_TYPE_TSAN) && GMX_SIMD_X86_AVX2_256
498 // This function is only implemented in this case
499 gatherLoadUTransposeSafe<align>(base, offset, v0, v1, v2);
501 gatherLoadUTranspose<align>(base, offset, v0, v1, v2);
505 /*! \brief Calculate the constraint distance vectors r to project on from x.
507 * Determine the right-hand side of the matrix equation using quantity f.
508 * This function only differs from calc_dr_x_xp_simd below in that
509 * no constraint length is subtracted and no PBC is used for f. */
510 static void gmx_simdcall
511 calc_dr_x_f_simd(int b0,
513 gmx::ArrayRef<const AtomPair> atoms,
514 const rvec * gmx_restrict x,
515 const rvec * gmx_restrict f,
516 const real * gmx_restrict blc,
517 const real * pbc_simd,
518 rvec * gmx_restrict r,
519 real * gmx_restrict rhs,
520 real * gmx_restrict sol)
522 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
524 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
526 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
531 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
533 SimdReal x0_S, y0_S, z0_S;
534 SimdReal x1_S, y1_S, z1_S;
535 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
536 SimdReal fx_S, fy_S, fz_S, ip_S, rhs_S;
537 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
538 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
540 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
542 offset0[i] = atoms[bs + i].index1;
543 offset1[i] = atoms[bs + i].index2;
546 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
547 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
552 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
554 n2_S = norm2(rx_S, ry_S, rz_S);
555 il_S = invsqrt(n2_S);
561 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
563 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(f), offset0, &x0_S, &y0_S, &z0_S);
564 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(f), offset1, &x1_S, &y1_S, &z1_S);
569 ip_S = iprod(rx_S, ry_S, rz_S, fx_S, fy_S, fz_S);
571 rhs_S = load<SimdReal>(blc + bs) * ip_S;
573 store(rhs + bs, rhs_S);
574 store(sol + bs, rhs_S);
577 #endif // GMX_SIMD_HAVE_REAL
579 /*! \brief LINCS projection, works on derivatives of the coordinates. */
580 static void do_lincsp(const rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
581 Lincs *lincsd, int th,
583 ConstraintVariable econq, bool bCalcDHDL,
584 bool bCalcVir, tensor rmdf)
586 const int b0 = lincsd->task[th].b0;
587 const int b1 = lincsd->task[th].b1;
589 gmx::ArrayRef<const AtomPair> atoms = lincsd->atoms;
590 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
591 gmx::ArrayRef<const int> blnr = lincsd->blnr;
592 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
594 gmx::ArrayRef<const real> blc;
595 gmx::ArrayRef<const real> blmf;
596 if (econq != ConstraintVariable::Force)
598 /* Use mass-weighted parameters */
604 /* Use non mass-weighted parameters */
606 blmf = lincsd->blmf1;
608 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
609 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
610 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
611 gmx::ArrayRef<real> sol = lincsd->tmp3;
613 #if GMX_SIMD_HAVE_REAL
614 /* This SIMD code does the same as the plain-C code after the #else.
615 * The only difference is that we always call pbc code, as with SIMD
616 * the overhead of pbc computation (when not needed) is small.
618 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
620 /* Convert the pbc struct for SIMD */
621 set_pbc_simd(pbc, pbc_simd);
623 /* Compute normalized x i-j vectors, store in r.
624 * Compute the inner product of r and xp i-j and store in rhs1.
626 calc_dr_x_f_simd(b0, b1, atoms, x, f, blc.data(),
628 as_rvec_array(r.data()), rhs1.data(), sol.data());
630 #else // GMX_SIMD_HAVE_REAL
632 /* Compute normalized i-j vectors */
635 for (int b = b0; b < b1; b++)
639 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
645 for (int b = b0; b < b1; b++)
649 rvec_sub(x[atoms[b].index1], x[atoms[b].index2], dx);
651 } /* 16 ncons flops */
654 for (int b = b0; b < b1; b++)
656 int i = atoms[b].index1;
657 int j = atoms[b].index2;
658 real mvb = blc[b]*(r[b][0]*(f[i][0] - f[j][0]) +
659 r[b][1]*(f[i][1] - f[j][1]) +
660 r[b][2]*(f[i][2] - f[j][2]));
666 #endif // GMX_SIMD_HAVE_REAL
668 if (lincsd->bTaskDep)
670 /* We need a barrier, since the matrix construction below
671 * can access entries in r of other threads.
676 /* Construct the (sparse) LINCS matrix */
677 for (int b = b0; b < b1; b++)
679 for (int n = blnr[b]; n < blnr[b+1]; n++)
681 blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
684 /* Together: 23*ncons + 6*nrtot flops */
686 lincs_matrix_expand(*lincsd, lincsd->task[th], blcc, rhs1, rhs2, sol);
687 /* nrec*(ncons+2*nrtot) flops */
689 if (econq == ConstraintVariable::Deriv_FlexCon)
691 /* We only want to constraint the flexible constraints,
692 * so we mask out the normal ones by setting sol to 0.
694 for (int b = b0; b < b1; b++)
696 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
703 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
704 for (int b = b0; b < b1; b++)
709 /* When constraining forces, we should not use mass weighting,
710 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
712 lincs_update_atoms(lincsd, th, 1.0, sol, r,
713 (econq != ConstraintVariable::Force) ? invmass : nullptr, fp);
718 for (int b = b0; b < b1; b++)
720 dhdlambda -= sol[b]*lincsd->ddist[b];
723 lincsd->task[th].dhdlambda = dhdlambda;
728 /* Constraint virial,
729 * determines sum r_bond x delta f,
730 * where delta f is the constraint correction
731 * of the quantity that is being constrained.
733 for (int b = b0; b < b1; b++)
735 const real mvb = lincsd->bllen[b]*sol[b];
736 for (int i = 0; i < DIM; i++)
738 const real tmp1 = mvb*r[b][i];
739 for (int j = 0; j < DIM; j++)
741 rmdf[i][j] += tmp1*r[b][j];
744 } /* 23 ncons flops */
748 #if GMX_SIMD_HAVE_REAL
750 /*! \brief Calculate the constraint distance vectors r to project on from x.
752 * Determine the right-hand side of the matrix equation using coordinates xp. */
753 static void gmx_simdcall
754 calc_dr_x_xp_simd(int b0,
756 gmx::ArrayRef<const AtomPair> atoms,
757 const rvec * gmx_restrict x,
758 const rvec * gmx_restrict xp,
759 const real * gmx_restrict bllen,
760 const real * gmx_restrict blc,
761 const real * pbc_simd,
762 rvec * gmx_restrict r,
763 real * gmx_restrict rhs,
764 real * gmx_restrict sol)
766 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
767 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
769 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
774 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
776 SimdReal x0_S, y0_S, z0_S;
777 SimdReal x1_S, y1_S, z1_S;
778 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
779 SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
780 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
781 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
783 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
785 offset0[i] = atoms[bs + i].index1;
786 offset1[i] = atoms[bs + i].index2;
789 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
790 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
795 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
797 n2_S = norm2(rx_S, ry_S, rz_S);
798 il_S = invsqrt(n2_S);
804 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
806 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(xp), offset0, &x0_S, &y0_S, &z0_S);
807 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(xp), offset1, &x1_S, &y1_S, &z1_S);
813 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
815 ip_S = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
817 rhs_S = load<SimdReal>(blc + bs) * (ip_S - load<SimdReal>(bllen + bs));
819 store(rhs + bs, rhs_S);
820 store(sol + bs, rhs_S);
823 #endif // GMX_SIMD_HAVE_REAL
825 /*! \brief Determine the distances and right-hand side for the next iteration. */
826 gmx_unused static void calc_dist_iter(
829 gmx::ArrayRef<const AtomPair> atoms,
830 const rvec * gmx_restrict xp,
831 const real * gmx_restrict bllen,
832 const real * gmx_restrict blc,
835 real * gmx_restrict rhs,
836 real * gmx_restrict sol,
839 for (int b = b0; b < b1; b++)
845 pbc_dx_aiuc(pbc, xp[atoms[b].index1], xp[atoms[b].index2], dx);
849 rvec_sub(xp[atoms[b].index1], xp[atoms[b].index2], dx);
852 real dlen2 = 2*len2 - ::norm2(dx);
853 if (dlen2 < wfac*len2)
855 /* not race free - see detailed comment in caller */
861 mvb = blc[b]*(len - dlen2*gmx::invsqrt(dlen2));
869 } /* 20*ncons flops */
872 #if GMX_SIMD_HAVE_REAL
873 /*! \brief As calc_dist_iter(), but using SIMD intrinsics. */
874 static void gmx_simdcall
875 calc_dist_iter_simd(int b0,
877 gmx::ArrayRef<const AtomPair> atoms,
878 const rvec * gmx_restrict x,
879 const real * gmx_restrict bllen,
880 const real * gmx_restrict blc,
881 const real * pbc_simd,
883 real * gmx_restrict rhs,
884 real * gmx_restrict sol,
887 SimdReal min_S(GMX_REAL_MIN);
889 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 (int 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] = atoms[bs + i].index1;
909 offset1[i] = atoms[bs + i].index2;
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 AtomPair> atoms = lincsd->atoms;
967 gmx::ArrayRef<gmx::RVec> r = lincsd->tmpv;
968 gmx::ArrayRef<const int> blnr = lincsd->blnr;
969 gmx::ArrayRef<const int> blbnb = lincsd->blbnb;
970 gmx::ArrayRef<const real> blc = lincsd->blc;
971 gmx::ArrayRef<const real> blmf = lincsd->blmf;
972 gmx::ArrayRef<const real> bllen = lincsd->bllen;
973 gmx::ArrayRef<real> blcc = lincsd->tmpncc;
974 gmx::ArrayRef<real> rhs1 = lincsd->tmp1;
975 gmx::ArrayRef<real> rhs2 = lincsd->tmp2;
976 gmx::ArrayRef<real> sol = lincsd->tmp3;
977 gmx::ArrayRef<real> blc_sol = lincsd->tmp4;
978 gmx::ArrayRef<real> mlambda = lincsd->mlambda;
979 gmx::ArrayRef<const int> nlocat = lincsd->nlocat;
981 #if GMX_SIMD_HAVE_REAL
983 /* This SIMD code does the same as the plain-C code after the #else.
984 * The only difference is that we always call pbc code, as with SIMD
985 * the overhead of pbc computation (when not needed) is small.
987 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
989 /* Convert the pbc struct for SIMD */
990 set_pbc_simd(pbc, pbc_simd);
992 /* Compute normalized x i-j vectors, store in r.
993 * Compute the inner product of r and xp i-j and store in rhs1.
995 calc_dr_x_xp_simd(b0, b1, atoms, x, xp, bllen.data(), blc.data(),
997 as_rvec_array(r.data()), rhs1.data(), sol.data());
999 #else // GMX_SIMD_HAVE_REAL
1003 /* Compute normalized i-j vectors */
1004 for (int b = b0; b < b1; b++)
1007 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
1010 pbc_dx_aiuc(pbc, xp[atoms[b].index1], xp[atoms[b].index2], dx);
1011 real mvb = blc[b]*(::iprod(r[b], dx) - bllen[b]);
1018 /* Compute normalized i-j vectors */
1019 for (int b = b0; b < b1; b++)
1021 int i = atoms[b].index1;
1022 int j = atoms[b].index2;
1023 real tmp0 = x[i][0] - x[j][0];
1024 real tmp1 = x[i][1] - x[j][1];
1025 real tmp2 = x[i][2] - x[j][2];
1026 real rlen = gmx::invsqrt(tmp0*tmp0 + tmp1*tmp1 + tmp2*tmp2);
1027 r[b][0] = rlen*tmp0;
1028 r[b][1] = rlen*tmp1;
1029 r[b][2] = rlen*tmp2;
1030 /* 16 ncons flops */
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, atoms,
1113 xp, bllen.data(), blc.data(), pbc_simd, wfac,
1114 rhs1.data(), sol.data(), bWarn);
1116 calc_dist_iter(b0, b1, atoms, 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)
1218 /* Construct the coupling coefficient matrix blmf */
1219 li_task->ntriangle = 0;
1221 *nCrossTaskTriangles = 0;
1222 for (int i = li_task->b0; i < li_task->b1; i++)
1224 const int a1 = li->atoms[i].index1;
1225 const int a2 = li->atoms[i].index2;
1226 for (int n = li->blnr[i]; n < li->blnr[i + 1]; n++)
1228 const int k = li->blbnb[n];
1230 /* If we are using multiple, independent tasks for LINCS,
1231 * the calls to check_assign_connected should have
1232 * put all connected constraints in our task.
1234 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1237 if (a1 == li->atoms[k].index1 || a2 == li->atoms[k].index2)
1247 if (a1 == li->atoms[k].index1 || a1 == li->atoms[k].index2)
1257 li->blmf[n] = sign*invmass[center]*li->blc[i]*li->blc[k];
1258 li->blmf1[n] = sign*0.5;
1259 if (li->ncg_triangle > 0)
1261 /* Look for constraint triangles */
1262 for (int nk = li->blnr[k]; nk < li->blnr[k + 1]; nk++)
1264 const int kk = li->blbnb[nk];
1265 if (kk != i && kk != k &&
1266 (li->atoms[kk].index1 == end ||
1267 li->atoms[kk].index2 == end))
1269 /* Check if the constraints in this triangle actually
1270 * belong to a different task. We still assign them
1271 * here, since it's convenient for the triangle
1272 * iterations, but we then need an extra barrier.
1274 if (k < li_task->b0 || k >= li_task->b1 ||
1275 kk < li_task->b0 || kk >= li_task->b1)
1277 (*nCrossTaskTriangles)++;
1280 if (li_task->ntriangle == 0 ||
1281 li_task->triangle[li_task->ntriangle - 1] < i)
1283 /* Add this constraint to the triangle list */
1284 li_task->triangle[li_task->ntriangle] = i;
1285 li_task->tri_bits[li_task->ntriangle] = 0;
1286 li_task->ntriangle++;
1287 if (li->blnr[i+1] - li->blnr[i] > static_cast<int>(sizeof(li_task->tri_bits[0])*8 - 1))
1289 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %zu allowed for constraints participating in triangles",
1290 li->blnr[i+1] - li->blnr[i],
1291 sizeof(li_task->tri_bits[0])*8-1);
1294 li_task->tri_bits[li_task->ntriangle-1] |= (1 << (n - li->blnr[i]));
1303 /*! \brief Sets the elements in the LINCS matrix. */
1304 static void set_lincs_matrix(Lincs *li, real *invmass, real lambda)
1306 const real invsqrt2 = 0.7071067811865475244;
1308 for (int i = 0; (i < li->nc); i++)
1310 const int a1 = li->atoms[i].index1;
1311 const int a2 = li->atoms[i].index2;
1312 li->blc[i] = gmx::invsqrt(invmass[a1] + invmass[a2]);
1313 li->blc1[i] = invsqrt2;
1316 /* Construct the coupling coefficient matrix blmf */
1317 int ntriangle = 0, ncc_triangle = 0, nCrossTaskTriangles = 0;
1318 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle, nCrossTaskTriangles) num_threads(li->ntask) schedule(static)
1319 for (int th = 0; th < li->ntask; th++)
1323 set_lincs_matrix_task(li, &li->task[th], invmass,
1324 &ncc_triangle, &nCrossTaskTriangles);
1325 ntriangle += li->task[th].ntriangle;
1327 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1329 li->ntriangle = ntriangle;
1330 li->ncc_triangle = ncc_triangle;
1331 li->bTaskDepTri = (nCrossTaskTriangles > 0);
1335 fprintf(debug, "The %d constraints participate in %d triangles\n",
1336 li->nc, li->ntriangle);
1337 fprintf(debug, "There are %d constraint couplings, of which %d in triangles\n",
1338 li->ncc, li->ncc_triangle);
1339 if (li->ntriangle > 0 && li->ntask > 1)
1341 fprintf(debug, "%d constraint triangles contain constraints assigned to different tasks\n",
1342 nCrossTaskTriangles);
1347 * so we know with which lambda value the masses have been set.
1349 li->matlam = lambda;
1352 //! Finds all triangles of atoms that share constraints to a central atom.
1353 static int count_triangle_constraints(const InteractionLists &ilist,
1354 const t_blocka &at2con)
1356 int ncon1, ncon_tot;
1357 int c0, n1, c1, ac1, n2, c2;
1360 ncon1 = ilist[F_CONSTR].size()/3;
1361 ncon_tot = ncon1 + ilist[F_CONSTRNC].size()/3;
1363 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1364 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1367 for (c0 = 0; c0 < ncon_tot; c0++)
1369 bool bTriangle = FALSE;
1370 const int *iap = constr_iatomptr(ia1, ia2, c0);
1371 const int a00 = iap[1];
1372 const int a01 = iap[2];
1373 for (n1 = at2con.index[a01]; n1 < at2con.index[a01+1]; n1++)
1378 const int *iap = constr_iatomptr(ia1, ia2, c1);
1379 const int a10 = iap[1];
1380 const int a11 = iap[2];
1389 for (n2 = at2con.index[ac1]; n2 < at2con.index[ac1+1]; n2++)
1392 if (c2 != c0 && c2 != c1)
1394 const int *iap = constr_iatomptr(ia1, ia2, c2);
1395 const int a20 = iap[1];
1396 const int a21 = iap[2];
1397 if (a20 == a00 || a21 == a00)
1411 return ncon_triangle;
1414 //! Finds sequences of sequential constraints.
1415 static bool more_than_two_sequential_constraints(const InteractionLists &ilist,
1416 const t_blocka &at2con)
1418 int ncon1, ncon_tot, c;
1419 bool bMoreThanTwoSequentialConstraints;
1421 ncon1 = ilist[F_CONSTR].size()/3;
1422 ncon_tot = ncon1 + ilist[F_CONSTRNC].size()/3;
1424 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1425 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1427 bMoreThanTwoSequentialConstraints = FALSE;
1428 for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
1430 const int *iap = constr_iatomptr(ia1, ia2, c);
1431 const int a1 = iap[1];
1432 const int a2 = iap[2];
1433 /* Check if this constraint has constraints connected at both atoms */
1434 if (at2con.index[a1+1] - at2con.index[a1] > 1 &&
1435 at2con.index[a2+1] - at2con.index[a2] > 1)
1437 bMoreThanTwoSequentialConstraints = TRUE;
1441 return bMoreThanTwoSequentialConstraints;
1444 Lincs *init_lincs(FILE *fplog, const gmx_mtop_t &mtop,
1445 int nflexcon_global, ArrayRef<const t_blocka> at2con,
1446 bool bPLINCS, int nIter, int nProjOrder)
1448 // TODO this should become a unique_ptr
1450 bool bMoreThanTwoSeq;
1454 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
1455 bPLINCS ? " Parallel" : "");
1461 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1462 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1463 li->ncg_flex = nflexcon_global;
1466 li->nOrder = nProjOrder;
1468 li->max_connect = 0;
1469 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1471 for (int a = 0; a < mtop.moltype[mt].atoms.nr; a++)
1473 li->max_connect = std::max(li->max_connect,
1474 at2con[mt].index[a + 1] - at2con[mt].index[a]);
1478 li->ncg_triangle = 0;
1479 bMoreThanTwoSeq = FALSE;
1480 for (const gmx_molblock_t &molb : mtop.molblock)
1482 const gmx_moltype_t &molt = mtop.moltype[molb.type];
1486 count_triangle_constraints(molt.ilist, at2con[molb.type]);
1488 if (!bMoreThanTwoSeq &&
1489 more_than_two_sequential_constraints(molt.ilist, at2con[molb.type]))
1491 bMoreThanTwoSeq = TRUE;
1495 /* Check if we need to communicate not only before LINCS,
1496 * but also before each iteration.
1497 * The check for only two sequential constraints is only
1498 * useful for the common case of H-bond only constraints.
1499 * With more effort we could also make it useful for small
1500 * molecules with nr. sequential constraints <= nOrder-1.
1502 li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1504 if (debug && bPLINCS)
1506 fprintf(debug, "PLINCS communication before each iteration: %d\n",
1507 static_cast<int>(li->bCommIter));
1510 /* LINCS can run on any number of threads.
1511 * Currently the number is fixed for the whole simulation,
1512 * but it could be set in set_lincs().
1513 * The current constraint to task assignment code can create independent
1514 * tasks only when not more than two constraints are connected sequentially.
1516 li->ntask = gmx_omp_nthreads_get(emntLINCS);
1517 li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1520 fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n",
1521 li->ntask, li->bTaskDep ? "" : "in");
1529 /* Allocate an extra elements for "task-overlap" constraints */
1530 li->task.resize(li->ntask + 1);
1533 if (bPLINCS || li->ncg_triangle > 0)
1535 please_cite(fplog, "Hess2008a");
1539 please_cite(fplog, "Hess97a");
1544 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1547 fprintf(fplog, "There are constraints between atoms in different decomposition domains,\n"
1548 "will communicate selected coordinates each lincs iteration\n");
1550 if (li->ncg_triangle > 0)
1553 "%d constraints are involved in constraint triangles,\n"
1554 "will apply an additional matrix expansion of order %d for couplings\n"
1555 "between constraints inside triangles\n",
1556 li->ncg_triangle, li->nOrder);
1563 void done_lincs(Lincs *li)
1568 /*! \brief Sets up the work division over the threads. */
1569 static void lincs_thread_setup(Lincs *li, int natoms)
1571 li->atf.resize(natoms);
1573 gmx::ArrayRef<gmx_bitmask_t> atf = li->atf;
1575 /* Clear the atom flags */
1576 for (gmx_bitmask_t &mask : atf)
1578 bitmask_clear(&mask);
1581 if (li->ntask > BITMASK_SIZE)
1583 gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1586 for (int th = 0; th < li->ntask; th++)
1588 const Task *li_task = &li->task[th];
1590 /* For each atom set a flag for constraints from each */
1591 for (int b = li_task->b0; b < li_task->b1; b++)
1593 bitmask_set_bit(&atf[li->atoms[b].index1], th);
1594 bitmask_set_bit(&atf[li->atoms[b].index2], th);
1598 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1599 for (int th = 0; th < li->ntask; th++)
1607 li_task = &li->task[th];
1609 bitmask_init_low_bits(&mask, th);
1611 li_task->ind.clear();
1612 li_task->ind_r.clear();
1613 for (b = li_task->b0; b < li_task->b1; b++)
1615 /* We let the constraint with the lowest thread index
1616 * operate on atoms with constraints from multiple threads.
1618 if (bitmask_is_disjoint(atf[li->atoms[b].index1], mask) &&
1619 bitmask_is_disjoint(atf[li->atoms[b].index2], mask))
1621 /* Add the constraint to the local atom update index */
1622 li_task->ind.push_back(b);
1626 /* Add the constraint to the rest block */
1627 li_task->ind_r.push_back(b);
1631 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1634 /* We need to copy all constraints which have not be assigned
1635 * to a thread to a separate list which will be handled by one thread.
1637 Task *li_m = &li->task[li->ntask];
1640 for (int th = 0; th < li->ntask; th++)
1642 const Task &li_task = li->task[th];
1644 for (int ind_r : li_task.ind_r)
1646 li_m->ind.push_back(ind_r);
1651 fprintf(debug, "LINCS thread %d: %zu constraints\n",
1652 th, li_task.ind.size());
1658 fprintf(debug, "LINCS thread r: %zu constraints\n",
1663 //! Assign a constraint.
1664 static void assign_constraint(Lincs *li,
1665 int constraint_index,
1667 real lenA, real lenB,
1668 const t_blocka *at2con)
1674 /* Make an mapping of local topology constraint index to LINCS index */
1675 li->con_index[constraint_index] = con;
1677 li->bllen0[con] = lenA;
1678 li->ddist[con] = lenB - lenA;
1679 /* Set the length to the topology A length */
1680 li->bllen[con] = lenA;
1681 li->atoms[con].index1 = a1;
1682 li->atoms[con].index2 = a2;
1684 /* Make space in the constraint connection matrix for constraints
1685 * connected to both end of the current constraint.
1688 at2con->index[a1 + 1] - at2con->index[a1] - 1 +
1689 at2con->index[a2 + 1] - at2con->index[a2] - 1;
1691 li->blnr[con + 1] = li->ncc;
1693 /* Increase the constraint count */
1697 /*! \brief Check if constraint with topology index constraint_index is connected
1698 * to other constraints, and if so add those connected constraints to our task. */
1699 static void check_assign_connected(Lincs *li,
1700 const t_iatom *iatom,
1704 const t_blocka *at2con)
1706 /* Currently this function only supports constraint groups
1707 * in which all constraints share at least one atom
1708 * (e.g. H-bond constraints).
1709 * Check both ends of the current constraint for
1710 * connected constraints. We need to assign those
1715 for (end = 0; end < 2; end++)
1719 a = (end == 0 ? a1 : a2);
1721 for (k = at2con->index[a]; k < at2con->index[a + 1]; k++)
1726 /* Check if constraint cc has not yet been assigned */
1727 if (li->con_index[cc] == -1)
1733 lenA = idef.iparams[type].constr.dA;
1734 lenB = idef.iparams[type].constr.dB;
1736 if (bDynamics || lenA != 0 || lenB != 0)
1738 assign_constraint(li, cc, iatom[3*cc + 1], iatom[3*cc + 2], lenA, lenB, at2con);
1745 /*! \brief Check if constraint with topology index constraint_index is involved
1746 * in a constraint triangle, and if so add the other two constraints
1747 * in the triangle to our task. */
1748 static void check_assign_triangle(Lincs *li,
1749 const t_iatom *iatom,
1752 int constraint_index,
1754 const t_blocka *at2con)
1756 int nca, cc[32], ca[32], k;
1757 int c_triangle[2] = { -1, -1 };
1760 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1765 if (c != constraint_index)
1769 aa1 = iatom[c*3 + 1];
1770 aa2 = iatom[c*3 + 2];
1786 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1791 if (c != constraint_index)
1795 aa1 = iatom[c*3 + 1];
1796 aa2 = iatom[c*3 + 2];
1799 for (i = 0; i < nca; i++)
1803 c_triangle[0] = cc[i];
1810 for (i = 0; i < nca; i++)
1814 c_triangle[0] = cc[i];
1822 if (c_triangle[0] >= 0)
1826 for (end = 0; end < 2; end++)
1828 /* Check if constraint c_triangle[end] has not yet been assigned */
1829 if (li->con_index[c_triangle[end]] == -1)
1834 i = c_triangle[end]*3;
1836 lenA = idef.iparams[type].constr.dA;
1837 lenB = idef.iparams[type].constr.dB;
1839 if (bDynamics || lenA != 0 || lenB != 0)
1841 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1848 //! Sets matrix indices.
1849 static void set_matrix_indices(Lincs *li,
1850 const Task &li_task,
1851 const t_blocka *at2con,
1854 for (int b = li_task.b0; b < li_task.b1; b++)
1856 const int a1 = li->atoms[b].index1;
1857 const int a2 = li->atoms[b].index2;
1859 int i = li->blnr[b];
1860 for (int k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1862 int concon = li->con_index[at2con->a[k]];
1865 li->blbnb[i++] = concon;
1868 for (int k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1870 int concon = li->con_index[at2con->a[k]];
1873 li->blbnb[i++] = concon;
1879 /* Order the blbnb matrix to optimize memory access */
1880 std::sort(&(li->blbnb[li->blnr[b]]), &(li->blbnb[li->blnr[b+1]]));
1885 void set_lincs(const t_idef &idef,
1886 const t_mdatoms &md,
1888 const t_commrec *cr,
1898 /* Zero the thread index ranges.
1899 * Otherwise without local constraints we could return with old ranges.
1901 for (int i = 0; i < li->ntask; i++)
1905 li->task[i].ind.clear();
1909 li->task[li->ntask].ind.clear();
1912 /* This is the local topology, so there are only F_CONSTR constraints */
1913 if (idef.il[F_CONSTR].nr == 0)
1915 /* There are no constraints,
1916 * we do not need to fill any data structures.
1923 fprintf(debug, "Building the LINCS connectivity\n");
1926 if (DOMAINDECOMP(cr))
1928 if (cr->dd->constraints)
1932 dd_get_constraint_range(cr->dd, &start, &natoms);
1936 natoms = dd_numHomeAtoms(*cr->dd);
1944 at2con = make_at2con(natoms, idef.il, idef.iparams,
1945 flexibleConstraintTreatment(bDynamics));
1947 const int ncon_tot = idef.il[F_CONSTR].nr/3;
1949 /* Ensure we have enough padding for aligned loads for each thread */
1950 const int numEntries = ncon_tot + li->ntask*simd_width;
1951 li->con_index.resize(numEntries);
1952 li->bllen0.resize(numEntries);
1953 li->ddist.resize(numEntries);
1954 li->atoms.resize(numEntries);
1955 li->blc.resize(numEntries);
1956 li->blc1.resize(numEntries);
1957 li->blnr.resize(numEntries + 1);
1958 li->bllen.resize(numEntries);
1959 li->tmpv.resizeWithPadding(numEntries);
1960 if (DOMAINDECOMP(cr))
1962 li->nlocat.resize(numEntries);
1964 li->tmp1.resize(numEntries);
1965 li->tmp2.resize(numEntries);
1966 li->tmp3.resize(numEntries);
1967 li->tmp4.resize(numEntries);
1968 li->mlambda.resize(numEntries);
1970 iatom = idef.il[F_CONSTR].iatoms;
1972 li->blnr[0] = li->ncc;
1974 /* Assign the constraints for li->ntask LINCS tasks.
1975 * We target a uniform distribution of constraints over the tasks.
1976 * Note that when flexible constraints are present, but are removed here
1977 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
1978 * happen during normal MD, that's ok.
1981 /* Determine the number of constraints we need to assign here */
1982 int ncon_assign = ncon_tot;
1985 /* With energy minimization, flexible constraints are ignored
1986 * (and thus minimized, as they should be).
1988 ncon_assign -= countFlexibleConstraints(idef.il, idef.iparams);
1991 /* Set the target constraint count per task to exactly uniform,
1992 * this might be overridden below.
1994 int ncon_target = (ncon_assign + li->ntask - 1)/li->ntask;
1996 /* Mark all constraints as unassigned by setting their index to -1 */
1997 for (int con = 0; con < ncon_tot; con++)
1999 li->con_index[con] = -1;
2003 for (int th = 0; th < li->ntask; th++)
2007 li_task = &li->task[th];
2009 #if GMX_SIMD_HAVE_REAL
2010 /* With indepedent tasks we likely have H-bond constraints or constraint
2011 * pairs. The connected constraints will be pulled into the task, so the
2012 * constraints per task will often exceed ncon_target.
2013 * Triangle constraints can also increase the count, but there are
2014 * relatively few of those, so we usually expect to get ncon_target.
2018 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
2019 * since otherwise a lot of operations can be wasted.
2020 * There are several ways to round here, we choose the one
2021 * that alternates block sizes, which helps with Intel HT.
2023 ncon_target = ((ncon_assign*(th + 1))/li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1) & ~(GMX_SIMD_REAL_WIDTH - 1);
2025 #endif // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
2027 /* Continue filling the arrays where we left off with the previous task,
2028 * including padding for SIMD.
2030 li_task->b0 = li->nc;
2032 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2034 if (li->con_index[con] == -1)
2039 type = iatom[3*con];
2040 a1 = iatom[3*con + 1];
2041 a2 = iatom[3*con + 2];
2042 lenA = idef.iparams[type].constr.dA;
2043 lenB = idef.iparams[type].constr.dB;
2044 /* Skip the flexible constraints when not doing dynamics */
2045 if (bDynamics || lenA != 0 || lenB != 0)
2047 assign_constraint(li, con, a1, a2, lenA, lenB, &at2con);
2049 if (li->ntask > 1 && !li->bTaskDep)
2051 /* We can generate independent tasks. Check if we
2052 * need to assign connected constraints to our task.
2054 check_assign_connected(li, iatom, idef, bDynamics,
2057 if (li->ntask > 1 && li->ncg_triangle > 0)
2059 /* Ensure constraints in one triangle are assigned
2062 check_assign_triangle(li, iatom, idef, bDynamics,
2063 con, a1, a2, &at2con);
2071 li_task->b1 = li->nc;
2075 /* Copy the last atom pair indices and lengths for constraints
2076 * up to a multiple of simd_width, such that we can do all
2077 * SIMD operations without having to worry about end effects.
2081 li->nc = ((li_task->b1 + simd_width - 1)/simd_width)*simd_width;
2082 last = li_task->b1 - 1;
2083 for (i = li_task->b1; i < li->nc; i++)
2085 li->atoms[i] = li->atoms[last];
2086 li->bllen0[i] = li->bllen0[last];
2087 li->ddist[i] = li->ddist[last];
2088 li->bllen[i] = li->bllen[last];
2089 li->blnr[i + 1] = li->blnr[last + 1];
2093 /* Keep track of how many constraints we assigned */
2094 li->nc_real += li_task->b1 - li_task->b0;
2098 fprintf(debug, "LINCS task %d constraints %d - %d\n",
2099 th, li_task->b0, li_task->b1);
2103 assert(li->nc_real == ncon_assign);
2107 /* Without DD we order the blbnb matrix to optimize memory access.
2108 * With DD the overhead of sorting is more than the gain during access.
2110 bSortMatrix = !DOMAINDECOMP(cr);
2112 li->blbnb.resize(li->ncc);
2114 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2115 for (int th = 0; th < li->ntask; th++)
2119 Task &li_task = li->task[th];
2121 if (li->ncg_triangle > 0)
2123 /* This is allocating too much, but it is difficult to improve */
2124 li_task.triangle.resize(li_task.b1 - li_task.b0);
2125 li_task.tri_bits.resize(li_task.b1 - li_task.b0);
2128 set_matrix_indices(li, li_task, &at2con, bSortMatrix);
2130 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2133 done_blocka(&at2con);
2135 if (cr->dd == nullptr)
2137 /* Since the matrix is static, we should free some memory */
2138 li->blbnb.resize(li->ncc);
2141 li->blmf.resize(li->ncc);
2142 li->blmf1.resize(li->ncc);
2143 li->tmpncc.resize(li->ncc);
2145 gmx::ArrayRef<const int> nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2146 if (!nlocat_dd.empty())
2148 /* Convert nlocat from local topology to LINCS constraint indexing */
2149 for (con = 0; con < ncon_tot; con++)
2151 li->nlocat[li->con_index[con]] = nlocat_dd[con];
2161 fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n",
2162 li->nc_real, li->nc, li->ncc);
2167 lincs_thread_setup(li, md.nr);
2170 set_lincs_matrix(li, md.invmass, md.lambda);
2173 //! Issues a warning when LINCS constraints cannot be satisfied.
2174 static void lincs_warning(gmx_domdec_t *dd, const rvec *x, rvec *xprime, t_pbc *pbc,
2175 int ncons, gmx::ArrayRef<const AtomPair> atoms,
2176 gmx::ArrayRef<const real> bllen, real wangle,
2177 int maxwarn, int *warncount)
2179 real wfac = std::cos(DEG2RAD*wangle);
2182 "bonds that rotated more than %g degrees:\n"
2183 " atom 1 atom 2 angle previous, current, constraint length\n",
2186 for (int b = 0; b < ncons; b++)
2188 const int i = atoms[b].index1;
2189 const int j = atoms[b].index2;
2194 pbc_dx_aiuc(pbc, x[i], x[j], v0);
2195 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2199 rvec_sub(x[i], x[j], v0);
2200 rvec_sub(xprime[i], xprime[j], v1);
2204 real cosine = ::iprod(v0, v1)/(d0*d1);
2208 " %6d %6d %5.1f %8.4f %8.4f %8.4f\n",
2209 ddglatnr(dd, i), ddglatnr(dd, j),
2210 RAD2DEG*std::acos(cosine), d0, d1, bllen[b]);
2211 if (!std::isfinite(d1))
2213 gmx_fatal(FARGS, "Bond length not finite.");
2219 if (*warncount > maxwarn)
2221 too_many_constraint_warnings(econtLINCS, *warncount);
2225 //! Determine how well the constraints have been satisfied.
2226 static void cconerr(const Lincs &lincsd,
2227 const rvec *x, const t_pbc *pbc,
2228 real *ncons_loc, real *ssd, real *max, int *imax)
2230 gmx::ArrayRef<const AtomPair> atoms = lincsd.atoms;
2231 gmx::ArrayRef<const real> bllen = lincsd.bllen;
2232 gmx::ArrayRef<const int> nlocat = lincsd.nlocat;
2238 for (int task = 0; task < lincsd.ntask; task++)
2240 for (int b = lincsd.task[task].b0; b < lincsd.task[task].b1; b++)
2245 pbc_dx_aiuc(pbc, x[atoms[b].index1], x[atoms[b].index2], dx);
2249 rvec_sub(x[atoms[b].index1], x[atoms[b].index2], dx);
2251 real r2 = ::norm2(dx);
2252 real len = r2*gmx::invsqrt(r2);
2253 real d = std::abs(len/bllen[b]-1);
2254 if (d > ma && (nlocat.empty() || nlocat[b]))
2266 ssd2 += nlocat[b]*d*d;
2272 *ncons_loc = (nlocat.empty() ? 1 : 0.5)*count;
2273 *ssd = (nlocat.empty() ? 1 : 0.5)*ssd2;
2278 bool constrain_lincs(bool computeRmsd,
2279 const t_inputrec &ir,
2281 Lincs *lincsd, const t_mdatoms &md,
2282 const t_commrec *cr,
2283 const gmx_multisim_t *ms,
2284 const rvec *x, rvec *xprime, rvec *min_proj,
2285 matrix box, t_pbc *pbc,
2286 real lambda, real *dvdlambda,
2287 real invdt, rvec *v,
2288 bool bCalcVir, tensor vir_r_m_dr,
2289 ConstraintVariable econq,
2291 int maxwarn, int *warncount)
2295 /* This boolean should be set by a flag passed to this routine.
2296 * We can also easily check if any constraint length is changed,
2297 * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2299 bool bCalcDHDL = (ir.efep != efepNO && dvdlambda != nullptr);
2301 if (lincsd->nc == 0 && cr->dd == nullptr)
2305 lincsd->rmsdData = {{0}};
2311 if (econq == ConstraintVariable::Positions)
2313 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2314 * also with efep!=fepNO.
2316 if (ir.efep != efepNO)
2318 if (md.nMassPerturbed && lincsd->matlam != md.lambda)
2320 set_lincs_matrix(lincsd, md.invmass, md.lambda);
2323 for (int i = 0; i < lincsd->nc; i++)
2325 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
2329 if (lincsd->ncg_flex)
2331 /* Set the flexible constraint lengths to the old lengths */
2334 for (int i = 0; i < lincsd->nc; i++)
2336 if (lincsd->bllen[i] == 0)
2339 pbc_dx_aiuc(pbc, x[lincsd->atoms[i].index1], x[lincsd->atoms[i].index2], dx);
2340 lincsd->bllen[i] = norm(dx);
2346 for (int i = 0; i < lincsd->nc; i++)
2348 if (lincsd->bllen[i] == 0)
2351 std::sqrt(distance2(x[lincsd->atoms[i].index1],
2352 x[lincsd->atoms[i].index2]));
2364 cconerr(*lincsd, xprime, pbc,
2365 &ncons_loc, &p_ssd, &p_max, &p_imax);
2368 /* This bWarn var can be updated by multiple threads
2369 * at the same time. But as we only need to detect
2370 * if a warning occurred or not, this is not an issue.
2374 /* The OpenMP parallel region of constrain_lincs for coords */
2375 #pragma omp parallel num_threads(lincsd->ntask)
2379 int th = gmx_omp_get_thread_num();
2381 clear_mat(lincsd->task[th].vir_r_m_dr);
2383 do_lincs(x, xprime, box, pbc, lincsd, th,
2386 ir.LincsWarnAngle, &bWarn,
2388 th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2390 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2393 if (debug && lincsd->nc > 0)
2395 fprintf(debug, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2396 fprintf(debug, " Before LINCS %.6f %.6f %6d %6d\n",
2397 std::sqrt(p_ssd/ncons_loc), p_max,
2398 ddglatnr(cr->dd, lincsd->atoms[p_imax].index1),
2399 ddglatnr(cr->dd, lincsd->atoms[p_imax].index2));
2401 if (computeRmsd || debug)
2403 cconerr(*lincsd, xprime, pbc,
2404 &ncons_loc, &p_ssd, &p_max, &p_imax);
2405 lincsd->rmsdData[0] = ncons_loc;
2406 lincsd->rmsdData[1] = p_ssd;
2410 lincsd->rmsdData = {{0}};
2412 if (debug && lincsd->nc > 0)
2415 " After LINCS %.6f %.6f %6d %6d\n\n",
2416 std::sqrt(p_ssd/ncons_loc), p_max,
2417 ddglatnr(cr->dd, lincsd->atoms[p_imax].index1),
2418 ddglatnr(cr->dd, lincsd->atoms[p_imax].index2));
2423 if (maxwarn < INT_MAX)
2425 cconerr(*lincsd, xprime, pbc,
2426 &ncons_loc, &p_ssd, &p_max, &p_imax);
2427 std::string simMesg;
2430 simMesg += gmx::formatString(" in simulation %d", ms->sim);
2433 "\nStep %" PRId64 ", time %g (ps) LINCS WARNING%s\n"
2434 "relative constraint deviation after LINCS:\n"
2435 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2436 step, ir.init_t+step*ir.delta_t,
2438 std::sqrt(p_ssd/ncons_loc), p_max,
2439 ddglatnr(cr->dd, lincsd->atoms[p_imax].index1),
2440 ddglatnr(cr->dd, lincsd->atoms[p_imax].index2));
2442 lincs_warning(cr->dd, x, xprime, pbc,
2443 lincsd->nc, lincsd->atoms, lincsd->bllen,
2444 ir.LincsWarnAngle, maxwarn, warncount);
2446 bOK = (p_max < 0.5);
2449 if (lincsd->ncg_flex)
2451 for (int i = 0; (i < lincsd->nc); i++)
2453 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2455 lincsd->bllen[i] = 0;
2462 /* The OpenMP parallel region of constrain_lincs for derivatives */
2463 #pragma omp parallel num_threads(lincsd->ntask)
2467 int th = gmx_omp_get_thread_num();
2469 do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
2470 md.invmass, econq, bCalcDHDL,
2471 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2473 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2479 /* Reduce the dH/dlambda contributions over the threads */
2484 for (th = 0; th < lincsd->ntask; th++)
2486 dhdlambda += lincsd->task[th].dhdlambda;
2488 if (econq == ConstraintVariable::Positions)
2490 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2491 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2492 dhdlambda /= ir.delta_t*ir.delta_t;
2494 *dvdlambda += dhdlambda;
2497 if (bCalcVir && lincsd->ntask > 1)
2499 for (int i = 1; i < lincsd->ntask; i++)
2501 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2505 /* count assuming nit=1 */
2506 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2507 inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
2508 if (lincsd->ntriangle > 0)
2510 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
2514 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real*2);
2518 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);