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 /*! \brief Calculate the constraint distance vectors r to project on from x.
505 * Determine the right-hand side of the matrix equation using quantity f.
506 * This function only differs from calc_dr_x_xp_simd below in that
507 * no constraint length is subtracted and no PBC is used for f. */
508 static void gmx_simdcall
509 calc_dr_x_f_simd(int b0,
512 const rvec * gmx_restrict x,
513 const rvec * gmx_restrict f,
514 const real * gmx_restrict blc,
515 const real * pbc_simd,
516 rvec * gmx_restrict r,
517 real * gmx_restrict rhs,
518 real * gmx_restrict sol)
520 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
522 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
524 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
529 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
531 SimdReal x0_S, y0_S, z0_S;
532 SimdReal x1_S, y1_S, z1_S;
533 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
534 SimdReal fx_S, fy_S, fz_S, ip_S, rhs_S;
535 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
536 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
538 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
540 offset0[i] = bla[bs*2 + i*2];
541 offset1[i] = bla[bs*2 + i*2 + 1];
544 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
545 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
550 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
552 n2_S = norm2(rx_S, ry_S, rz_S);
553 il_S = invsqrt(n2_S);
559 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
561 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(f), offset0, &x0_S, &y0_S, &z0_S);
562 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(f), offset1, &x1_S, &y1_S, &z1_S);
567 ip_S = iprod(rx_S, ry_S, rz_S, fx_S, fy_S, fz_S);
569 rhs_S = load<SimdReal>(blc + bs) * ip_S;
571 store(rhs + bs, rhs_S);
572 store(sol + bs, rhs_S);
575 #endif // GMX_SIMD_HAVE_REAL
577 /*! \brief LINCS projection, works on derivatives of the coordinates. */
578 static void do_lincsp(const rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
579 Lincs *lincsd, int th,
581 ConstraintVariable econq, bool bCalcDHDL,
582 bool bCalcVir, tensor rmdf)
585 int *bla, *blnr, *blbnb;
587 real *blc, *blmf, *blcc, *rhs1, *rhs2, *sol;
589 b0 = lincsd->task[th].b0;
590 b1 = lincsd->task[th].b1;
595 blbnb = lincsd->blbnb;
596 if (econq != ConstraintVariable::Force)
598 /* Use mass-weighted parameters */
604 /* Use non mass-weighted parameters */
606 blmf = lincsd->blmf1;
608 blcc = lincsd->tmpncc;
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, bla, x, f, blc,
630 #else // GMX_SIMD_HAVE_REAL
632 /* Compute normalized i-j vectors */
635 for (b = b0; b < b1; b++)
639 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
645 for (b = b0; b < b1; b++)
649 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
651 } /* 16 ncons flops */
654 for (b = b0; b < b1; b++)
661 mvb = blc[b]*(r[b][0]*(f[i][0] - f[j][0]) +
662 r[b][1]*(f[i][1] - f[j][1]) +
663 r[b][2]*(f[i][2] - f[j][2]));
669 #endif // GMX_SIMD_HAVE_REAL
671 if (lincsd->bTaskDep)
673 /* We need a barrier, since the matrix construction below
674 * can access entries in r of other threads.
679 /* Construct the (sparse) LINCS matrix */
680 for (b = b0; b < b1; b++)
684 for (n = blnr[b]; n < blnr[b+1]; n++)
686 blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
689 /* Together: 23*ncons + 6*nrtot flops */
691 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
692 /* nrec*(ncons+2*nrtot) flops */
694 if (econq == ConstraintVariable::Deriv_FlexCon)
696 /* We only want to constraint the flexible constraints,
697 * so we mask out the normal ones by setting sol to 0.
699 for (b = b0; b < b1; b++)
701 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
708 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
709 for (b = b0; b < b1; b++)
714 /* When constraining forces, we should not use mass weighting,
715 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
717 lincs_update_atoms(lincsd, th, 1.0, sol, r,
718 (econq != ConstraintVariable::Force) ? invmass : nullptr, fp);
725 for (b = b0; b < b1; b++)
727 dhdlambda -= sol[b]*lincsd->ddist[b];
730 lincsd->task[th].dhdlambda = dhdlambda;
735 /* Constraint virial,
736 * determines sum r_bond x delta f,
737 * where delta f is the constraint correction
738 * of the quantity that is being constrained.
740 for (b = b0; b < b1; b++)
745 mvb = lincsd->bllen[b]*sol[b];
746 for (i = 0; i < DIM; i++)
749 for (j = 0; j < DIM; j++)
751 rmdf[i][j] += tmp1*r[b][j];
754 } /* 23 ncons flops */
758 #if GMX_SIMD_HAVE_REAL
759 /*! \brief Calculate the constraint distance vectors r to project on from x.
761 * Determine the right-hand side of the matrix equation using coordinates xp. */
762 static void gmx_simdcall
763 calc_dr_x_xp_simd(int b0,
766 const rvec * gmx_restrict x,
767 const rvec * gmx_restrict xp,
768 const real * gmx_restrict bllen,
769 const real * gmx_restrict blc,
770 const real * pbc_simd,
771 rvec * gmx_restrict r,
772 real * gmx_restrict rhs,
773 real * gmx_restrict sol)
775 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
776 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
778 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
783 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
785 SimdReal x0_S, y0_S, z0_S;
786 SimdReal x1_S, y1_S, z1_S;
787 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
788 SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
789 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
790 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
792 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
794 offset0[i] = bla[bs*2 + i*2];
795 offset1[i] = bla[bs*2 + i*2 + 1];
798 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
799 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
804 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
806 n2_S = norm2(rx_S, ry_S, rz_S);
807 il_S = invsqrt(n2_S);
813 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
815 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(xp), offset0, &x0_S, &y0_S, &z0_S);
816 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(xp), offset1, &x1_S, &y1_S, &z1_S);
821 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
823 ip_S = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
825 rhs_S = load<SimdReal>(blc + bs) * (ip_S - load<SimdReal>(bllen + bs));
827 store(rhs + bs, rhs_S);
828 store(sol + bs, rhs_S);
831 #endif // GMX_SIMD_HAVE_REAL
833 /*! \brief Determine the distances and right-hand side for the next iteration. */
834 gmx_unused static void calc_dist_iter(
838 const rvec * gmx_restrict xp,
839 const real * gmx_restrict bllen,
840 const real * gmx_restrict blc,
843 real * gmx_restrict rhs,
844 real * gmx_restrict sol,
849 for (b = b0; b < b1; b++)
851 real len, len2, dlen2, mvb;
857 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
861 rvec_sub(xp[bla[2*b]], xp[bla[2*b+1]], dx);
864 dlen2 = 2*len2 - ::norm2(dx);
865 if (dlen2 < wfac*len2)
867 /* not race free - see detailed comment in caller */
872 mvb = blc[b]*(len - dlen2*gmx::invsqrt(dlen2));
880 } /* 20*ncons flops */
883 #if GMX_SIMD_HAVE_REAL
884 /*! \brief As calc_dist_iter(), but using SIMD intrinsics. */
885 static void gmx_simdcall
886 calc_dist_iter_simd(int b0,
889 const rvec * gmx_restrict x,
890 const real * gmx_restrict bllen,
891 const real * gmx_restrict blc,
892 const real * pbc_simd,
894 real * gmx_restrict rhs,
895 real * gmx_restrict sol,
898 SimdReal min_S(GMX_REAL_MIN);
900 SimdReal wfac_S(wfac);
905 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
907 /* Initialize all to FALSE */
908 warn_S = (two_S < setZero());
910 for (bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
912 SimdReal x0_S, y0_S, z0_S;
913 SimdReal x1_S, y1_S, z1_S;
914 SimdReal rx_S, ry_S, rz_S, n2_S;
915 SimdReal len_S, len2_S, dlen2_S, lc_S, blc_S;
916 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
917 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
919 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
921 offset0[i] = bla[bs*2 + i*2];
922 offset1[i] = bla[bs*2 + i*2 + 1];
925 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
926 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
931 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
933 n2_S = norm2(rx_S, ry_S, rz_S);
935 len_S = load<SimdReal>(bllen + bs);
936 len2_S = len_S * len_S;
938 dlen2_S = fms(two_S, len2_S, n2_S);
940 warn_S = warn_S || (dlen2_S < (wfac_S * len2_S));
942 /* Avoid 1/0 by taking the max with REAL_MIN.
943 * Note: when dlen2 is close to zero (90 degree constraint rotation),
944 * the accuracy of the algorithm is no longer relevant.
946 dlen2_S = max(dlen2_S, min_S);
948 lc_S = fnma(dlen2_S, invsqrt(dlen2_S), len_S);
950 blc_S = load<SimdReal>(blc + bs);
954 store(rhs + bs, lc_S);
955 store(sol + bs, lc_S);
963 #endif // GMX_SIMD_HAVE_REAL
965 //! Implements LINCS constraining.
966 static void do_lincs(const rvec *x, rvec *xp, matrix box, t_pbc *pbc,
967 Lincs *lincsd, int th,
971 real wangle, bool *bWarn,
972 real invdt, rvec * gmx_restrict v,
973 bool bCalcVir, tensor vir_r_m_dr)
975 int b0, b1, b, i, j, n, iter;
976 int *bla, *blnr, *blbnb;
978 real *blc, *blmf, *bllen, *blcc, *rhs1, *rhs2, *sol, *blc_sol, *mlambda;
981 b0 = lincsd->task[th].b0;
982 b1 = lincsd->task[th].b1;
987 blbnb = lincsd->blbnb;
990 bllen = lincsd->bllen;
991 blcc = lincsd->tmpncc;
995 blc_sol = lincsd->tmp4;
996 mlambda = lincsd->mlambda;
997 nlocat = lincsd->nlocat;
999 #if GMX_SIMD_HAVE_REAL
1001 /* This SIMD code does the same as the plain-C code after the #else.
1002 * The only difference is that we always call pbc code, as with SIMD
1003 * the overhead of pbc computation (when not needed) is small.
1005 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
1007 /* Convert the pbc struct for SIMD */
1008 set_pbc_simd(pbc, pbc_simd);
1010 /* Compute normalized x i-j vectors, store in r.
1011 * Compute the inner product of r and xp i-j and store in rhs1.
1013 calc_dr_x_xp_simd(b0, b1, bla, x, xp, bllen, blc,
1017 #else // GMX_SIMD_HAVE_REAL
1021 /* Compute normalized i-j vectors */
1022 for (b = b0; b < b1; b++)
1027 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
1030 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
1031 mvb = blc[b]*(::iprod(r[b], dx) - bllen[b]);
1038 /* Compute normalized i-j vectors */
1039 for (b = b0; b < b1; b++)
1041 real tmp0, tmp1, tmp2, rlen, mvb;
1045 tmp0 = x[i][0] - x[j][0];
1046 tmp1 = x[i][1] - x[j][1];
1047 tmp2 = x[i][2] - x[j][2];
1048 rlen = gmx::invsqrt(tmp0*tmp0 + tmp1*tmp1 + tmp2*tmp2);
1049 r[b][0] = rlen*tmp0;
1050 r[b][1] = rlen*tmp1;
1051 r[b][2] = rlen*tmp2;
1052 /* 16 ncons flops */
1056 mvb = blc[b]*(r[b][0]*(xp[i][0] - xp[j][0]) +
1057 r[b][1]*(xp[i][1] - xp[j][1]) +
1058 r[b][2]*(xp[i][2] - xp[j][2]) - bllen[b]);
1063 /* Together: 26*ncons + 6*nrtot flops */
1066 #endif // GMX_SIMD_HAVE_REAL
1068 if (lincsd->bTaskDep)
1070 /* We need a barrier, since the matrix construction below
1071 * can access entries in r of other threads.
1076 /* Construct the (sparse) LINCS matrix */
1077 for (b = b0; b < b1; b++)
1079 for (n = blnr[b]; n < blnr[b+1]; n++)
1081 blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
1084 /* Together: 26*ncons + 6*nrtot flops */
1086 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1087 /* nrec*(ncons+2*nrtot) flops */
1089 #if GMX_SIMD_HAVE_REAL
1090 for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1092 SimdReal t1 = load<SimdReal>(blc + b);
1093 SimdReal t2 = load<SimdReal>(sol + b);
1094 store(mlambda + b, t1 * t2);
1097 for (b = b0; b < b1; b++)
1099 mlambda[b] = blc[b]*sol[b];
1101 #endif // GMX_SIMD_HAVE_REAL
1103 /* Update the coordinates */
1104 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1107 ******** Correction for centripetal effects ********
1112 wfac = std::cos(DEG2RAD*wangle);
1115 for (iter = 0; iter < lincsd->nIter; iter++)
1117 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
1122 /* Communicate the corrected non-local coordinates */
1123 if (DOMAINDECOMP(cr))
1125 dd_move_x_constraints(cr->dd, box, xp, nullptr, FALSE);
1130 else if (lincsd->bTaskDep)
1135 #if GMX_SIMD_HAVE_REAL
1136 calc_dist_iter_simd(b0, b1, bla, xp, bllen, blc, pbc_simd, wfac,
1139 calc_dist_iter(b0, b1, bla, xp, bllen, blc, pbc, wfac,
1141 /* 20*ncons flops */
1142 #endif // GMX_SIMD_HAVE_REAL
1144 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1145 /* nrec*(ncons+2*nrtot) flops */
1147 #if GMX_SIMD_HAVE_REAL
1148 for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1150 SimdReal t1 = load<SimdReal>(blc + b);
1151 SimdReal t2 = load<SimdReal>(sol + b);
1152 SimdReal mvb = t1 * t2;
1153 store(blc_sol + b, mvb);
1154 store(mlambda + b, load<SimdReal>(mlambda + b) + mvb);
1157 for (b = b0; b < b1; b++)
1161 mvb = blc[b]*sol[b];
1165 #endif // GMX_SIMD_HAVE_REAL
1167 /* Update the coordinates */
1168 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1170 /* nit*ncons*(37+9*nrec) flops */
1174 /* Update the velocities */
1175 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1176 /* 16 ncons flops */
1179 if (nlocat != nullptr && (bCalcDHDL || bCalcVir))
1181 if (lincsd->bTaskDep)
1183 /* In lincs_update_atoms threads might cross-read mlambda */
1187 /* Only account for local atoms */
1188 for (b = b0; b < b1; b++)
1190 mlambda[b] *= 0.5*nlocat[b];
1199 for (b = b0; b < b1; b++)
1201 /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1202 * later after the contributions are reduced over the threads.
1204 dhdl -= lincsd->mlambda[b]*lincsd->ddist[b];
1206 lincsd->task[th].dhdlambda = dhdl;
1211 /* Constraint virial */
1212 for (b = b0; b < b1; b++)
1216 tmp0 = -bllen[b]*mlambda[b];
1217 for (i = 0; i < DIM; i++)
1219 tmp1 = tmp0*r[b][i];
1220 for (j = 0; j < DIM; j++)
1222 vir_r_m_dr[i][j] -= tmp1*r[b][j];
1225 } /* 22 ncons flops */
1229 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1230 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1232 * (26+nrec)*ncons + (6+2*nrec)*nrtot
1233 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1235 * (63+nrec)*ncons + (6+4*nrec)*nrtot
1239 /*! \brief Sets the elements in the LINCS matrix for task task. */
1240 static void set_lincs_matrix_task(Lincs *li,
1242 const real *invmass,
1244 int *nCrossTaskTriangles)
1248 /* Construct the coupling coefficient matrix blmf */
1249 li_task->ntriangle = 0;
1251 *nCrossTaskTriangles = 0;
1252 for (i = li_task->b0; i < li_task->b1; i++)
1257 a2 = li->bla[2*i+1];
1258 for (n = li->blnr[i]; (n < li->blnr[i+1]); n++)
1260 int k, sign, center, end;
1264 /* If we are using multiple, independent tasks for LINCS,
1265 * the calls to check_assign_connected should have
1266 * put all connected constraints in our task.
1268 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1270 if (a1 == li->bla[2*k] || a2 == li->bla[2*k+1])
1278 if (a1 == li->bla[2*k] || a1 == li->bla[2*k+1])
1288 li->blmf[n] = sign*invmass[center]*li->blc[i]*li->blc[k];
1289 li->blmf1[n] = sign*0.5;
1290 if (li->ncg_triangle > 0)
1294 /* Look for constraint triangles */
1295 for (nk = li->blnr[k]; (nk < li->blnr[k+1]); nk++)
1298 if (kk != i && kk != k &&
1299 (li->bla[2*kk] == end || li->bla[2*kk+1] == end))
1301 /* Check if the constraints in this triangle actually
1302 * belong to a different task. We still assign them
1303 * here, since it's convenient for the triangle
1304 * iterations, but we then need an extra barrier.
1306 if (k < li_task->b0 || k >= li_task->b1 ||
1307 kk < li_task->b0 || kk >= li_task->b1)
1309 (*nCrossTaskTriangles)++;
1312 if (li_task->ntriangle == 0 ||
1313 li_task->triangle[li_task->ntriangle - 1] < i)
1315 /* Add this constraint to the triangle list */
1316 li_task->triangle[li_task->ntriangle] = i;
1317 li_task->tri_bits[li_task->ntriangle] = 0;
1318 li_task->ntriangle++;
1319 if (li->blnr[i+1] - li->blnr[i] > static_cast<int>(sizeof(li_task->tri_bits[0])*8 - 1))
1321 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %zu allowed for constraints participating in triangles",
1322 li->blnr[i+1] - li->blnr[i],
1323 sizeof(li_task->tri_bits[0])*8-1);
1326 li_task->tri_bits[li_task->ntriangle-1] |= (1 << (n - li->blnr[i]));
1335 /*! \brief Sets the elements in the LINCS matrix. */
1336 static void set_lincs_matrix(Lincs *li, real *invmass, real lambda)
1339 const real invsqrt2 = 0.7071067811865475244;
1341 for (i = 0; (i < li->nc); i++)
1346 a2 = li->bla[2*i+1];
1347 li->blc[i] = gmx::invsqrt(invmass[a1] + invmass[a2]);
1348 li->blc1[i] = invsqrt2;
1351 /* Construct the coupling coefficient matrix blmf */
1352 int th, ntriangle = 0, ncc_triangle = 0, nCrossTaskTriangles = 0;
1353 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle, nCrossTaskTriangles) num_threads(li->ntask) schedule(static)
1354 for (th = 0; th < li->ntask; th++)
1358 set_lincs_matrix_task(li, &li->task[th], invmass,
1359 &ncc_triangle, &nCrossTaskTriangles);
1360 ntriangle = li->task[th].ntriangle;
1362 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1364 li->ntriangle = ntriangle;
1365 li->ncc_triangle = ncc_triangle;
1366 li->bTaskDepTri = (nCrossTaskTriangles > 0);
1370 fprintf(debug, "The %d constraints participate in %d triangles\n",
1371 li->nc, li->ntriangle);
1372 fprintf(debug, "There are %d constraint couplings, of which %d in triangles\n",
1373 li->ncc, li->ncc_triangle);
1374 if (li->ntriangle > 0 && li->ntask > 1)
1376 fprintf(debug, "%d constraint triangles contain constraints assigned to different tasks\n",
1377 nCrossTaskTriangles);
1382 * so we know with which lambda value the masses have been set.
1384 li->matlam = lambda;
1387 //! Finds all triangles of atoms that share constraints to a central atom.
1388 static int count_triangle_constraints(const InteractionList *ilist,
1389 const t_blocka &at2con)
1391 int ncon1, ncon_tot;
1392 int c0, n1, c1, ac1, n2, c2;
1395 ncon1 = ilist[F_CONSTR].size()/3;
1396 ncon_tot = ncon1 + ilist[F_CONSTRNC].size()/3;
1398 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1399 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1402 for (c0 = 0; c0 < ncon_tot; c0++)
1404 bool bTriangle = FALSE;
1405 const int *iap = constr_iatomptr(ia1, ia2, c0);
1406 const int a00 = iap[1];
1407 const int a01 = iap[2];
1408 for (n1 = at2con.index[a01]; n1 < at2con.index[a01+1]; n1++)
1413 const int *iap = constr_iatomptr(ia1, ia2, c1);
1414 const int a10 = iap[1];
1415 const int a11 = iap[2];
1424 for (n2 = at2con.index[ac1]; n2 < at2con.index[ac1+1]; n2++)
1427 if (c2 != c0 && c2 != c1)
1429 const int *iap = constr_iatomptr(ia1, ia2, c2);
1430 const int a20 = iap[1];
1431 const int a21 = iap[2];
1432 if (a20 == a00 || a21 == a00)
1446 return ncon_triangle;
1449 //! Finds sequences of sequential constraints.
1450 static bool more_than_two_sequential_constraints(const InteractionList *ilist,
1451 const t_blocka &at2con)
1453 int ncon1, ncon_tot, c;
1454 bool bMoreThanTwoSequentialConstraints;
1456 ncon1 = ilist[F_CONSTR].size()/3;
1457 ncon_tot = ncon1 + ilist[F_CONSTRNC].size()/3;
1459 gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
1460 gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
1462 bMoreThanTwoSequentialConstraints = FALSE;
1463 for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
1465 const int *iap = constr_iatomptr(ia1, ia2, c);
1466 const int a1 = iap[1];
1467 const int a2 = iap[2];
1468 /* Check if this constraint has constraints connected at both atoms */
1469 if (at2con.index[a1+1] - at2con.index[a1] > 1 &&
1470 at2con.index[a2+1] - at2con.index[a2] > 1)
1472 bMoreThanTwoSequentialConstraints = TRUE;
1476 return bMoreThanTwoSequentialConstraints;
1479 Lincs *init_lincs(FILE *fplog, const gmx_mtop_t &mtop,
1480 int nflexcon_global, ArrayRef<const t_blocka> at2con,
1481 bool bPLINCS, int nIter, int nProjOrder)
1483 // TODO this should become a unique_ptr
1485 bool bMoreThanTwoSeq;
1489 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
1490 bPLINCS ? " Parallel" : "");
1496 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1497 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1498 li->ncg_flex = nflexcon_global;
1501 li->nOrder = nProjOrder;
1503 li->max_connect = 0;
1504 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1506 for (int a = 0; a < mtop.moltype[mt].atoms.nr; a++)
1508 li->max_connect = std::max(li->max_connect,
1509 at2con[mt].index[a + 1] - at2con[mt].index[a]);
1513 li->ncg_triangle = 0;
1514 bMoreThanTwoSeq = FALSE;
1515 for (const gmx_molblock_t &molb : mtop.molblock)
1517 const gmx_moltype_t &molt = mtop.moltype[molb.type];
1521 count_triangle_constraints(molt.ilist, at2con[molb.type]);
1523 if (!bMoreThanTwoSeq &&
1524 more_than_two_sequential_constraints(molt.ilist, at2con[molb.type]))
1526 bMoreThanTwoSeq = TRUE;
1530 /* Check if we need to communicate not only before LINCS,
1531 * but also before each iteration.
1532 * The check for only two sequential constraints is only
1533 * useful for the common case of H-bond only constraints.
1534 * With more effort we could also make it useful for small
1535 * molecules with nr. sequential constraints <= nOrder-1.
1537 li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1539 if (debug && bPLINCS)
1541 fprintf(debug, "PLINCS communication before each iteration: %d\n",
1542 static_cast<int>(li->bCommIter));
1545 /* LINCS can run on any number of threads.
1546 * Currently the number is fixed for the whole simulation,
1547 * but it could be set in set_lincs().
1548 * The current constraint to task assignment code can create independent
1549 * tasks only when not more than two constraints are connected sequentially.
1551 li->ntask = gmx_omp_nthreads_get(emntLINCS);
1552 li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1555 fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n",
1556 li->ntask, li->bTaskDep ? "" : "in");
1564 /* Allocate an extra elements for "task-overlap" constraints */
1565 li->task.resize(li->ntask + 1);
1568 if (bPLINCS || li->ncg_triangle > 0)
1570 please_cite(fplog, "Hess2008a");
1574 please_cite(fplog, "Hess97a");
1579 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1582 fprintf(fplog, "There are inter charge-group constraints,\n"
1583 "will communicate selected coordinates each lincs iteration\n");
1585 if (li->ncg_triangle > 0)
1588 "%d constraints are involved in constraint triangles,\n"
1589 "will apply an additional matrix expansion of order %d for couplings\n"
1590 "between constraints inside triangles\n",
1591 li->ncg_triangle, li->nOrder);
1598 void done_lincs(Lincs *li)
1603 /*! \brief Sets up the work division over the threads. */
1604 static void lincs_thread_setup(Lincs *li, int natoms)
1611 if (natoms > li->atf_nalloc)
1613 li->atf_nalloc = over_alloc_large(natoms);
1614 srenew(li->atf, li->atf_nalloc);
1618 /* Clear the atom flags */
1619 for (a = 0; a < natoms; a++)
1621 bitmask_clear(&atf[a]);
1624 if (li->ntask > BITMASK_SIZE)
1626 gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1629 for (th = 0; th < li->ntask; th++)
1634 li_task = &li->task[th];
1636 /* For each atom set a flag for constraints from each */
1637 for (b = li_task->b0; b < li_task->b1; b++)
1639 bitmask_set_bit(&atf[li->bla[b*2 ]], th);
1640 bitmask_set_bit(&atf[li->bla[b*2 + 1]], th);
1644 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1645 for (th = 0; th < li->ntask; th++)
1653 li_task = &li->task[th];
1655 if (li_task->b1 - li_task->b0 > li_task->ind_nalloc)
1657 li_task->ind_nalloc = over_alloc_large(li_task->b1-li_task->b0);
1658 srenew(li_task->ind, li_task->ind_nalloc);
1659 srenew(li_task->ind_r, li_task->ind_nalloc);
1662 bitmask_init_low_bits(&mask, th);
1665 li_task->nind_r = 0;
1666 for (b = li_task->b0; b < li_task->b1; b++)
1668 /* We let the constraint with the lowest thread index
1669 * operate on atoms with constraints from multiple threads.
1671 if (bitmask_is_disjoint(atf[li->bla[b*2]], mask) &&
1672 bitmask_is_disjoint(atf[li->bla[b*2+1]], mask))
1674 /* Add the constraint to the local atom update index */
1675 li_task->ind[li_task->nind++] = b;
1679 /* Add the constraint to the rest block */
1680 li_task->ind_r[li_task->nind_r++] = b;
1684 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1687 /* We need to copy all constraints which have not be assigned
1688 * to a thread to a separate list which will be handled by one thread.
1690 li_m = &li->task[li->ntask];
1693 for (th = 0; th < li->ntask; th++)
1698 li_task = &li->task[th];
1700 if (li_m->nind + li_task->nind_r > li_m->ind_nalloc)
1702 li_m->ind_nalloc = over_alloc_large(li_m->nind+li_task->nind_r);
1703 srenew(li_m->ind, li_m->ind_nalloc);
1706 for (b = 0; b < li_task->nind_r; b++)
1708 li_m->ind[li_m->nind++] = li_task->ind_r[b];
1713 fprintf(debug, "LINCS thread %d: %d constraints\n",
1720 fprintf(debug, "LINCS thread r: %d constraints\n",
1725 /*! \brief There is no realloc with alignment, so here we make one for reals.
1726 * Note that this function does not preserve the contents of the memory.
1728 static void resize_real_aligned(real **ptr, int nelem)
1730 sfree_aligned(*ptr);
1731 snew_aligned(*ptr, nelem, align_bytes);
1734 //! Assign a constraint.
1735 static void assign_constraint(Lincs *li,
1736 int constraint_index,
1738 real lenA, real lenB,
1739 const t_blocka *at2con)
1745 /* Make an mapping of local topology constraint index to LINCS index */
1746 li->con_index[constraint_index] = con;
1748 li->bllen0[con] = lenA;
1749 li->ddist[con] = lenB - lenA;
1750 /* Set the length to the topology A length */
1751 li->bllen[con] = lenA;
1752 li->bla[2*con] = a1;
1753 li->bla[2*con+1] = a2;
1755 /* Make space in the constraint connection matrix for constraints
1756 * connected to both end of the current constraint.
1759 at2con->index[a1 + 1] - at2con->index[a1] - 1 +
1760 at2con->index[a2 + 1] - at2con->index[a2] - 1;
1762 li->blnr[con + 1] = li->ncc;
1764 /* Increase the constraint count */
1768 /*! \brief Check if constraint with topology index constraint_index is connected
1769 * to other constraints, and if so add those connected constraints to our task. */
1770 static void check_assign_connected(Lincs *li,
1771 const t_iatom *iatom,
1775 const t_blocka *at2con)
1777 /* Currently this function only supports constraint groups
1778 * in which all constraints share at least one atom
1779 * (e.g. H-bond constraints).
1780 * Check both ends of the current constraint for
1781 * connected constraints. We need to assign those
1786 for (end = 0; end < 2; end++)
1790 a = (end == 0 ? a1 : a2);
1792 for (k = at2con->index[a]; k < at2con->index[a + 1]; k++)
1797 /* Check if constraint cc has not yet been assigned */
1798 if (li->con_index[cc] == -1)
1804 lenA = idef.iparams[type].constr.dA;
1805 lenB = idef.iparams[type].constr.dB;
1807 if (bDynamics || lenA != 0 || lenB != 0)
1809 assign_constraint(li, cc, iatom[3*cc + 1], iatom[3*cc + 2], lenA, lenB, at2con);
1816 /*! \brief Check if constraint with topology index constraint_index is involved
1817 * in a constraint triangle, and if so add the other two constraints
1818 * in the triangle to our task. */
1819 static void check_assign_triangle(Lincs *li,
1820 const t_iatom *iatom,
1823 int constraint_index,
1825 const t_blocka *at2con)
1827 int nca, cc[32], ca[32], k;
1828 int c_triangle[2] = { -1, -1 };
1831 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1836 if (c != constraint_index)
1840 aa1 = iatom[c*3 + 1];
1841 aa2 = iatom[c*3 + 2];
1857 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1862 if (c != constraint_index)
1866 aa1 = iatom[c*3 + 1];
1867 aa2 = iatom[c*3 + 2];
1870 for (i = 0; i < nca; i++)
1874 c_triangle[0] = cc[i];
1881 for (i = 0; i < nca; i++)
1885 c_triangle[0] = cc[i];
1893 if (c_triangle[0] >= 0)
1897 for (end = 0; end < 2; end++)
1899 /* Check if constraint c_triangle[end] has not yet been assigned */
1900 if (li->con_index[c_triangle[end]] == -1)
1905 i = c_triangle[end]*3;
1907 lenA = idef.iparams[type].constr.dA;
1908 lenB = idef.iparams[type].constr.dB;
1910 if (bDynamics || lenA != 0 || lenB != 0)
1912 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1919 //! Sets matrix indices.
1920 static void set_matrix_indices(Lincs *li,
1921 const Task *li_task,
1922 const t_blocka *at2con,
1927 for (b = li_task->b0; b < li_task->b1; b++)
1932 a2 = li->bla[b*2 + 1];
1935 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1939 concon = li->con_index[at2con->a[k]];
1942 li->blbnb[i++] = concon;
1945 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1949 concon = li->con_index[at2con->a[k]];
1952 li->blbnb[i++] = concon;
1958 /* Order the blbnb matrix to optimize memory access */
1959 std::sort(&(li->blbnb[li->blnr[b]]), &(li->blbnb[li->blnr[b+1]]));
1964 void set_lincs(const t_idef &idef,
1965 const t_mdatoms &md,
1967 const t_commrec *cr,
1973 int i, ncc_alloc_old, ncon_tot;
1978 /* Zero the thread index ranges.
1979 * Otherwise without local constraints we could return with old ranges.
1981 for (i = 0; i < li->ntask; i++)
1985 li->task[i].nind = 0;
1989 li->task[li->ntask].nind = 0;
1992 /* This is the local topology, so there are only F_CONSTR constraints */
1993 if (idef.il[F_CONSTR].nr == 0)
1995 /* There are no constraints,
1996 * we do not need to fill any data structures.
2003 fprintf(debug, "Building the LINCS connectivity\n");
2006 if (DOMAINDECOMP(cr))
2008 if (cr->dd->constraints)
2012 dd_get_constraint_range(cr->dd, &start, &natoms);
2016 natoms = dd_numHomeAtoms(*cr->dd);
2024 at2con = make_at2con(natoms, idef.il, idef.iparams,
2025 flexibleConstraintTreatment(bDynamics));
2027 ncon_tot = idef.il[F_CONSTR].nr/3;
2029 /* Ensure we have enough padding for aligned loads for each thread */
2030 if (ncon_tot + li->ntask*simd_width > li->nc_alloc || li->nc_alloc == 0)
2032 li->nc_alloc = over_alloc_dd(ncon_tot + li->ntask*simd_width);
2033 srenew(li->con_index, li->nc_alloc);
2034 resize_real_aligned(&li->bllen0, li->nc_alloc);
2035 resize_real_aligned(&li->ddist, li->nc_alloc);
2036 srenew(li->bla, 2*li->nc_alloc);
2037 resize_real_aligned(&li->blc, li->nc_alloc);
2038 resize_real_aligned(&li->blc1, li->nc_alloc);
2039 srenew(li->blnr, li->nc_alloc + 1);
2040 resize_real_aligned(&li->bllen, li->nc_alloc);
2041 srenew(li->tmpv, li->nc_alloc);
2042 if (DOMAINDECOMP(cr))
2044 srenew(li->nlocat, li->nc_alloc);
2046 resize_real_aligned(&li->tmp1, li->nc_alloc);
2047 resize_real_aligned(&li->tmp2, li->nc_alloc);
2048 resize_real_aligned(&li->tmp3, li->nc_alloc);
2049 resize_real_aligned(&li->tmp4, li->nc_alloc);
2050 resize_real_aligned(&li->mlambda, li->nc_alloc);
2053 iatom = idef.il[F_CONSTR].iatoms;
2055 ncc_alloc_old = li->ncc_alloc;
2056 li->blnr[0] = li->ncc;
2058 /* Assign the constraints for li->ntask LINCS tasks.
2059 * We target a uniform distribution of constraints over the tasks.
2060 * Note that when flexible constraints are present, but are removed here
2061 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
2062 * happen during normal MD, that's ok.
2064 int ncon_assign, ncon_target, con, th;
2066 /* Determine the number of constraints we need to assign here */
2067 ncon_assign = ncon_tot;
2070 /* With energy minimization, flexible constraints are ignored
2071 * (and thus minimized, as they should be).
2073 ncon_assign -= countFlexibleConstraints(idef.il, idef.iparams);
2076 /* Set the target constraint count per task to exactly uniform,
2077 * this might be overridden below.
2079 ncon_target = (ncon_assign + li->ntask - 1)/li->ntask;
2081 /* Mark all constraints as unassigned by setting their index to -1 */
2082 for (con = 0; con < ncon_tot; con++)
2084 li->con_index[con] = -1;
2088 for (th = 0; th < li->ntask; th++)
2092 li_task = &li->task[th];
2094 #if GMX_SIMD_HAVE_REAL
2095 /* With indepedent tasks we likely have H-bond constraints or constraint
2096 * pairs. The connected constraints will be pulled into the task, so the
2097 * constraints per task will often exceed ncon_target.
2098 * Triangle constraints can also increase the count, but there are
2099 * relatively few of those, so we usually expect to get ncon_target.
2103 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
2104 * since otherwise a lot of operations can be wasted.
2105 * There are several ways to round here, we choose the one
2106 * that alternates block sizes, which helps with Intel HT.
2108 ncon_target = ((ncon_assign*(th + 1))/li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1) & ~(GMX_SIMD_REAL_WIDTH - 1);
2110 #endif // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
2112 /* Continue filling the arrays where we left off with the previous task,
2113 * including padding for SIMD.
2115 li_task->b0 = li->nc;
2117 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2119 if (li->con_index[con] == -1)
2124 type = iatom[3*con];
2125 a1 = iatom[3*con + 1];
2126 a2 = iatom[3*con + 2];
2127 lenA = idef.iparams[type].constr.dA;
2128 lenB = idef.iparams[type].constr.dB;
2129 /* Skip the flexible constraints when not doing dynamics */
2130 if (bDynamics || lenA != 0 || lenB != 0)
2132 assign_constraint(li, con, a1, a2, lenA, lenB, &at2con);
2134 if (li->ntask > 1 && !li->bTaskDep)
2136 /* We can generate independent tasks. Check if we
2137 * need to assign connected constraints to our task.
2139 check_assign_connected(li, iatom, idef, bDynamics,
2142 if (li->ntask > 1 && li->ncg_triangle > 0)
2144 /* Ensure constraints in one triangle are assigned
2147 check_assign_triangle(li, iatom, idef, bDynamics,
2148 con, a1, a2, &at2con);
2156 li_task->b1 = li->nc;
2160 /* Copy the last atom pair indices and lengths for constraints
2161 * up to a multiple of simd_width, such that we can do all
2162 * SIMD operations without having to worry about end effects.
2166 li->nc = ((li_task->b1 + simd_width - 1)/simd_width)*simd_width;
2167 last = li_task->b1 - 1;
2168 for (i = li_task->b1; i < li->nc; i++)
2170 li->bla[i*2 ] = li->bla[last*2 ];
2171 li->bla[i*2 + 1] = li->bla[last*2 + 1];
2172 li->bllen0[i] = li->bllen0[last];
2173 li->ddist[i] = li->ddist[last];
2174 li->bllen[i] = li->bllen[last];
2175 li->blnr[i + 1] = li->blnr[last + 1];
2179 /* Keep track of how many constraints we assigned */
2180 li->nc_real += li_task->b1 - li_task->b0;
2184 fprintf(debug, "LINCS task %d constraints %d - %d\n",
2185 th, li_task->b0, li_task->b1);
2189 assert(li->nc_real == ncon_assign);
2193 /* Without DD we order the blbnb matrix to optimize memory access.
2194 * With DD the overhead of sorting is more than the gain during access.
2196 bSortMatrix = !DOMAINDECOMP(cr);
2198 if (li->ncc > li->ncc_alloc)
2200 li->ncc_alloc = over_alloc_small(li->ncc);
2201 srenew(li->blbnb, li->ncc_alloc);
2204 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2205 for (th = 0; th < li->ntask; th++)
2211 li_task = &li->task[th];
2213 if (li->ncg_triangle > 0 &&
2214 li_task->b1 - li_task->b0 > li_task->tri_alloc)
2216 /* This is allocating too much, but it is difficult to improve */
2217 li_task->tri_alloc = over_alloc_dd(li_task->b1 - li_task->b0);
2218 srenew(li_task->triangle, li_task->tri_alloc);
2219 srenew(li_task->tri_bits, li_task->tri_alloc);
2222 set_matrix_indices(li, li_task, &at2con, bSortMatrix);
2224 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2227 done_blocka(&at2con);
2229 if (cr->dd == nullptr)
2231 /* Since the matrix is static, we should free some memory */
2232 li->ncc_alloc = li->ncc;
2233 srenew(li->blbnb, li->ncc_alloc);
2236 if (li->ncc_alloc > ncc_alloc_old)
2238 srenew(li->blmf, li->ncc_alloc);
2239 srenew(li->blmf1, li->ncc_alloc);
2240 srenew(li->tmpncc, li->ncc_alloc);
2243 gmx::ArrayRef<const int> nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2244 if (!nlocat_dd.empty())
2246 /* Convert nlocat from local topology to LINCS constraint indexing */
2247 for (con = 0; con < ncon_tot; con++)
2249 li->nlocat[li->con_index[con]] = nlocat_dd[con];
2254 li->nlocat = nullptr;
2259 fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n",
2260 li->nc_real, li->nc, li->ncc);
2265 lincs_thread_setup(li, md.nr);
2268 set_lincs_matrix(li, md.invmass, md.lambda);
2271 //! Issues a warning when LINCS constraints cannot be satisfied.
2272 static void lincs_warning(gmx_domdec_t *dd, const rvec *x, rvec *xprime, t_pbc *pbc,
2273 int ncons, const int *bla, real *bllen, real wangle,
2274 int maxwarn, int *warncount)
2278 real wfac, d0, d1, cosine;
2280 wfac = std::cos(DEG2RAD*wangle);
2283 "bonds that rotated more than %g degrees:\n"
2284 " atom 1 atom 2 angle previous, current, constraint length\n",
2287 for (b = 0; b < ncons; b++)
2293 pbc_dx_aiuc(pbc, x[i], x[j], v0);
2294 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2298 rvec_sub(x[i], x[j], v0);
2299 rvec_sub(xprime[i], xprime[j], v1);
2303 cosine = ::iprod(v0, v1)/(d0*d1);
2307 " %6d %6d %5.1f %8.4f %8.4f %8.4f\n",
2308 ddglatnr(dd, i), ddglatnr(dd, j),
2309 RAD2DEG*std::acos(cosine), d0, d1, bllen[b]);
2310 if (!std::isfinite(d1))
2312 gmx_fatal(FARGS, "Bond length not finite.");
2318 if (*warncount > maxwarn)
2320 too_many_constraint_warnings(econtLINCS, *warncount);
2324 //! Determine how well the constraints have been satisfied.
2325 static void cconerr(const Lincs *lincsd,
2326 rvec *x, t_pbc *pbc,
2327 real *ncons_loc, real *ssd, real *max, int *imax)
2329 const int *bla, *nlocat;
2332 int count, im, task;
2335 bllen = lincsd->bllen;
2336 nlocat = lincsd->nlocat;
2342 for (task = 0; task < lincsd->ntask; task++)
2346 for (b = lincsd->task[task].b0; b < lincsd->task[task].b1; b++)
2353 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
2357 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
2360 len = r2*gmx::invsqrt(r2);
2361 d = std::abs(len/bllen[b]-1);
2362 if (d > ma && (nlocat == nullptr || nlocat[b]))
2367 if (nlocat == nullptr)
2374 ssd2 += nlocat[b]*d*d;
2380 *ncons_loc = (nlocat ? 0.5 : 1)*count;
2381 *ssd = (nlocat ? 0.5 : 1)*ssd2;
2386 bool constrain_lincs(bool computeRmsd,
2387 const t_inputrec &ir,
2389 Lincs *lincsd, const t_mdatoms &md,
2390 const t_commrec *cr,
2391 const gmx_multisim_t &ms,
2392 const rvec *x, rvec *xprime, rvec *min_proj,
2393 matrix box, t_pbc *pbc,
2394 real lambda, real *dvdlambda,
2395 real invdt, rvec *v,
2396 bool bCalcVir, tensor vir_r_m_dr,
2397 ConstraintVariable econq,
2399 int maxwarn, int *warncount)
2402 char buf2[22], buf3[STRLEN];
2404 real ncons_loc, p_ssd, p_max = 0;
2410 /* This boolean should be set by a flag passed to this routine.
2411 * We can also easily check if any constraint length is changed,
2412 * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2414 bCalcDHDL = (ir.efep != efepNO && dvdlambda != nullptr);
2416 if (lincsd->nc == 0 && cr->dd == nullptr)
2420 lincsd->rmsdData = {{0}};
2426 if (econq == ConstraintVariable::Positions)
2428 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2429 * also with efep!=fepNO.
2431 if (ir.efep != efepNO)
2433 if (md.nMassPerturbed && lincsd->matlam != md.lambda)
2435 set_lincs_matrix(lincsd, md.invmass, md.lambda);
2438 for (i = 0; i < lincsd->nc; i++)
2440 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
2444 if (lincsd->ncg_flex)
2446 /* Set the flexible constraint lengths to the old lengths */
2449 for (i = 0; i < lincsd->nc; i++)
2451 if (lincsd->bllen[i] == 0)
2453 pbc_dx_aiuc(pbc, x[lincsd->bla[2*i]], x[lincsd->bla[2*i+1]], dx);
2454 lincsd->bllen[i] = norm(dx);
2460 for (i = 0; i < lincsd->nc; i++)
2462 if (lincsd->bllen[i] == 0)
2465 std::sqrt(distance2(x[lincsd->bla[2*i]],
2466 x[lincsd->bla[2*i+1]]));
2474 cconerr(lincsd, xprime, pbc,
2475 &ncons_loc, &p_ssd, &p_max, &p_imax);
2478 /* This bWarn var can be updated by multiple threads
2479 * at the same time. But as we only need to detect
2480 * if a warning occurred or not, this is not an issue.
2484 /* The OpenMP parallel region of constrain_lincs for coords */
2485 #pragma omp parallel num_threads(lincsd->ntask)
2489 int th = gmx_omp_get_thread_num();
2491 clear_mat(lincsd->task[th].vir_r_m_dr);
2493 do_lincs(x, xprime, box, pbc, lincsd, th,
2496 ir.LincsWarnAngle, &bWarn,
2498 th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2500 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2503 if (debug && lincsd->nc > 0)
2505 fprintf(debug, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2506 fprintf(debug, " Before LINCS %.6f %.6f %6d %6d\n",
2507 std::sqrt(p_ssd/ncons_loc), p_max,
2508 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2509 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2511 if (computeRmsd || debug)
2513 cconerr(lincsd, xprime, pbc,
2514 &ncons_loc, &p_ssd, &p_max, &p_imax);
2515 lincsd->rmsdData[0] = ncons_loc;
2516 lincsd->rmsdData[1] = p_ssd;
2520 lincsd->rmsdData = {{0}};
2522 if (debug && lincsd->nc > 0)
2525 " After LINCS %.6f %.6f %6d %6d\n\n",
2526 std::sqrt(p_ssd/ncons_loc), p_max,
2527 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2528 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2533 if (maxwarn < INT_MAX)
2535 cconerr(lincsd, xprime, pbc,
2536 &ncons_loc, &p_ssd, &p_max, &p_imax);
2537 if (isMultiSim(&ms))
2539 sprintf(buf3, " in simulation %d", ms.sim);
2546 "\nStep %s, time %g (ps) LINCS WARNING%s\n"
2547 "relative constraint deviation after LINCS:\n"
2548 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2549 gmx_step_str(step, buf2), ir.init_t+step*ir.delta_t,
2551 std::sqrt(p_ssd/ncons_loc), p_max,
2552 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2553 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2555 lincs_warning(cr->dd, x, xprime, pbc,
2556 lincsd->nc, lincsd->bla, lincsd->bllen,
2557 ir.LincsWarnAngle, maxwarn, warncount);
2559 bOK = (p_max < 0.5);
2562 if (lincsd->ncg_flex)
2564 for (i = 0; (i < lincsd->nc); i++)
2566 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2568 lincsd->bllen[i] = 0;
2575 /* The OpenMP parallel region of constrain_lincs for derivatives */
2576 #pragma omp parallel num_threads(lincsd->ntask)
2580 int th = gmx_omp_get_thread_num();
2582 do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
2583 md.invmass, econq, bCalcDHDL,
2584 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2586 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2592 /* Reduce the dH/dlambda contributions over the threads */
2597 for (th = 0; th < lincsd->ntask; th++)
2599 dhdlambda += lincsd->task[th].dhdlambda;
2601 if (econq == ConstraintVariable::Positions)
2603 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2604 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2605 dhdlambda /= ir.delta_t*ir.delta_t;
2607 *dvdlambda += dhdlambda;
2610 if (bCalcVir && lincsd->ntask > 1)
2612 for (i = 1; i < lincsd->ntask; i++)
2614 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2618 /* count assuming nit=1 */
2619 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2620 inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
2621 if (lincsd->ntriangle > 0)
2623 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
2627 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real*2);
2631 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);