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/units.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/constr.h"
64 #include "gromacs/mdlib/gmx_omp_nthreads.h"
65 #include "gromacs/mdlib/mdrun.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/arrayref.h"
78 #include "gromacs/utility/basedefinitions.h"
79 #include "gromacs/utility/bitmask.h"
80 #include "gromacs/utility/cstringutil.h"
81 #include "gromacs/utility/exceptions.h"
82 #include "gromacs/utility/fatalerror.h"
83 #include "gromacs/utility/gmxomp.h"
84 #include "gromacs/utility/pleasecite.h"
85 #include "gromacs/utility/smalloc.h"
87 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
92 //! Unit of work within LINCS.
95 //! First constraint for this task.
97 //! b1-1 is the last constraint for this task.
99 //! The number of constraints in triangles.
101 //! The list of triangle constraints.
102 int *triangle = nullptr;
103 //! The bits tell if the matrix element should be used.
104 int *tri_bits = nullptr;
105 //! Allocation size of triangle and tri_bits.
107 //! Number of indices.
109 //! Constraint index for updating atom data.
111 //! Number of indices.
113 //! Constraint index for updating atom data.
114 int *ind_r = nullptr;
115 //! Allocation size of ind and ind_r.
117 //! Temporary variable for virial calculation.
118 tensor vir_r_m_dr = {{0}};
119 //! Temporary variable for lambda derivative.
123 /*! \brief Data for LINCS algorithm.
128 //! The global number of constraints.
130 //! The global number of flexible constraints.
132 //! The global number of constraints in triangles.
133 int ncg_triangle = 0;
134 //! The number of iterations.
136 //! The order of the matrix expansion.
138 //! The maximum number of constraints connected to a single atom.
141 //! The number of real constraints.
143 //! The number of constraints including padding for SIMD.
145 //! The number we allocated memory for.
147 //! The number of constraint connections.
149 //! The number we allocated memory for.
151 //! The FE lambda value used for filling blc and blmf.
153 //! mapping from topology to LINCS constraints.
154 int *con_index = nullptr;
155 //! The reference distance in topology A.
156 real *bllen0 = nullptr;
157 //! The reference distance in top B - the r.d. in top A.
158 real *ddist = nullptr;
159 //! The atom pairs involved in the constraints.
161 //! 1/sqrt(invmass1 invmass2).
163 //! As blc, but with all masses 1.
164 real *blc1 = nullptr;
165 //! Index into blbnb and blmf.
167 //! List of constraint connections.
168 int *blbnb = nullptr;
169 //! The local number of constraints in triangles.
171 //! The number of constraint connections in triangles.
172 int ncc_triangle = 0;
173 //! Communicate before each LINCS interation.
174 bool bCommIter = false;
175 //! Matrix of mass factors for constraint connections.
176 real *blmf = nullptr;
177 //! As blmf, but with all masses 1.
178 real *blmf1 = nullptr;
179 //! The reference bond length.
180 real *bllen = nullptr;
181 //! The local atom count per constraint, can be NULL.
182 int *nlocat = nullptr;
184 /*! \brief The number of tasks used for LINCS work.
186 * \todo This is mostly used to loop over \c task, which would
187 * be nicer to do with range-based for loops, but the thread
188 * index is used for constructing bit masks and organizing the
189 * virial output buffer, so other things need to change,
192 /*! \brief LINCS thread division */
193 std::vector<Task> task;
194 //! Atom flags for thread parallelization.
195 gmx_bitmask_t *atf = nullptr;
196 //! Allocation size of 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 rvec *tmpv = nullptr;
205 real *tmpncc = nullptr;
206 real *tmp1 = nullptr;
207 real *tmp2 = nullptr;
208 real *tmp3 = nullptr;
209 real *tmp4 = nullptr;
211 //! The Lagrange multipliers times -1.
212 real *mlambda = nullptr;
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 /*! \brief Align to 128 bytes, consistent with the current implementation of
225 AlignedAllocator, which currently forces 128 byte alignment. */
226 static const int align_bytes = 128;
228 ArrayRef<real> lincs_rmsdData(Lincs *lincsd)
230 return lincsd->rmsdData;
233 real lincs_rmsd(const Lincs *lincsd)
235 if (lincsd->rmsdData[0] > 0)
237 return std::sqrt(lincsd->rmsdData[1]/lincsd->rmsdData[0]);
245 /*! \brief Do a set of nrec LINCS matrix multiplications.
247 * This function will return with up to date thread-local
248 * constraint data, without an OpenMP barrier.
250 static void lincs_matrix_expand(const Lincs *lincsd,
253 real *rhs1, real *rhs2, real *sol)
255 int b0, b1, nrec, rec;
256 const int *blnr = lincsd->blnr;
257 const int *blbnb = lincsd->blbnb;
261 nrec = lincsd->nOrder;
263 for (rec = 0; rec < nrec; rec++)
267 if (lincsd->bTaskDep)
271 for (b = b0; b < b1; b++)
277 for (n = blnr[b]; n < blnr[b+1]; n++)
279 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
282 sol[b] = sol[b] + mvb;
290 } /* nrec*(ncons+2*nrtot) flops */
292 if (lincsd->ntriangle > 0)
294 /* Perform an extra nrec recursions for only the constraints
295 * involved in rigid triangles.
296 * In this way their accuracy should come close to those of the other
297 * constraints, since traingles of constraints can produce eigenvalues
298 * around 0.7, while the effective eigenvalue for bond constraints
299 * is around 0.4 (and 0.7*0.7=0.5).
302 if (lincsd->bTaskDep)
304 /* We need a barrier here, since other threads might still be
305 * reading the contents of rhs1 and/o rhs2.
306 * We could avoid this barrier by introducing two extra rhs
307 * arrays for the triangle constraints only.
312 /* Constraints involved in a triangle are ensured to be in the same
313 * LINCS task. This means no barriers are required during the extra
314 * iterations for the triangle constraints.
316 const int *triangle = li_task->triangle;
317 const int *tri_bits = li_task->tri_bits;
319 for (rec = 0; rec < nrec; rec++)
323 for (tb = 0; tb < li_task->ntriangle; tb++)
325 int b, bits, nr0, nr1, n;
333 for (n = nr0; n < nr1; n++)
335 if (bits & (1 << (n - nr0)))
337 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
341 sol[b] = sol[b] + mvb;
349 } /* nrec*(ntriangle + ncc_triangle*2) flops */
351 if (lincsd->bTaskDepTri)
353 /* The constraints triangles are decoupled from each other,
354 * but constraints in one triangle cross thread task borders.
355 * We could probably avoid this with more advanced setup code.
362 //! Update atomic coordinates when an index is not required.
363 static void lincs_update_atoms_noind(int ncons, const int *bla,
365 const real *fac, rvec *r,
370 real mvb, im1, im2, tmp0, tmp1, tmp2;
372 if (invmass != nullptr)
374 for (b = 0; b < ncons; b++)
390 } /* 16 ncons flops */
394 for (b = 0; b < ncons; b++)
412 //! Update atomic coordinates when an index is required.
413 static void lincs_update_atoms_ind(int ncons, const int *ind, const int *bla,
415 const real *fac, rvec *r,
420 real mvb, im1, im2, tmp0, tmp1, tmp2;
422 if (invmass != nullptr)
424 for (bi = 0; bi < ncons; bi++)
441 } /* 16 ncons flops */
445 for (bi = 0; bi < ncons; bi++)
460 } /* 16 ncons flops */
464 //! Update coordinates for atoms.
465 static void lincs_update_atoms(Lincs *li, int th,
467 const real *fac, rvec *r,
473 /* Single thread, we simply update for all constraints */
474 lincs_update_atoms_noind(li->nc_real,
475 li->bla, prefac, fac, r, invmass, x);
479 /* Update the atom vector components for our thread local
480 * constraints that only access our local atom range.
481 * This can be done without a barrier.
483 lincs_update_atoms_ind(li->task[th].nind, li->task[th].ind,
484 li->bla, prefac, fac, r, invmass, x);
486 if (li->task[li->ntask].nind > 0)
488 /* Update the constraints that operate on atoms
489 * in multiple thread atom blocks on the master thread.
494 lincs_update_atoms_ind(li->task[li->ntask].nind,
495 li->task[li->ntask].ind,
496 li->bla, prefac, fac, r, invmass, x);
502 #if GMX_SIMD_HAVE_REAL
503 //! Helper function so that we can run TSAN with SIMD support (where implemented).
505 static inline void gmx_simdcall
506 gatherLoadUTransposeTSANSafe(const real *base,
507 const std::int32_t *offset,
512 #if (CMAKE_BUILD_TYPE == CMAKE_BUILD_TYPE_TSAN) && GMX_SIMD_X86_AVX2_256
513 // This function is only implemented in this case
514 gatherLoadUTransposeSafe<align>(base, offset, v0, v1, v2);
516 gatherLoadUTranspose<align>(base, offset, v0, v1, v2);
520 /*! \brief Calculate the constraint distance vectors r to project on from x.
522 * Determine the right-hand side of the matrix equation using quantity f.
523 * This function only differs from calc_dr_x_xp_simd below in that
524 * no constraint length is subtracted and no PBC is used for f. */
525 static void gmx_simdcall
526 calc_dr_x_f_simd(int b0,
529 const rvec * gmx_restrict x,
530 const rvec * gmx_restrict f,
531 const real * gmx_restrict blc,
532 const real * pbc_simd,
533 rvec * gmx_restrict r,
534 real * gmx_restrict rhs,
535 real * gmx_restrict sol)
537 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
539 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
541 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
546 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
548 SimdReal x0_S, y0_S, z0_S;
549 SimdReal x1_S, y1_S, z1_S;
550 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
551 SimdReal fx_S, fy_S, fz_S, ip_S, rhs_S;
552 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
553 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
555 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
557 offset0[i] = bla[bs*2 + i*2];
558 offset1[i] = bla[bs*2 + i*2 + 1];
561 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
562 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
567 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
569 n2_S = norm2(rx_S, ry_S, rz_S);
570 il_S = invsqrt(n2_S);
576 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
578 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(f), offset0, &x0_S, &y0_S, &z0_S);
579 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(f), offset1, &x1_S, &y1_S, &z1_S);
584 ip_S = iprod(rx_S, ry_S, rz_S, fx_S, fy_S, fz_S);
586 rhs_S = load<SimdReal>(blc + bs) * ip_S;
588 store(rhs + bs, rhs_S);
589 store(sol + bs, rhs_S);
592 #endif // GMX_SIMD_HAVE_REAL
594 /*! \brief LINCS projection, works on derivatives of the coordinates. */
595 static void do_lincsp(const rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
596 Lincs *lincsd, int th,
598 ConstraintVariable econq, bool bCalcDHDL,
599 bool bCalcVir, tensor rmdf)
602 int *bla, *blnr, *blbnb;
604 real *blc, *blmf, *blcc, *rhs1, *rhs2, *sol;
606 b0 = lincsd->task[th].b0;
607 b1 = lincsd->task[th].b1;
612 blbnb = lincsd->blbnb;
613 if (econq != ConstraintVariable::Force)
615 /* Use mass-weighted parameters */
621 /* Use non mass-weighted parameters */
623 blmf = lincsd->blmf1;
625 blcc = lincsd->tmpncc;
630 #if GMX_SIMD_HAVE_REAL
631 /* This SIMD code does the same as the plain-C code after the #else.
632 * The only difference is that we always call pbc code, as with SIMD
633 * the overhead of pbc computation (when not needed) is small.
635 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
637 /* Convert the pbc struct for SIMD */
638 set_pbc_simd(pbc, pbc_simd);
640 /* Compute normalized x i-j vectors, store in r.
641 * Compute the inner product of r and xp i-j and store in rhs1.
643 calc_dr_x_f_simd(b0, b1, bla, x, f, blc,
647 #else // GMX_SIMD_HAVE_REAL
649 /* Compute normalized i-j vectors */
652 for (b = b0; b < b1; b++)
656 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
662 for (b = b0; b < b1; b++)
666 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
668 } /* 16 ncons flops */
671 for (b = b0; b < b1; b++)
678 mvb = blc[b]*(r[b][0]*(f[i][0] - f[j][0]) +
679 r[b][1]*(f[i][1] - f[j][1]) +
680 r[b][2]*(f[i][2] - f[j][2]));
686 #endif // GMX_SIMD_HAVE_REAL
688 if (lincsd->bTaskDep)
690 /* We need a barrier, since the matrix construction below
691 * can access entries in r of other threads.
696 /* Construct the (sparse) LINCS matrix */
697 for (b = b0; b < b1; b++)
701 for (n = blnr[b]; n < blnr[b+1]; n++)
703 blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
706 /* Together: 23*ncons + 6*nrtot flops */
708 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
709 /* nrec*(ncons+2*nrtot) flops */
711 if (econq == ConstraintVariable::Deriv_FlexCon)
713 /* We only want to constraint the flexible constraints,
714 * so we mask out the normal ones by setting sol to 0.
716 for (b = b0; b < b1; b++)
718 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
725 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
726 for (b = b0; b < b1; b++)
731 /* When constraining forces, we should not use mass weighting,
732 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
734 lincs_update_atoms(lincsd, th, 1.0, sol, r,
735 (econq != ConstraintVariable::Force) ? invmass : nullptr, fp);
742 for (b = b0; b < b1; b++)
744 dhdlambda -= sol[b]*lincsd->ddist[b];
747 lincsd->task[th].dhdlambda = dhdlambda;
752 /* Constraint virial,
753 * determines sum r_bond x delta f,
754 * where delta f is the constraint correction
755 * of the quantity that is being constrained.
757 for (b = b0; b < b1; b++)
762 mvb = lincsd->bllen[b]*sol[b];
763 for (i = 0; i < DIM; i++)
766 for (j = 0; j < DIM; j++)
768 rmdf[i][j] += tmp1*r[b][j];
771 } /* 23 ncons flops */
775 #if GMX_SIMD_HAVE_REAL
777 /*! \brief Calculate the constraint distance vectors r to project on from x.
779 * Determine the right-hand side of the matrix equation using coordinates xp. */
780 static void gmx_simdcall
781 calc_dr_x_xp_simd(int b0,
784 const rvec * gmx_restrict x,
785 const rvec * gmx_restrict xp,
786 const real * gmx_restrict bllen,
787 const real * gmx_restrict blc,
788 const real * pbc_simd,
789 rvec * gmx_restrict r,
790 real * gmx_restrict rhs,
791 real * gmx_restrict sol)
793 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
794 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
796 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
801 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
803 SimdReal x0_S, y0_S, z0_S;
804 SimdReal x1_S, y1_S, z1_S;
805 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
806 SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
807 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
808 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
810 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
812 offset0[i] = bla[bs*2 + i*2];
813 offset1[i] = bla[bs*2 + i*2 + 1];
816 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
817 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
822 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
824 n2_S = norm2(rx_S, ry_S, rz_S);
825 il_S = invsqrt(n2_S);
831 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
833 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(xp), offset0, &x0_S, &y0_S, &z0_S);
834 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(xp), offset1, &x1_S, &y1_S, &z1_S);
840 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
842 ip_S = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
844 rhs_S = load<SimdReal>(blc + bs) * (ip_S - load<SimdReal>(bllen + bs));
846 store(rhs + bs, rhs_S);
847 store(sol + bs, rhs_S);
850 #endif // GMX_SIMD_HAVE_REAL
852 /*! \brief Determine the distances and right-hand side for the next iteration. */
853 gmx_unused static void calc_dist_iter(
857 const rvec * gmx_restrict xp,
858 const real * gmx_restrict bllen,
859 const real * gmx_restrict blc,
862 real * gmx_restrict rhs,
863 real * gmx_restrict sol,
868 for (b = b0; b < b1; b++)
870 real len, len2, dlen2, mvb;
876 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
880 rvec_sub(xp[bla[2*b]], xp[bla[2*b+1]], dx);
883 dlen2 = 2*len2 - ::norm2(dx);
884 if (dlen2 < wfac*len2)
886 /* not race free - see detailed comment in caller */
891 mvb = blc[b]*(len - dlen2*gmx::invsqrt(dlen2));
899 } /* 20*ncons flops */
902 #if GMX_SIMD_HAVE_REAL
903 /*! \brief As calc_dist_iter(), but using SIMD intrinsics. */
904 static void gmx_simdcall
905 calc_dist_iter_simd(int b0,
908 const rvec * gmx_restrict x,
909 const real * gmx_restrict bllen,
910 const real * gmx_restrict blc,
911 const real * pbc_simd,
913 real * gmx_restrict rhs,
914 real * gmx_restrict sol,
917 SimdReal min_S(GMX_REAL_MIN);
919 SimdReal wfac_S(wfac);
924 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
926 /* Initialize all to FALSE */
927 warn_S = (two_S < setZero());
929 for (bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
931 SimdReal x0_S, y0_S, z0_S;
932 SimdReal x1_S, y1_S, z1_S;
933 SimdReal rx_S, ry_S, rz_S, n2_S;
934 SimdReal len_S, len2_S, dlen2_S, lc_S, blc_S;
935 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
936 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
938 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
940 offset0[i] = bla[bs*2 + i*2];
941 offset1[i] = bla[bs*2 + i*2 + 1];
944 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
945 gatherLoadUTransposeTSANSafe<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
951 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
953 n2_S = norm2(rx_S, ry_S, rz_S);
955 len_S = load<SimdReal>(bllen + bs);
956 len2_S = len_S * len_S;
958 dlen2_S = fms(two_S, len2_S, n2_S);
960 warn_S = warn_S || (dlen2_S < (wfac_S * len2_S));
962 /* Avoid 1/0 by taking the max with REAL_MIN.
963 * Note: when dlen2 is close to zero (90 degree constraint rotation),
964 * the accuracy of the algorithm is no longer relevant.
966 dlen2_S = max(dlen2_S, min_S);
968 lc_S = fnma(dlen2_S, invsqrt(dlen2_S), len_S);
970 blc_S = load<SimdReal>(blc + bs);
974 store(rhs + bs, lc_S);
975 store(sol + bs, lc_S);
983 #endif // GMX_SIMD_HAVE_REAL
985 //! Implements LINCS constraining.
986 static void do_lincs(const rvec *x, rvec *xp, matrix box, t_pbc *pbc,
987 Lincs *lincsd, int th,
991 real wangle, bool *bWarn,
992 real invdt, rvec * gmx_restrict v,
993 bool bCalcVir, tensor vir_r_m_dr)
995 int b0, b1, b, i, j, n, iter;
996 int *bla, *blnr, *blbnb;
998 real *blc, *blmf, *bllen, *blcc, *rhs1, *rhs2, *sol, *blc_sol, *mlambda;
1001 b0 = lincsd->task[th].b0;
1002 b1 = lincsd->task[th].b1;
1006 blnr = lincsd->blnr;
1007 blbnb = lincsd->blbnb;
1009 blmf = lincsd->blmf;
1010 bllen = lincsd->bllen;
1011 blcc = lincsd->tmpncc;
1012 rhs1 = lincsd->tmp1;
1013 rhs2 = lincsd->tmp2;
1015 blc_sol = lincsd->tmp4;
1016 mlambda = lincsd->mlambda;
1017 nlocat = lincsd->nlocat;
1019 #if GMX_SIMD_HAVE_REAL
1021 /* This SIMD code does the same as the plain-C code after the #else.
1022 * The only difference is that we always call pbc code, as with SIMD
1023 * the overhead of pbc computation (when not needed) is small.
1025 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
1027 /* Convert the pbc struct for SIMD */
1028 set_pbc_simd(pbc, pbc_simd);
1030 /* Compute normalized x i-j vectors, store in r.
1031 * Compute the inner product of r and xp i-j and store in rhs1.
1033 calc_dr_x_xp_simd(b0, b1, bla, x, xp, bllen, blc,
1037 #else // GMX_SIMD_HAVE_REAL
1041 /* Compute normalized i-j vectors */
1042 for (b = b0; b < b1; b++)
1047 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
1050 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
1051 mvb = blc[b]*(::iprod(r[b], dx) - bllen[b]);
1058 /* Compute normalized i-j vectors */
1059 for (b = b0; b < b1; b++)
1061 real tmp0, tmp1, tmp2, rlen, mvb;
1065 tmp0 = x[i][0] - x[j][0];
1066 tmp1 = x[i][1] - x[j][1];
1067 tmp2 = x[i][2] - x[j][2];
1068 rlen = gmx::invsqrt(tmp0*tmp0 + tmp1*tmp1 + tmp2*tmp2);
1069 r[b][0] = rlen*tmp0;
1070 r[b][1] = rlen*tmp1;
1071 r[b][2] = rlen*tmp2;
1072 /* 16 ncons flops */
1076 mvb = blc[b]*(r[b][0]*(xp[i][0] - xp[j][0]) +
1077 r[b][1]*(xp[i][1] - xp[j][1]) +
1078 r[b][2]*(xp[i][2] - xp[j][2]) - bllen[b]);
1083 /* Together: 26*ncons + 6*nrtot flops */
1086 #endif // GMX_SIMD_HAVE_REAL
1088 if (lincsd->bTaskDep)
1090 /* We need a barrier, since the matrix construction below
1091 * can access entries in r of other threads.
1096 /* Construct the (sparse) LINCS matrix */
1097 for (b = b0; b < b1; b++)
1099 for (n = blnr[b]; n < blnr[b+1]; n++)
1101 blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
1104 /* Together: 26*ncons + 6*nrtot flops */
1106 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1107 /* nrec*(ncons+2*nrtot) flops */
1109 #if GMX_SIMD_HAVE_REAL
1110 for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1112 SimdReal t1 = load<SimdReal>(blc + b);
1113 SimdReal t2 = load<SimdReal>(sol + b);
1114 store(mlambda + b, t1 * t2);
1117 for (b = b0; b < b1; b++)
1119 mlambda[b] = blc[b]*sol[b];
1121 #endif // GMX_SIMD_HAVE_REAL
1123 /* Update the coordinates */
1124 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1127 ******** Correction for centripetal effects ********
1132 wfac = std::cos(DEG2RAD*wangle);
1135 for (iter = 0; iter < lincsd->nIter; iter++)
1137 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
1142 /* Communicate the corrected non-local coordinates */
1143 if (DOMAINDECOMP(cr))
1145 dd_move_x_constraints(cr->dd, box, xp, nullptr, FALSE);
1150 else if (lincsd->bTaskDep)
1155 #if GMX_SIMD_HAVE_REAL
1156 calc_dist_iter_simd(b0, b1, bla, xp, bllen, blc, pbc_simd, wfac,
1159 calc_dist_iter(b0, b1, bla, xp, bllen, blc, pbc, wfac,
1161 /* 20*ncons flops */
1162 #endif // GMX_SIMD_HAVE_REAL
1164 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1165 /* nrec*(ncons+2*nrtot) flops */
1167 #if GMX_SIMD_HAVE_REAL
1168 for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1170 SimdReal t1 = load<SimdReal>(blc + b);
1171 SimdReal t2 = load<SimdReal>(sol + b);
1172 SimdReal mvb = t1 * t2;
1173 store(blc_sol + b, mvb);
1174 store(mlambda + b, load<SimdReal>(mlambda + b) + mvb);
1177 for (b = b0; b < b1; b++)
1181 mvb = blc[b]*sol[b];
1185 #endif // GMX_SIMD_HAVE_REAL
1187 /* Update the coordinates */
1188 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1190 /* nit*ncons*(37+9*nrec) flops */
1194 /* Update the velocities */
1195 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1196 /* 16 ncons flops */
1199 if (nlocat != nullptr && (bCalcDHDL || bCalcVir))
1201 if (lincsd->bTaskDep)
1203 /* In lincs_update_atoms threads might cross-read mlambda */
1207 /* Only account for local atoms */
1208 for (b = b0; b < b1; b++)
1210 mlambda[b] *= 0.5*nlocat[b];
1219 for (b = b0; b < b1; b++)
1221 /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1222 * later after the contributions are reduced over the threads.
1224 dhdl -= lincsd->mlambda[b]*lincsd->ddist[b];
1226 lincsd->task[th].dhdlambda = dhdl;
1231 /* Constraint virial */
1232 for (b = b0; b < b1; b++)
1236 tmp0 = -bllen[b]*mlambda[b];
1237 for (i = 0; i < DIM; i++)
1239 tmp1 = tmp0*r[b][i];
1240 for (j = 0; j < DIM; j++)
1242 vir_r_m_dr[i][j] -= tmp1*r[b][j];
1245 } /* 22 ncons flops */
1249 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1250 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1252 * (26+nrec)*ncons + (6+2*nrec)*nrtot
1253 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1255 * (63+nrec)*ncons + (6+4*nrec)*nrtot
1259 /*! \brief Sets the elements in the LINCS matrix for task task. */
1260 static void set_lincs_matrix_task(Lincs *li,
1262 const real *invmass,
1264 int *nCrossTaskTriangles)
1268 /* Construct the coupling coefficient matrix blmf */
1269 li_task->ntriangle = 0;
1271 *nCrossTaskTriangles = 0;
1272 for (i = li_task->b0; i < li_task->b1; i++)
1277 a2 = li->bla[2*i+1];
1278 for (n = li->blnr[i]; (n < li->blnr[i+1]); n++)
1280 int k, sign, center, end;
1284 /* If we are using multiple, independent tasks for LINCS,
1285 * the calls to check_assign_connected should have
1286 * put all connected constraints in our task.
1288 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1290 if (a1 == li->bla[2*k] || a2 == li->bla[2*k+1])
1298 if (a1 == li->bla[2*k] || a1 == li->bla[2*k+1])
1308 li->blmf[n] = sign*invmass[center]*li->blc[i]*li->blc[k];
1309 li->blmf1[n] = sign*0.5;
1310 if (li->ncg_triangle > 0)
1314 /* Look for constraint triangles */
1315 for (nk = li->blnr[k]; (nk < li->blnr[k+1]); nk++)
1318 if (kk != i && kk != k &&
1319 (li->bla[2*kk] == end || li->bla[2*kk+1] == end))
1321 /* Check if the constraints in this triangle actually
1322 * belong to a different task. We still assign them
1323 * here, since it's convenient for the triangle
1324 * iterations, but we then need an extra barrier.
1326 if (k < li_task->b0 || k >= li_task->b1 ||
1327 kk < li_task->b0 || kk >= li_task->b1)
1329 (*nCrossTaskTriangles)++;
1332 if (li_task->ntriangle == 0 ||
1333 li_task->triangle[li_task->ntriangle - 1] < i)
1335 /* Add this constraint to the triangle list */
1336 li_task->triangle[li_task->ntriangle] = i;
1337 li_task->tri_bits[li_task->ntriangle] = 0;
1338 li_task->ntriangle++;
1339 if (li->blnr[i+1] - li->blnr[i] > static_cast<int>(sizeof(li_task->tri_bits[0])*8 - 1))
1341 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %zu allowed for constraints participating in triangles",
1342 li->blnr[i+1] - li->blnr[i],
1343 sizeof(li_task->tri_bits[0])*8-1);
1346 li_task->tri_bits[li_task->ntriangle-1] |= (1 << (n - li->blnr[i]));
1355 /*! \brief Sets the elements in the LINCS matrix. */
1356 static void set_lincs_matrix(Lincs *li, real *invmass, real lambda)
1359 const real invsqrt2 = 0.7071067811865475244;
1361 for (i = 0; (i < li->nc); i++)
1366 a2 = li->bla[2*i+1];
1367 li->blc[i] = gmx::invsqrt(invmass[a1] + invmass[a2]);
1368 li->blc1[i] = invsqrt2;
1371 /* Construct the coupling coefficient matrix blmf */
1372 int th, ntriangle = 0, ncc_triangle = 0, nCrossTaskTriangles = 0;
1373 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle, nCrossTaskTriangles) num_threads(li->ntask) schedule(static)
1374 for (th = 0; th < li->ntask; th++)
1378 set_lincs_matrix_task(li, &li->task[th], invmass,
1379 &ncc_triangle, &nCrossTaskTriangles);
1380 ntriangle += li->task[th].ntriangle;
1382 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1384 li->ntriangle = ntriangle;
1385 li->ncc_triangle = ncc_triangle;
1386 li->bTaskDepTri = (nCrossTaskTriangles > 0);
1390 fprintf(debug, "The %d constraints participate in %d triangles\n",
1391 li->nc, li->ntriangle);
1392 fprintf(debug, "There are %d constraint couplings, of which %d in triangles\n",
1393 li->ncc, li->ncc_triangle);
1394 if (li->ntriangle > 0 && li->ntask > 1)
1396 fprintf(debug, "%d constraint triangles contain constraints assigned to different tasks\n",
1397 nCrossTaskTriangles);
1402 * so we know with which lambda value the masses have been set.
1404 li->matlam = lambda;
1407 //! Finds all triangles of atoms that share constraints to a central atom.
1408 static int count_triangle_constraints(const InteractionLists &ilist,
1409 const t_blocka &at2con)
1411 int ncon1, ncon_tot;
1412 int c0, n1, c1, ac1, n2, c2;
1415 ncon1 = ilist[F_CONSTR].size()/3;
1416 ncon_tot = ncon1 + ilist[F_CONSTRNC].size()/3;
1418 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1419 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1422 for (c0 = 0; c0 < ncon_tot; c0++)
1424 bool bTriangle = FALSE;
1425 const int *iap = constr_iatomptr(ia1, ia2, c0);
1426 const int a00 = iap[1];
1427 const int a01 = iap[2];
1428 for (n1 = at2con.index[a01]; n1 < at2con.index[a01+1]; n1++)
1433 const int *iap = constr_iatomptr(ia1, ia2, c1);
1434 const int a10 = iap[1];
1435 const int a11 = iap[2];
1444 for (n2 = at2con.index[ac1]; n2 < at2con.index[ac1+1]; n2++)
1447 if (c2 != c0 && c2 != c1)
1449 const int *iap = constr_iatomptr(ia1, ia2, c2);
1450 const int a20 = iap[1];
1451 const int a21 = iap[2];
1452 if (a20 == a00 || a21 == a00)
1466 return ncon_triangle;
1469 //! Finds sequences of sequential constraints.
1470 static bool more_than_two_sequential_constraints(const InteractionLists &ilist,
1471 const t_blocka &at2con)
1473 int ncon1, ncon_tot, c;
1474 bool bMoreThanTwoSequentialConstraints;
1476 ncon1 = ilist[F_CONSTR].size()/3;
1477 ncon_tot = ncon1 + ilist[F_CONSTRNC].size()/3;
1479 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1480 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1482 bMoreThanTwoSequentialConstraints = FALSE;
1483 for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
1485 const int *iap = constr_iatomptr(ia1, ia2, c);
1486 const int a1 = iap[1];
1487 const int a2 = iap[2];
1488 /* Check if this constraint has constraints connected at both atoms */
1489 if (at2con.index[a1+1] - at2con.index[a1] > 1 &&
1490 at2con.index[a2+1] - at2con.index[a2] > 1)
1492 bMoreThanTwoSequentialConstraints = TRUE;
1496 return bMoreThanTwoSequentialConstraints;
1499 Lincs *init_lincs(FILE *fplog, const gmx_mtop_t &mtop,
1500 int nflexcon_global, ArrayRef<const t_blocka> at2con,
1501 bool bPLINCS, int nIter, int nProjOrder)
1503 // TODO this should become a unique_ptr
1505 bool bMoreThanTwoSeq;
1509 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
1510 bPLINCS ? " Parallel" : "");
1516 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1517 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1518 li->ncg_flex = nflexcon_global;
1521 li->nOrder = nProjOrder;
1523 li->max_connect = 0;
1524 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1526 for (int a = 0; a < mtop.moltype[mt].atoms.nr; a++)
1528 li->max_connect = std::max(li->max_connect,
1529 at2con[mt].index[a + 1] - at2con[mt].index[a]);
1533 li->ncg_triangle = 0;
1534 bMoreThanTwoSeq = FALSE;
1535 for (const gmx_molblock_t &molb : mtop.molblock)
1537 const gmx_moltype_t &molt = mtop.moltype[molb.type];
1541 count_triangle_constraints(molt.ilist, at2con[molb.type]);
1543 if (!bMoreThanTwoSeq &&
1544 more_than_two_sequential_constraints(molt.ilist, at2con[molb.type]))
1546 bMoreThanTwoSeq = TRUE;
1550 /* Check if we need to communicate not only before LINCS,
1551 * but also before each iteration.
1552 * The check for only two sequential constraints is only
1553 * useful for the common case of H-bond only constraints.
1554 * With more effort we could also make it useful for small
1555 * molecules with nr. sequential constraints <= nOrder-1.
1557 li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1559 if (debug && bPLINCS)
1561 fprintf(debug, "PLINCS communication before each iteration: %d\n",
1562 static_cast<int>(li->bCommIter));
1565 /* LINCS can run on any number of threads.
1566 * Currently the number is fixed for the whole simulation,
1567 * but it could be set in set_lincs().
1568 * The current constraint to task assignment code can create independent
1569 * tasks only when not more than two constraints are connected sequentially.
1571 li->ntask = gmx_omp_nthreads_get(emntLINCS);
1572 li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1575 fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n",
1576 li->ntask, li->bTaskDep ? "" : "in");
1584 /* Allocate an extra elements for "task-overlap" constraints */
1585 li->task.resize(li->ntask + 1);
1588 if (bPLINCS || li->ncg_triangle > 0)
1590 please_cite(fplog, "Hess2008a");
1594 please_cite(fplog, "Hess97a");
1599 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1602 fprintf(fplog, "There are constraints between atoms in different decomposition domains,\n"
1603 "will communicate selected coordinates each lincs iteration\n");
1605 if (li->ncg_triangle > 0)
1608 "%d constraints are involved in constraint triangles,\n"
1609 "will apply an additional matrix expansion of order %d for couplings\n"
1610 "between constraints inside triangles\n",
1611 li->ncg_triangle, li->nOrder);
1618 void done_lincs(Lincs *li)
1623 /*! \brief Sets up the work division over the threads. */
1624 static void lincs_thread_setup(Lincs *li, int natoms)
1631 if (natoms > li->atf_nalloc)
1633 li->atf_nalloc = over_alloc_large(natoms);
1634 srenew(li->atf, li->atf_nalloc);
1638 /* Clear the atom flags */
1639 for (a = 0; a < natoms; a++)
1641 bitmask_clear(&atf[a]);
1644 if (li->ntask > BITMASK_SIZE)
1646 gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1649 for (th = 0; th < li->ntask; th++)
1654 li_task = &li->task[th];
1656 /* For each atom set a flag for constraints from each */
1657 for (b = li_task->b0; b < li_task->b1; b++)
1659 bitmask_set_bit(&atf[li->bla[b*2 ]], th);
1660 bitmask_set_bit(&atf[li->bla[b*2 + 1]], th);
1664 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1665 for (th = 0; th < li->ntask; th++)
1673 li_task = &li->task[th];
1675 if (li_task->b1 - li_task->b0 > li_task->ind_nalloc)
1677 li_task->ind_nalloc = over_alloc_large(li_task->b1-li_task->b0);
1678 srenew(li_task->ind, li_task->ind_nalloc);
1679 srenew(li_task->ind_r, li_task->ind_nalloc);
1682 bitmask_init_low_bits(&mask, th);
1685 li_task->nind_r = 0;
1686 for (b = li_task->b0; b < li_task->b1; b++)
1688 /* We let the constraint with the lowest thread index
1689 * operate on atoms with constraints from multiple threads.
1691 if (bitmask_is_disjoint(atf[li->bla[b*2]], mask) &&
1692 bitmask_is_disjoint(atf[li->bla[b*2+1]], mask))
1694 /* Add the constraint to the local atom update index */
1695 li_task->ind[li_task->nind++] = b;
1699 /* Add the constraint to the rest block */
1700 li_task->ind_r[li_task->nind_r++] = b;
1704 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1707 /* We need to copy all constraints which have not be assigned
1708 * to a thread to a separate list which will be handled by one thread.
1710 li_m = &li->task[li->ntask];
1713 for (th = 0; th < li->ntask; th++)
1718 li_task = &li->task[th];
1720 if (li_m->nind + li_task->nind_r > li_m->ind_nalloc)
1722 li_m->ind_nalloc = over_alloc_large(li_m->nind+li_task->nind_r);
1723 srenew(li_m->ind, li_m->ind_nalloc);
1726 for (b = 0; b < li_task->nind_r; b++)
1728 li_m->ind[li_m->nind++] = li_task->ind_r[b];
1733 fprintf(debug, "LINCS thread %d: %d constraints\n",
1740 fprintf(debug, "LINCS thread r: %d constraints\n",
1745 /*! \brief There is no realloc with alignment, so here we make one for reals.
1746 * Note that this function does not preserve the contents of the memory.
1748 static void resize_real_aligned(real **ptr, int nelem)
1750 sfree_aligned(*ptr);
1751 snew_aligned(*ptr, nelem, align_bytes);
1754 //! Assign a constraint.
1755 static void assign_constraint(Lincs *li,
1756 int constraint_index,
1758 real lenA, real lenB,
1759 const t_blocka *at2con)
1765 /* Make an mapping of local topology constraint index to LINCS index */
1766 li->con_index[constraint_index] = con;
1768 li->bllen0[con] = lenA;
1769 li->ddist[con] = lenB - lenA;
1770 /* Set the length to the topology A length */
1771 li->bllen[con] = lenA;
1772 li->bla[2*con] = a1;
1773 li->bla[2*con+1] = a2;
1775 /* Make space in the constraint connection matrix for constraints
1776 * connected to both end of the current constraint.
1779 at2con->index[a1 + 1] - at2con->index[a1] - 1 +
1780 at2con->index[a2 + 1] - at2con->index[a2] - 1;
1782 li->blnr[con + 1] = li->ncc;
1784 /* Increase the constraint count */
1788 /*! \brief Check if constraint with topology index constraint_index is connected
1789 * to other constraints, and if so add those connected constraints to our task. */
1790 static void check_assign_connected(Lincs *li,
1791 const t_iatom *iatom,
1795 const t_blocka *at2con)
1797 /* Currently this function only supports constraint groups
1798 * in which all constraints share at least one atom
1799 * (e.g. H-bond constraints).
1800 * Check both ends of the current constraint for
1801 * connected constraints. We need to assign those
1806 for (end = 0; end < 2; end++)
1810 a = (end == 0 ? a1 : a2);
1812 for (k = at2con->index[a]; k < at2con->index[a + 1]; k++)
1817 /* Check if constraint cc has not yet been assigned */
1818 if (li->con_index[cc] == -1)
1824 lenA = idef.iparams[type].constr.dA;
1825 lenB = idef.iparams[type].constr.dB;
1827 if (bDynamics || lenA != 0 || lenB != 0)
1829 assign_constraint(li, cc, iatom[3*cc + 1], iatom[3*cc + 2], lenA, lenB, at2con);
1836 /*! \brief Check if constraint with topology index constraint_index is involved
1837 * in a constraint triangle, and if so add the other two constraints
1838 * in the triangle to our task. */
1839 static void check_assign_triangle(Lincs *li,
1840 const t_iatom *iatom,
1843 int constraint_index,
1845 const t_blocka *at2con)
1847 int nca, cc[32], ca[32], k;
1848 int c_triangle[2] = { -1, -1 };
1851 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1856 if (c != constraint_index)
1860 aa1 = iatom[c*3 + 1];
1861 aa2 = iatom[c*3 + 2];
1877 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1882 if (c != constraint_index)
1886 aa1 = iatom[c*3 + 1];
1887 aa2 = iatom[c*3 + 2];
1890 for (i = 0; i < nca; i++)
1894 c_triangle[0] = cc[i];
1901 for (i = 0; i < nca; i++)
1905 c_triangle[0] = cc[i];
1913 if (c_triangle[0] >= 0)
1917 for (end = 0; end < 2; end++)
1919 /* Check if constraint c_triangle[end] has not yet been assigned */
1920 if (li->con_index[c_triangle[end]] == -1)
1925 i = c_triangle[end]*3;
1927 lenA = idef.iparams[type].constr.dA;
1928 lenB = idef.iparams[type].constr.dB;
1930 if (bDynamics || lenA != 0 || lenB != 0)
1932 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1939 //! Sets matrix indices.
1940 static void set_matrix_indices(Lincs *li,
1941 const Task *li_task,
1942 const t_blocka *at2con,
1947 for (b = li_task->b0; b < li_task->b1; b++)
1952 a2 = li->bla[b*2 + 1];
1955 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1959 concon = li->con_index[at2con->a[k]];
1962 li->blbnb[i++] = concon;
1965 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1969 concon = li->con_index[at2con->a[k]];
1972 li->blbnb[i++] = concon;
1978 /* Order the blbnb matrix to optimize memory access */
1979 std::sort(&(li->blbnb[li->blnr[b]]), &(li->blbnb[li->blnr[b+1]]));
1984 void set_lincs(const t_idef &idef,
1985 const t_mdatoms &md,
1987 const t_commrec *cr,
1993 int i, ncc_alloc_old, ncon_tot;
1998 /* Zero the thread index ranges.
1999 * Otherwise without local constraints we could return with old ranges.
2001 for (i = 0; i < li->ntask; i++)
2005 li->task[i].nind = 0;
2009 li->task[li->ntask].nind = 0;
2012 /* This is the local topology, so there are only F_CONSTR constraints */
2013 if (idef.il[F_CONSTR].nr == 0)
2015 /* There are no constraints,
2016 * we do not need to fill any data structures.
2023 fprintf(debug, "Building the LINCS connectivity\n");
2026 if (DOMAINDECOMP(cr))
2028 if (cr->dd->constraints)
2032 dd_get_constraint_range(cr->dd, &start, &natoms);
2036 natoms = dd_numHomeAtoms(*cr->dd);
2044 at2con = make_at2con(natoms, idef.il, idef.iparams,
2045 flexibleConstraintTreatment(bDynamics));
2047 ncon_tot = idef.il[F_CONSTR].nr/3;
2049 /* Ensure we have enough padding for aligned loads for each thread */
2050 if (ncon_tot + li->ntask*simd_width > li->nc_alloc || li->nc_alloc == 0)
2052 int old_alloc = li->nc_alloc;
2053 li->nc_alloc = over_alloc_dd(ncon_tot + li->ntask*simd_width);
2054 srenew(li->con_index, li->nc_alloc);
2055 resize_real_aligned(&li->bllen0, li->nc_alloc);
2056 resize_real_aligned(&li->ddist, li->nc_alloc);
2057 srenew(li->bla, 2*li->nc_alloc);
2058 resize_real_aligned(&li->blc, li->nc_alloc);
2059 resize_real_aligned(&li->blc1, li->nc_alloc);
2060 srenew(li->blnr, li->nc_alloc + 1);
2061 resize_real_aligned(&li->bllen, li->nc_alloc);
2062 srenew(li->tmpv, li->nc_alloc);
2063 /* Need to clear the SIMD padding */
2064 for (int i = old_alloc; i < li->nc_alloc; i++)
2066 clear_rvec(li->tmpv[i]);
2068 if (DOMAINDECOMP(cr))
2070 srenew(li->nlocat, li->nc_alloc);
2072 resize_real_aligned(&li->tmp1, li->nc_alloc);
2073 resize_real_aligned(&li->tmp2, li->nc_alloc);
2074 resize_real_aligned(&li->tmp3, li->nc_alloc);
2075 resize_real_aligned(&li->tmp4, li->nc_alloc);
2076 resize_real_aligned(&li->mlambda, li->nc_alloc);
2079 iatom = idef.il[F_CONSTR].iatoms;
2081 ncc_alloc_old = li->ncc_alloc;
2082 li->blnr[0] = li->ncc;
2084 /* Assign the constraints for li->ntask LINCS tasks.
2085 * We target a uniform distribution of constraints over the tasks.
2086 * Note that when flexible constraints are present, but are removed here
2087 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
2088 * happen during normal MD, that's ok.
2090 int ncon_assign, ncon_target, con, th;
2092 /* Determine the number of constraints we need to assign here */
2093 ncon_assign = ncon_tot;
2096 /* With energy minimization, flexible constraints are ignored
2097 * (and thus minimized, as they should be).
2099 ncon_assign -= countFlexibleConstraints(idef.il, idef.iparams);
2102 /* Set the target constraint count per task to exactly uniform,
2103 * this might be overridden below.
2105 ncon_target = (ncon_assign + li->ntask - 1)/li->ntask;
2107 /* Mark all constraints as unassigned by setting their index to -1 */
2108 for (con = 0; con < ncon_tot; con++)
2110 li->con_index[con] = -1;
2114 for (th = 0; th < li->ntask; th++)
2118 li_task = &li->task[th];
2120 #if GMX_SIMD_HAVE_REAL
2121 /* With indepedent tasks we likely have H-bond constraints or constraint
2122 * pairs. The connected constraints will be pulled into the task, so the
2123 * constraints per task will often exceed ncon_target.
2124 * Triangle constraints can also increase the count, but there are
2125 * relatively few of those, so we usually expect to get ncon_target.
2129 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
2130 * since otherwise a lot of operations can be wasted.
2131 * There are several ways to round here, we choose the one
2132 * that alternates block sizes, which helps with Intel HT.
2134 ncon_target = ((ncon_assign*(th + 1))/li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1) & ~(GMX_SIMD_REAL_WIDTH - 1);
2136 #endif // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
2138 /* Continue filling the arrays where we left off with the previous task,
2139 * including padding for SIMD.
2141 li_task->b0 = li->nc;
2143 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2145 if (li->con_index[con] == -1)
2150 type = iatom[3*con];
2151 a1 = iatom[3*con + 1];
2152 a2 = iatom[3*con + 2];
2153 lenA = idef.iparams[type].constr.dA;
2154 lenB = idef.iparams[type].constr.dB;
2155 /* Skip the flexible constraints when not doing dynamics */
2156 if (bDynamics || lenA != 0 || lenB != 0)
2158 assign_constraint(li, con, a1, a2, lenA, lenB, &at2con);
2160 if (li->ntask > 1 && !li->bTaskDep)
2162 /* We can generate independent tasks. Check if we
2163 * need to assign connected constraints to our task.
2165 check_assign_connected(li, iatom, idef, bDynamics,
2168 if (li->ntask > 1 && li->ncg_triangle > 0)
2170 /* Ensure constraints in one triangle are assigned
2173 check_assign_triangle(li, iatom, idef, bDynamics,
2174 con, a1, a2, &at2con);
2182 li_task->b1 = li->nc;
2186 /* Copy the last atom pair indices and lengths for constraints
2187 * up to a multiple of simd_width, such that we can do all
2188 * SIMD operations without having to worry about end effects.
2192 li->nc = ((li_task->b1 + simd_width - 1)/simd_width)*simd_width;
2193 last = li_task->b1 - 1;
2194 for (i = li_task->b1; i < li->nc; i++)
2196 li->bla[i*2 ] = li->bla[last*2 ];
2197 li->bla[i*2 + 1] = li->bla[last*2 + 1];
2198 li->bllen0[i] = li->bllen0[last];
2199 li->ddist[i] = li->ddist[last];
2200 li->bllen[i] = li->bllen[last];
2201 li->blnr[i + 1] = li->blnr[last + 1];
2205 /* Keep track of how many constraints we assigned */
2206 li->nc_real += li_task->b1 - li_task->b0;
2210 fprintf(debug, "LINCS task %d constraints %d - %d\n",
2211 th, li_task->b0, li_task->b1);
2215 assert(li->nc_real == ncon_assign);
2219 /* Without DD we order the blbnb matrix to optimize memory access.
2220 * With DD the overhead of sorting is more than the gain during access.
2222 bSortMatrix = !DOMAINDECOMP(cr);
2224 if (li->ncc > li->ncc_alloc)
2226 li->ncc_alloc = over_alloc_small(li->ncc);
2227 srenew(li->blbnb, li->ncc_alloc);
2230 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2231 for (th = 0; th < li->ntask; th++)
2237 li_task = &li->task[th];
2239 if (li->ncg_triangle > 0 &&
2240 li_task->b1 - li_task->b0 > li_task->tri_alloc)
2242 /* This is allocating too much, but it is difficult to improve */
2243 li_task->tri_alloc = over_alloc_dd(li_task->b1 - li_task->b0);
2244 srenew(li_task->triangle, li_task->tri_alloc);
2245 srenew(li_task->tri_bits, li_task->tri_alloc);
2248 set_matrix_indices(li, li_task, &at2con, bSortMatrix);
2250 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2253 done_blocka(&at2con);
2255 if (cr->dd == nullptr)
2257 /* Since the matrix is static, we should free some memory */
2258 li->ncc_alloc = li->ncc;
2259 srenew(li->blbnb, li->ncc_alloc);
2262 if (li->ncc_alloc > ncc_alloc_old)
2264 srenew(li->blmf, li->ncc_alloc);
2265 srenew(li->blmf1, li->ncc_alloc);
2266 srenew(li->tmpncc, li->ncc_alloc);
2269 gmx::ArrayRef<const int> nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2270 if (!nlocat_dd.empty())
2272 /* Convert nlocat from local topology to LINCS constraint indexing */
2273 for (con = 0; con < ncon_tot; con++)
2275 li->nlocat[li->con_index[con]] = nlocat_dd[con];
2280 li->nlocat = nullptr;
2285 fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n",
2286 li->nc_real, li->nc, li->ncc);
2291 lincs_thread_setup(li, md.nr);
2294 set_lincs_matrix(li, md.invmass, md.lambda);
2297 //! Issues a warning when LINCS constraints cannot be satisfied.
2298 static void lincs_warning(gmx_domdec_t *dd, const rvec *x, rvec *xprime, t_pbc *pbc,
2299 int ncons, const int *bla, real *bllen, real wangle,
2300 int maxwarn, int *warncount)
2304 real wfac, d0, d1, cosine;
2306 wfac = std::cos(DEG2RAD*wangle);
2309 "bonds that rotated more than %g degrees:\n"
2310 " atom 1 atom 2 angle previous, current, constraint length\n",
2313 for (b = 0; b < ncons; b++)
2319 pbc_dx_aiuc(pbc, x[i], x[j], v0);
2320 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2324 rvec_sub(x[i], x[j], v0);
2325 rvec_sub(xprime[i], xprime[j], v1);
2329 cosine = ::iprod(v0, v1)/(d0*d1);
2333 " %6d %6d %5.1f %8.4f %8.4f %8.4f\n",
2334 ddglatnr(dd, i), ddglatnr(dd, j),
2335 RAD2DEG*std::acos(cosine), d0, d1, bllen[b]);
2336 if (!std::isfinite(d1))
2338 gmx_fatal(FARGS, "Bond length not finite.");
2344 if (*warncount > maxwarn)
2346 too_many_constraint_warnings(econtLINCS, *warncount);
2350 //! Determine how well the constraints have been satisfied.
2351 static void cconerr(const Lincs *lincsd,
2352 rvec *x, t_pbc *pbc,
2353 real *ncons_loc, real *ssd, real *max, int *imax)
2355 const int *bla, *nlocat;
2358 int count, im, task;
2361 bllen = lincsd->bllen;
2362 nlocat = lincsd->nlocat;
2368 for (task = 0; task < lincsd->ntask; task++)
2372 for (b = lincsd->task[task].b0; b < lincsd->task[task].b1; b++)
2379 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
2383 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
2386 len = r2*gmx::invsqrt(r2);
2387 d = std::abs(len/bllen[b]-1);
2388 if (d > ma && (nlocat == nullptr || nlocat[b]))
2393 if (nlocat == nullptr)
2400 ssd2 += nlocat[b]*d*d;
2406 *ncons_loc = (nlocat ? 0.5 : 1)*count;
2407 *ssd = (nlocat ? 0.5 : 1)*ssd2;
2412 bool constrain_lincs(bool computeRmsd,
2413 const t_inputrec &ir,
2415 Lincs *lincsd, const t_mdatoms &md,
2416 const t_commrec *cr,
2417 const gmx_multisim_t &ms,
2418 const rvec *x, rvec *xprime, rvec *min_proj,
2419 matrix box, t_pbc *pbc,
2420 real lambda, real *dvdlambda,
2421 real invdt, rvec *v,
2422 bool bCalcVir, tensor vir_r_m_dr,
2423 ConstraintVariable econq,
2425 int maxwarn, int *warncount)
2428 char buf2[22], buf3[STRLEN];
2430 real ncons_loc, p_ssd, p_max = 0;
2436 /* This boolean should be set by a flag passed to this routine.
2437 * We can also easily check if any constraint length is changed,
2438 * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2440 bCalcDHDL = (ir.efep != efepNO && dvdlambda != nullptr);
2442 if (lincsd->nc == 0 && cr->dd == nullptr)
2446 lincsd->rmsdData = {{0}};
2452 if (econq == ConstraintVariable::Positions)
2454 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2455 * also with efep!=fepNO.
2457 if (ir.efep != efepNO)
2459 if (md.nMassPerturbed && lincsd->matlam != md.lambda)
2461 set_lincs_matrix(lincsd, md.invmass, md.lambda);
2464 for (i = 0; i < lincsd->nc; i++)
2466 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
2470 if (lincsd->ncg_flex)
2472 /* Set the flexible constraint lengths to the old lengths */
2475 for (i = 0; i < lincsd->nc; i++)
2477 if (lincsd->bllen[i] == 0)
2479 pbc_dx_aiuc(pbc, x[lincsd->bla[2*i]], x[lincsd->bla[2*i+1]], dx);
2480 lincsd->bllen[i] = norm(dx);
2486 for (i = 0; i < lincsd->nc; i++)
2488 if (lincsd->bllen[i] == 0)
2491 std::sqrt(distance2(x[lincsd->bla[2*i]],
2492 x[lincsd->bla[2*i+1]]));
2500 cconerr(lincsd, xprime, pbc,
2501 &ncons_loc, &p_ssd, &p_max, &p_imax);
2504 /* This bWarn var can be updated by multiple threads
2505 * at the same time. But as we only need to detect
2506 * if a warning occurred or not, this is not an issue.
2510 /* The OpenMP parallel region of constrain_lincs for coords */
2511 #pragma omp parallel num_threads(lincsd->ntask)
2515 int th = gmx_omp_get_thread_num();
2517 clear_mat(lincsd->task[th].vir_r_m_dr);
2519 do_lincs(x, xprime, box, pbc, lincsd, th,
2522 ir.LincsWarnAngle, &bWarn,
2524 th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2526 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2529 if (debug && lincsd->nc > 0)
2531 fprintf(debug, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2532 fprintf(debug, " Before LINCS %.6f %.6f %6d %6d\n",
2533 std::sqrt(p_ssd/ncons_loc), p_max,
2534 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2535 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2537 if (computeRmsd || debug)
2539 cconerr(lincsd, xprime, pbc,
2540 &ncons_loc, &p_ssd, &p_max, &p_imax);
2541 lincsd->rmsdData[0] = ncons_loc;
2542 lincsd->rmsdData[1] = p_ssd;
2546 lincsd->rmsdData = {{0}};
2548 if (debug && lincsd->nc > 0)
2551 " After LINCS %.6f %.6f %6d %6d\n\n",
2552 std::sqrt(p_ssd/ncons_loc), p_max,
2553 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2554 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2559 if (maxwarn < INT_MAX)
2561 cconerr(lincsd, xprime, pbc,
2562 &ncons_loc, &p_ssd, &p_max, &p_imax);
2563 if (isMultiSim(&ms))
2565 sprintf(buf3, " in simulation %d", ms.sim);
2572 "\nStep %s, time %g (ps) LINCS WARNING%s\n"
2573 "relative constraint deviation after LINCS:\n"
2574 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2575 gmx_step_str(step, buf2), ir.init_t+step*ir.delta_t,
2577 std::sqrt(p_ssd/ncons_loc), p_max,
2578 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2579 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2581 lincs_warning(cr->dd, x, xprime, pbc,
2582 lincsd->nc, lincsd->bla, lincsd->bllen,
2583 ir.LincsWarnAngle, maxwarn, warncount);
2585 bOK = (p_max < 0.5);
2588 if (lincsd->ncg_flex)
2590 for (i = 0; (i < lincsd->nc); i++)
2592 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2594 lincsd->bllen[i] = 0;
2601 /* The OpenMP parallel region of constrain_lincs for derivatives */
2602 #pragma omp parallel num_threads(lincsd->ntask)
2606 int th = gmx_omp_get_thread_num();
2608 do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
2609 md.invmass, econq, bCalcDHDL,
2610 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2612 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2618 /* Reduce the dH/dlambda contributions over the threads */
2623 for (th = 0; th < lincsd->ntask; th++)
2625 dhdlambda += lincsd->task[th].dhdlambda;
2627 if (econq == ConstraintVariable::Positions)
2629 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2630 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2631 dhdlambda /= ir.delta_t*ir.delta_t;
2633 *dvdlambda += dhdlambda;
2636 if (bCalcVir && lincsd->ntask > 1)
2638 for (i = 1; i < lincsd->ntask; i++)
2640 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2644 /* count assuming nit=1 */
2645 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2646 inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
2647 if (lincsd->ntriangle > 0)
2649 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
2653 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real*2);
2657 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);