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
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/domdec/domdec_struct.h"
60 #include "gromacs/gmxlib/nrnb.h"
61 #include "gromacs/math/functions.h"
62 #include "gromacs/math/units.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/mdlib/constr.h"
65 #include "gromacs/mdlib/gmx_omp_nthreads.h"
66 #include "gromacs/mdlib/mdrun.h"
67 #include "gromacs/mdtypes/commrec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/mdatom.h"
71 #include "gromacs/pbcutil/pbc.h"
72 #include "gromacs/pbcutil/pbc-simd.h"
73 #include "gromacs/simd/simd.h"
74 #include "gromacs/simd/simd_math.h"
75 #include "gromacs/simd/vector_operations.h"
76 #include "gromacs/topology/block.h"
77 #include "gromacs/topology/mtop_util.h"
78 #include "gromacs/utility/arrayref.h"
79 #include "gromacs/utility/basedefinitions.h"
80 #include "gromacs/utility/bitmask.h"
81 #include "gromacs/utility/cstringutil.h"
82 #include "gromacs/utility/exceptions.h"
83 #include "gromacs/utility/fatalerror.h"
84 #include "gromacs/utility/gmxomp.h"
85 #include "gromacs/utility/pleasecite.h"
86 #include "gromacs/utility/smalloc.h"
88 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
93 //! Unit of work within LINCS.
96 //! First constraint for this task.
98 //! b1-1 is the last constraint for this task.
100 //! The number of constraints in triangles.
102 //! The list of triangle constraints.
103 int *triangle = nullptr;
104 //! The bits tell if the matrix element should be used.
105 int *tri_bits = nullptr;
106 //! Allocation size of triangle and tri_bits.
108 //! Number of indices.
110 //! Constraint index for updating atom data.
112 //! Number of indices.
114 //! Constraint index for updating atom data.
115 int *ind_r = nullptr;
116 //! Allocation size of ind and ind_r.
118 //! Temporary variable for virial calculation.
119 tensor vir_r_m_dr = {{0}};
120 //! Temporary variable for lambda derivative.
124 /*! \brief Data for LINCS algorithm.
129 //! The global number of constraints.
131 //! The global number of flexible constraints.
133 //! The global number of constraints in triangles.
134 int ncg_triangle = 0;
135 //! The number of iterations.
137 //! The order of the matrix expansion.
139 //! The maximum number of constraints connected to a single atom.
142 //! The number of real constraints.
144 //! The number of constraints including padding for SIMD.
146 //! The number we allocated memory for.
148 //! The number of constraint connections.
150 //! The number we allocated memory for.
152 //! The FE lambda value used for filling blc and blmf.
154 //! mapping from topology to LINCS constraints.
155 int *con_index = nullptr;
156 //! The reference distance in topology A.
157 real *bllen0 = nullptr;
158 //! The reference distance in top B - the r.d. in top A.
159 real *ddist = nullptr;
160 //! The atom pairs involved in the constraints.
162 //! 1/sqrt(invmass1 invmass2).
164 //! As blc, but with all masses 1.
165 real *blc1 = nullptr;
166 //! Index into blbnb and blmf.
168 //! List of constraint connections.
169 int *blbnb = nullptr;
170 //! The local number of constraints in triangles.
172 //! The number of constraint connections in triangles.
173 int ncc_triangle = 0;
174 //! Communicate before each LINCS interation.
175 bool bCommIter = false;
176 //! Matrix of mass factors for constraint connections.
177 real *blmf = nullptr;
178 //! As blmf, but with all masses 1.
179 real *blmf1 = nullptr;
180 //! The reference bond length.
181 real *bllen = nullptr;
182 //! The local atom count per constraint, can be NULL.
183 int *nlocat = nullptr;
185 /*! \brief The number of tasks used for LINCS work.
187 * \todo This is mostly used to loop over \c task, which would
188 * be nicer to do with range-based for loops, but the thread
189 * index is used for constructing bit masks and organizing the
190 * virial output buffer, so other things need to change,
193 /*! \brief LINCS thread division */
194 std::vector<Task> task;
195 //! Atom flags for thread parallelization.
196 gmx_bitmask_t *atf = nullptr;
197 //! Allocation size of atf
199 //! Are the LINCS tasks interdependent?
200 bool bTaskDep = false;
201 //! Are there triangle constraints that cross task borders?
202 bool bTaskDepTri = false;
203 //! Arrays for temporary storage in the LINCS algorithm.
205 rvec *tmpv = nullptr;
206 real *tmpncc = nullptr;
207 real *tmp1 = nullptr;
208 real *tmp2 = nullptr;
209 real *tmp3 = nullptr;
210 real *tmp4 = nullptr;
212 //! The Lagrange multipliers times -1.
213 real *mlambda = nullptr;
214 //! Storage for the constraint RMS relative deviation output.
215 std::array<real, 2> rmsdData = {{0}};
218 /*! \brief Define simd_width for memory allocation used for SIMD code */
219 #if GMX_SIMD_HAVE_REAL
220 static const int simd_width = GMX_SIMD_REAL_WIDTH;
222 static const int simd_width = 1;
225 /*! \brief Align to 128 bytes, consistent with the current implementation of
226 AlignedAllocator, which currently forces 128 byte alignment. */
227 static const int align_bytes = 128;
229 ArrayRef<real> lincs_rmsdData(Lincs *lincsd)
231 return lincsd->rmsdData;
234 real lincs_rmsd(const Lincs *lincsd)
236 if (lincsd->rmsdData[0] > 0)
238 return std::sqrt(lincsd->rmsdData[1]/lincsd->rmsdData[0]);
246 /*! \brief Do a set of nrec LINCS matrix multiplications.
248 * This function will return with up to date thread-local
249 * constraint data, without an OpenMP barrier.
251 static void lincs_matrix_expand(const Lincs *lincsd,
254 real *rhs1, real *rhs2, real *sol)
256 int b0, b1, nrec, rec;
257 const int *blnr = lincsd->blnr;
258 const int *blbnb = lincsd->blbnb;
262 nrec = lincsd->nOrder;
264 for (rec = 0; rec < nrec; rec++)
268 if (lincsd->bTaskDep)
272 for (b = b0; b < b1; b++)
278 for (n = blnr[b]; n < blnr[b+1]; n++)
280 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
283 sol[b] = sol[b] + mvb;
291 } /* nrec*(ncons+2*nrtot) flops */
293 if (lincsd->ntriangle > 0)
295 /* Perform an extra nrec recursions for only the constraints
296 * involved in rigid triangles.
297 * In this way their accuracy should come close to those of the other
298 * constraints, since traingles of constraints can produce eigenvalues
299 * around 0.7, while the effective eigenvalue for bond constraints
300 * is around 0.4 (and 0.7*0.7=0.5).
303 if (lincsd->bTaskDep)
305 /* We need a barrier here, since other threads might still be
306 * reading the contents of rhs1 and/o rhs2.
307 * We could avoid this barrier by introducing two extra rhs
308 * arrays for the triangle constraints only.
313 /* Constraints involved in a triangle are ensured to be in the same
314 * LINCS task. This means no barriers are required during the extra
315 * iterations for the triangle constraints.
317 const int *triangle = li_task->triangle;
318 const int *tri_bits = li_task->tri_bits;
320 for (rec = 0; rec < nrec; rec++)
324 for (tb = 0; tb < li_task->ntriangle; tb++)
326 int b, bits, nr0, nr1, n;
334 for (n = nr0; n < nr1; n++)
336 if (bits & (1 << (n - nr0)))
338 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
342 sol[b] = sol[b] + mvb;
350 } /* nrec*(ntriangle + ncc_triangle*2) flops */
352 if (lincsd->bTaskDepTri)
354 /* The constraints triangles are decoupled from each other,
355 * but constraints in one triangle cross thread task borders.
356 * We could probably avoid this with more advanced setup code.
363 //! Update atomic coordinates when an index is not required.
364 static void lincs_update_atoms_noind(int ncons, const int *bla,
366 const real *fac, rvec *r,
371 real mvb, im1, im2, tmp0, tmp1, tmp2;
373 if (invmass != nullptr)
375 for (b = 0; b < ncons; b++)
391 } /* 16 ncons flops */
395 for (b = 0; b < ncons; b++)
413 //! Update atomic coordinates when an index is required.
414 static void lincs_update_atoms_ind(int ncons, const int *ind, const int *bla,
416 const real *fac, rvec *r,
421 real mvb, im1, im2, tmp0, tmp1, tmp2;
423 if (invmass != nullptr)
425 for (bi = 0; bi < ncons; bi++)
442 } /* 16 ncons flops */
446 for (bi = 0; bi < ncons; bi++)
461 } /* 16 ncons flops */
465 //! Update coordinates for atoms.
466 static void lincs_update_atoms(Lincs *li, int th,
468 const real *fac, rvec *r,
474 /* Single thread, we simply update for all constraints */
475 lincs_update_atoms_noind(li->nc_real,
476 li->bla, prefac, fac, r, invmass, x);
480 /* Update the atom vector components for our thread local
481 * constraints that only access our local atom range.
482 * This can be done without a barrier.
484 lincs_update_atoms_ind(li->task[th].nind, li->task[th].ind,
485 li->bla, prefac, fac, r, invmass, x);
487 if (li->task[li->ntask].nind > 0)
489 /* Update the constraints that operate on atoms
490 * in multiple thread atom blocks on the master thread.
495 lincs_update_atoms_ind(li->task[li->ntask].nind,
496 li->task[li->ntask].ind,
497 li->bla, prefac, fac, r, invmass, x);
503 #if GMX_SIMD_HAVE_REAL
504 /*! \brief Calculate the constraint distance vectors r to project on from x.
506 * Determine the right-hand side of the matrix equation using quantity f.
507 * This function only differs from calc_dr_x_xp_simd below in that
508 * no constraint length is subtracted and no PBC is used for f. */
509 static void gmx_simdcall
510 calc_dr_x_f_simd(int b0,
513 const rvec * gmx_restrict x,
514 const rvec * gmx_restrict f,
515 const real * gmx_restrict blc,
516 const real * pbc_simd,
517 rvec * gmx_restrict r,
518 real * gmx_restrict rhs,
519 real * gmx_restrict sol)
521 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
523 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
525 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
530 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
532 SimdReal x0_S, y0_S, z0_S;
533 SimdReal x1_S, y1_S, z1_S;
534 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
535 SimdReal fx_S, fy_S, fz_S, ip_S, rhs_S;
536 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
537 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
539 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
541 offset0[i] = bla[bs*2 + i*2];
542 offset1[i] = bla[bs*2 + i*2 + 1];
545 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
546 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
551 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
553 n2_S = norm2(rx_S, ry_S, rz_S);
554 il_S = invsqrt(n2_S);
560 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
562 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(f), offset0, &x0_S, &y0_S, &z0_S);
563 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(f), offset1, &x1_S, &y1_S, &z1_S);
568 ip_S = iprod(rx_S, ry_S, rz_S, fx_S, fy_S, fz_S);
570 rhs_S = load<SimdReal>(blc + bs) * ip_S;
572 store(rhs + bs, rhs_S);
573 store(sol + bs, rhs_S);
576 #endif // GMX_SIMD_HAVE_REAL
578 /*! \brief LINCS projection, works on derivatives of the coordinates. */
579 static void do_lincsp(const rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
580 Lincs *lincsd, int th,
582 ConstraintVariable econq, bool bCalcDHDL,
583 bool bCalcVir, tensor rmdf)
586 int *bla, *blnr, *blbnb;
588 real *blc, *blmf, *blcc, *rhs1, *rhs2, *sol;
590 b0 = lincsd->task[th].b0;
591 b1 = lincsd->task[th].b1;
596 blbnb = lincsd->blbnb;
597 if (econq != ConstraintVariable::Force)
599 /* Use mass-weighted parameters */
605 /* Use non mass-weighted parameters */
607 blmf = lincsd->blmf1;
609 blcc = lincsd->tmpncc;
614 #if GMX_SIMD_HAVE_REAL
615 /* This SIMD code does the same as the plain-C code after the #else.
616 * The only difference is that we always call pbc code, as with SIMD
617 * the overhead of pbc computation (when not needed) is small.
619 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
621 /* Convert the pbc struct for SIMD */
622 set_pbc_simd(pbc, pbc_simd);
624 /* Compute normalized x i-j vectors, store in r.
625 * Compute the inner product of r and xp i-j and store in rhs1.
627 calc_dr_x_f_simd(b0, b1, bla, x, f, blc,
631 #else // GMX_SIMD_HAVE_REAL
633 /* Compute normalized i-j vectors */
636 for (b = b0; b < b1; b++)
640 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
646 for (b = b0; b < b1; b++)
650 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
652 } /* 16 ncons flops */
655 for (b = b0; b < b1; b++)
662 mvb = blc[b]*(r[b][0]*(f[i][0] - f[j][0]) +
663 r[b][1]*(f[i][1] - f[j][1]) +
664 r[b][2]*(f[i][2] - f[j][2]));
670 #endif // GMX_SIMD_HAVE_REAL
672 if (lincsd->bTaskDep)
674 /* We need a barrier, since the matrix construction below
675 * can access entries in r of other threads.
680 /* Construct the (sparse) LINCS matrix */
681 for (b = b0; b < b1; b++)
685 for (n = blnr[b]; n < blnr[b+1]; n++)
687 blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
690 /* Together: 23*ncons + 6*nrtot flops */
692 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
693 /* nrec*(ncons+2*nrtot) flops */
695 if (econq == ConstraintVariable::Deriv_FlexCon)
697 /* We only want to constraint the flexible constraints,
698 * so we mask out the normal ones by setting sol to 0.
700 for (b = b0; b < b1; b++)
702 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
709 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
710 for (b = b0; b < b1; b++)
715 /* When constraining forces, we should not use mass weighting,
716 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
718 lincs_update_atoms(lincsd, th, 1.0, sol, r,
719 (econq != ConstraintVariable::Force) ? invmass : nullptr, fp);
726 for (b = b0; b < b1; b++)
728 dhdlambda -= sol[b]*lincsd->ddist[b];
731 lincsd->task[th].dhdlambda = dhdlambda;
736 /* Constraint virial,
737 * determines sum r_bond x delta f,
738 * where delta f is the constraint correction
739 * of the quantity that is being constrained.
741 for (b = b0; b < b1; b++)
746 mvb = lincsd->bllen[b]*sol[b];
747 for (i = 0; i < DIM; i++)
750 for (j = 0; j < DIM; j++)
752 rmdf[i][j] += tmp1*r[b][j];
755 } /* 23 ncons flops */
759 #if GMX_SIMD_HAVE_REAL
760 /*! \brief Calculate the constraint distance vectors r to project on from x.
762 * Determine the right-hand side of the matrix equation using coordinates xp. */
763 static void gmx_simdcall
764 calc_dr_x_xp_simd(int b0,
767 const rvec * gmx_restrict x,
768 const rvec * gmx_restrict xp,
769 const real * gmx_restrict bllen,
770 const real * gmx_restrict blc,
771 const real * pbc_simd,
772 rvec * gmx_restrict r,
773 real * gmx_restrict rhs,
774 real * gmx_restrict sol)
776 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
777 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset2[GMX_SIMD_REAL_WIDTH];
779 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
784 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
786 SimdReal x0_S, y0_S, z0_S;
787 SimdReal x1_S, y1_S, z1_S;
788 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
789 SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
790 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
791 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
793 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
795 offset0[i] = bla[bs*2 + i*2];
796 offset1[i] = bla[bs*2 + i*2 + 1];
799 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
800 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
805 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
807 n2_S = norm2(rx_S, ry_S, rz_S);
808 il_S = invsqrt(n2_S);
814 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
816 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(xp), offset0, &x0_S, &y0_S, &z0_S);
817 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(xp), offset1, &x1_S, &y1_S, &z1_S);
822 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
824 ip_S = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
826 rhs_S = load<SimdReal>(blc + bs) * (ip_S - load<SimdReal>(bllen + bs));
828 store(rhs + bs, rhs_S);
829 store(sol + bs, rhs_S);
832 #endif // GMX_SIMD_HAVE_REAL
834 /*! \brief Determine the distances and right-hand side for the next iteration. */
835 gmx_unused static void calc_dist_iter(
839 const rvec * gmx_restrict xp,
840 const real * gmx_restrict bllen,
841 const real * gmx_restrict blc,
844 real * gmx_restrict rhs,
845 real * gmx_restrict sol,
850 for (b = b0; b < b1; b++)
852 real len, len2, dlen2, mvb;
858 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
862 rvec_sub(xp[bla[2*b]], xp[bla[2*b+1]], dx);
865 dlen2 = 2*len2 - ::norm2(dx);
866 if (dlen2 < wfac*len2)
868 /* not race free - see detailed comment in caller */
873 mvb = blc[b]*(len - dlen2*gmx::invsqrt(dlen2));
881 } /* 20*ncons flops */
884 #if GMX_SIMD_HAVE_REAL
885 /*! \brief As calc_dist_iter(), but using SIMD intrinsics. */
886 static void gmx_simdcall
887 calc_dist_iter_simd(int b0,
890 const rvec * gmx_restrict x,
891 const real * gmx_restrict bllen,
892 const real * gmx_restrict blc,
893 const real * pbc_simd,
895 real * gmx_restrict rhs,
896 real * gmx_restrict sol,
899 SimdReal min_S(GMX_REAL_MIN);
901 SimdReal wfac_S(wfac);
906 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
908 /* Initialize all to FALSE */
909 warn_S = (two_S < setZero());
911 for (bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
913 SimdReal x0_S, y0_S, z0_S;
914 SimdReal x1_S, y1_S, z1_S;
915 SimdReal rx_S, ry_S, rz_S, n2_S;
916 SimdReal len_S, len2_S, dlen2_S, lc_S, blc_S;
917 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset0[GMX_SIMD_REAL_WIDTH];
918 alignas(GMX_SIMD_ALIGNMENT) std::int32_t offset1[GMX_SIMD_REAL_WIDTH];
920 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
922 offset0[i] = bla[bs*2 + i*2];
923 offset1[i] = bla[bs*2 + i*2 + 1];
926 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
927 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
932 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
934 n2_S = norm2(rx_S, ry_S, rz_S);
936 len_S = load<SimdReal>(bllen + bs);
937 len2_S = len_S * len_S;
939 dlen2_S = fms(two_S, len2_S, n2_S);
941 warn_S = warn_S || (dlen2_S < (wfac_S * len2_S));
943 /* Avoid 1/0 by taking the max with REAL_MIN.
944 * Note: when dlen2 is close to zero (90 degree constraint rotation),
945 * the accuracy of the algorithm is no longer relevant.
947 dlen2_S = max(dlen2_S, min_S);
949 lc_S = fnma(dlen2_S, invsqrt(dlen2_S), len_S);
951 blc_S = load<SimdReal>(blc + bs);
955 store(rhs + bs, lc_S);
956 store(sol + bs, lc_S);
964 #endif // GMX_SIMD_HAVE_REAL
966 //! Implements LINCS constraining.
967 static void do_lincs(const rvec *x, rvec *xp, matrix box, t_pbc *pbc,
968 Lincs *lincsd, int th,
972 real wangle, bool *bWarn,
973 real invdt, rvec * gmx_restrict v,
974 bool bCalcVir, tensor vir_r_m_dr)
976 int b0, b1, b, i, j, n, iter;
977 int *bla, *blnr, *blbnb;
979 real *blc, *blmf, *bllen, *blcc, *rhs1, *rhs2, *sol, *blc_sol, *mlambda;
982 b0 = lincsd->task[th].b0;
983 b1 = lincsd->task[th].b1;
988 blbnb = lincsd->blbnb;
991 bllen = lincsd->bllen;
992 blcc = lincsd->tmpncc;
996 blc_sol = lincsd->tmp4;
997 mlambda = lincsd->mlambda;
998 nlocat = lincsd->nlocat;
1000 #if GMX_SIMD_HAVE_REAL
1002 /* This SIMD code does the same as the plain-C code after the #else.
1003 * The only difference is that we always call pbc code, as with SIMD
1004 * the overhead of pbc computation (when not needed) is small.
1006 alignas(GMX_SIMD_ALIGNMENT) real pbc_simd[9*GMX_SIMD_REAL_WIDTH];
1008 /* Convert the pbc struct for SIMD */
1009 set_pbc_simd(pbc, pbc_simd);
1011 /* Compute normalized x i-j vectors, store in r.
1012 * Compute the inner product of r and xp i-j and store in rhs1.
1014 calc_dr_x_xp_simd(b0, b1, bla, x, xp, bllen, blc,
1018 #else // GMX_SIMD_HAVE_REAL
1022 /* Compute normalized i-j vectors */
1023 for (b = b0; b < b1; b++)
1028 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
1031 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
1032 mvb = blc[b]*(::iprod(r[b], dx) - bllen[b]);
1039 /* Compute normalized i-j vectors */
1040 for (b = b0; b < b1; b++)
1042 real tmp0, tmp1, tmp2, rlen, mvb;
1046 tmp0 = x[i][0] - x[j][0];
1047 tmp1 = x[i][1] - x[j][1];
1048 tmp2 = x[i][2] - x[j][2];
1049 rlen = gmx::invsqrt(tmp0*tmp0 + tmp1*tmp1 + tmp2*tmp2);
1050 r[b][0] = rlen*tmp0;
1051 r[b][1] = rlen*tmp1;
1052 r[b][2] = rlen*tmp2;
1053 /* 16 ncons flops */
1057 mvb = blc[b]*(r[b][0]*(xp[i][0] - xp[j][0]) +
1058 r[b][1]*(xp[i][1] - xp[j][1]) +
1059 r[b][2]*(xp[i][2] - xp[j][2]) - bllen[b]);
1064 /* Together: 26*ncons + 6*nrtot flops */
1067 #endif // GMX_SIMD_HAVE_REAL
1069 if (lincsd->bTaskDep)
1071 /* We need a barrier, since the matrix construction below
1072 * can access entries in r of other threads.
1077 /* Construct the (sparse) LINCS matrix */
1078 for (b = b0; b < b1; b++)
1080 for (n = blnr[b]; n < blnr[b+1]; n++)
1082 blcc[n] = blmf[n]*::iprod(r[b], r[blbnb[n]]);
1085 /* Together: 26*ncons + 6*nrtot flops */
1087 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1088 /* nrec*(ncons+2*nrtot) flops */
1090 #if GMX_SIMD_HAVE_REAL
1091 for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1093 SimdReal t1 = load<SimdReal>(blc + b);
1094 SimdReal t2 = load<SimdReal>(sol + b);
1095 store(mlambda + b, t1 * t2);
1098 for (b = b0; b < b1; b++)
1100 mlambda[b] = blc[b]*sol[b];
1102 #endif // GMX_SIMD_HAVE_REAL
1104 /* Update the coordinates */
1105 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1108 ******** Correction for centripetal effects ********
1113 wfac = std::cos(DEG2RAD*wangle);
1116 for (iter = 0; iter < lincsd->nIter; iter++)
1118 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
1123 /* Communicate the corrected non-local coordinates */
1124 if (DOMAINDECOMP(cr))
1126 dd_move_x_constraints(cr->dd, box, xp, nullptr, FALSE);
1131 else if (lincsd->bTaskDep)
1136 #if GMX_SIMD_HAVE_REAL
1137 calc_dist_iter_simd(b0, b1, bla, xp, bllen, blc, pbc_simd, wfac,
1140 calc_dist_iter(b0, b1, bla, xp, bllen, blc, pbc, wfac,
1142 /* 20*ncons flops */
1143 #endif // GMX_SIMD_HAVE_REAL
1145 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1146 /* nrec*(ncons+2*nrtot) flops */
1148 #if GMX_SIMD_HAVE_REAL
1149 for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1151 SimdReal t1 = load<SimdReal>(blc + b);
1152 SimdReal t2 = load<SimdReal>(sol + b);
1153 SimdReal mvb = t1 * t2;
1154 store(blc_sol + b, mvb);
1155 store(mlambda + b, load<SimdReal>(mlambda + b) + mvb);
1158 for (b = b0; b < b1; b++)
1162 mvb = blc[b]*sol[b];
1166 #endif // GMX_SIMD_HAVE_REAL
1168 /* Update the coordinates */
1169 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1171 /* nit*ncons*(37+9*nrec) flops */
1175 /* Update the velocities */
1176 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1177 /* 16 ncons flops */
1180 if (nlocat != nullptr && (bCalcDHDL || bCalcVir))
1182 if (lincsd->bTaskDep)
1184 /* In lincs_update_atoms threads might cross-read mlambda */
1188 /* Only account for local atoms */
1189 for (b = b0; b < b1; b++)
1191 mlambda[b] *= 0.5*nlocat[b];
1200 for (b = b0; b < b1; b++)
1202 /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1203 * later after the contributions are reduced over the threads.
1205 dhdl -= lincsd->mlambda[b]*lincsd->ddist[b];
1207 lincsd->task[th].dhdlambda = dhdl;
1212 /* Constraint virial */
1213 for (b = b0; b < b1; b++)
1217 tmp0 = -bllen[b]*mlambda[b];
1218 for (i = 0; i < DIM; i++)
1220 tmp1 = tmp0*r[b][i];
1221 for (j = 0; j < DIM; j++)
1223 vir_r_m_dr[i][j] -= tmp1*r[b][j];
1226 } /* 22 ncons flops */
1230 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1231 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1233 * (26+nrec)*ncons + (6+2*nrec)*nrtot
1234 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1236 * (63+nrec)*ncons + (6+4*nrec)*nrtot
1240 /*! \brief Sets the elements in the LINCS matrix for task task. */
1241 static void set_lincs_matrix_task(Lincs *li,
1243 const real *invmass,
1245 int *nCrossTaskTriangles)
1249 /* Construct the coupling coefficient matrix blmf */
1250 li_task->ntriangle = 0;
1252 *nCrossTaskTriangles = 0;
1253 for (i = li_task->b0; i < li_task->b1; i++)
1258 a2 = li->bla[2*i+1];
1259 for (n = li->blnr[i]; (n < li->blnr[i+1]); n++)
1261 int k, sign, center, end;
1265 /* If we are using multiple, independent tasks for LINCS,
1266 * the calls to check_assign_connected should have
1267 * put all connected constraints in our task.
1269 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1271 if (a1 == li->bla[2*k] || a2 == li->bla[2*k+1])
1279 if (a1 == li->bla[2*k] || a1 == li->bla[2*k+1])
1289 li->blmf[n] = sign*invmass[center]*li->blc[i]*li->blc[k];
1290 li->blmf1[n] = sign*0.5;
1291 if (li->ncg_triangle > 0)
1295 /* Look for constraint triangles */
1296 for (nk = li->blnr[k]; (nk < li->blnr[k+1]); nk++)
1299 if (kk != i && kk != k &&
1300 (li->bla[2*kk] == end || li->bla[2*kk+1] == end))
1302 /* Check if the constraints in this triangle actually
1303 * belong to a different task. We still assign them
1304 * here, since it's convenient for the triangle
1305 * iterations, but we then need an extra barrier.
1307 if (k < li_task->b0 || k >= li_task->b1 ||
1308 kk < li_task->b0 || kk >= li_task->b1)
1310 (*nCrossTaskTriangles)++;
1313 if (li_task->ntriangle == 0 ||
1314 li_task->triangle[li_task->ntriangle - 1] < i)
1316 /* Add this constraint to the triangle list */
1317 li_task->triangle[li_task->ntriangle] = i;
1318 li_task->tri_bits[li_task->ntriangle] = 0;
1319 li_task->ntriangle++;
1320 if (li->blnr[i+1] - li->blnr[i] > static_cast<int>(sizeof(li_task->tri_bits[0])*8 - 1))
1322 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %d allowed for constraints participating in triangles",
1323 li->blnr[i+1] - li->blnr[i],
1324 sizeof(li_task->tri_bits[0])*8-1);
1327 li_task->tri_bits[li_task->ntriangle-1] |= (1 << (n - li->blnr[i]));
1336 /*! \brief Sets the elements in the LINCS matrix. */
1337 static void set_lincs_matrix(Lincs *li, real *invmass, real lambda)
1340 const real invsqrt2 = 0.7071067811865475244;
1342 for (i = 0; (i < li->nc); i++)
1347 a2 = li->bla[2*i+1];
1348 li->blc[i] = gmx::invsqrt(invmass[a1] + invmass[a2]);
1349 li->blc1[i] = invsqrt2;
1352 /* Construct the coupling coefficient matrix blmf */
1353 int th, ntriangle = 0, ncc_triangle = 0, nCrossTaskTriangles = 0;
1354 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle, nCrossTaskTriangles) num_threads(li->ntask) schedule(static)
1355 for (th = 0; th < li->ntask; th++)
1359 set_lincs_matrix_task(li, &li->task[th], invmass,
1360 &ncc_triangle, &nCrossTaskTriangles);
1361 ntriangle = li->task[th].ntriangle;
1363 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1365 li->ntriangle = ntriangle;
1366 li->ncc_triangle = ncc_triangle;
1367 li->bTaskDepTri = (nCrossTaskTriangles > 0);
1371 fprintf(debug, "The %d constraints participate in %d triangles\n",
1372 li->nc, li->ntriangle);
1373 fprintf(debug, "There are %d constraint couplings, of which %d in triangles\n",
1374 li->ncc, li->ncc_triangle);
1375 if (li->ntriangle > 0 && li->ntask > 1)
1377 fprintf(debug, "%d constraint triangles contain constraints assigned to different tasks\n",
1378 nCrossTaskTriangles);
1383 * so we know with which lambda value the masses have been set.
1385 li->matlam = lambda;
1388 //! Finds all triangles of atoms that share constraints to a central atom.
1389 static int count_triangle_constraints(const t_ilist *ilist,
1390 const t_blocka &at2con)
1392 int ncon1, ncon_tot;
1393 int c0, a00, a01, n1, c1, a10, a11, ac1, n2, c2, a20, a21;
1396 t_iatom *ia1, *ia2, *iap;
1398 ncon1 = ilist[F_CONSTR].nr/3;
1399 ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
1401 ia1 = ilist[F_CONSTR].iatoms;
1402 ia2 = ilist[F_CONSTRNC].iatoms;
1405 for (c0 = 0; c0 < ncon_tot; c0++)
1408 iap = constr_iatomptr(ncon1, ia1, ia2, c0);
1411 for (n1 = at2con.index[a01]; n1 < at2con.index[a01+1]; n1++)
1416 iap = constr_iatomptr(ncon1, ia1, ia2, c1);
1427 for (n2 = at2con.index[ac1]; n2 < at2con.index[ac1+1]; n2++)
1430 if (c2 != c0 && c2 != c1)
1432 iap = constr_iatomptr(ncon1, ia1, ia2, c2);
1435 if (a20 == a00 || a21 == a00)
1449 return ncon_triangle;
1452 //! Finds sequences of sequential constraints.
1453 static bool more_than_two_sequential_constraints(const t_ilist *ilist,
1454 const t_blocka &at2con)
1456 t_iatom *ia1, *ia2, *iap;
1457 int ncon1, ncon_tot, c;
1459 bool bMoreThanTwoSequentialConstraints;
1461 ncon1 = ilist[F_CONSTR].nr/3;
1462 ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
1464 ia1 = ilist[F_CONSTR].iatoms;
1465 ia2 = ilist[F_CONSTRNC].iatoms;
1467 bMoreThanTwoSequentialConstraints = FALSE;
1468 for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
1470 iap = constr_iatomptr(ncon1, ia1, ia2, c);
1473 /* Check if this constraint has constraints connected at both atoms */
1474 if (at2con.index[a1+1] - at2con.index[a1] > 1 &&
1475 at2con.index[a2+1] - at2con.index[a2] > 1)
1477 bMoreThanTwoSequentialConstraints = TRUE;
1481 return bMoreThanTwoSequentialConstraints;
1484 //! Sorting helper function to compare two integers.
1485 static int int_comp(const void *a, const void *b)
1487 return (*(int *)a) - (*(int *)b);
1490 Lincs *init_lincs(FILE *fplog, const gmx_mtop_t &mtop,
1491 int nflexcon_global, ArrayRef<const t_blocka> at2con,
1492 bool bPLINCS, int nIter, int nProjOrder)
1494 // TODO this should become a unique_ptr
1496 bool bMoreThanTwoSeq;
1500 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
1501 bPLINCS ? " Parallel" : "");
1507 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1508 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1509 li->ncg_flex = nflexcon_global;
1512 li->nOrder = nProjOrder;
1514 li->max_connect = 0;
1515 for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1517 for (int a = 0; a < mtop.moltype[mt].atoms.nr; a++)
1519 li->max_connect = std::max(li->max_connect,
1520 at2con[mt].index[a + 1] - at2con[mt].index[a]);
1524 li->ncg_triangle = 0;
1525 bMoreThanTwoSeq = FALSE;
1526 for (const gmx_molblock_t &molb : mtop.molblock)
1528 const gmx_moltype_t &molt = mtop.moltype[molb.type];
1532 count_triangle_constraints(molt.ilist, at2con[molb.type]);
1534 if (!bMoreThanTwoSeq &&
1535 more_than_two_sequential_constraints(molt.ilist, at2con[molb.type]))
1537 bMoreThanTwoSeq = TRUE;
1541 /* Check if we need to communicate not only before LINCS,
1542 * but also before each iteration.
1543 * The check for only two sequential constraints is only
1544 * useful for the common case of H-bond only constraints.
1545 * With more effort we could also make it useful for small
1546 * molecules with nr. sequential constraints <= nOrder-1.
1548 li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1550 if (debug && bPLINCS)
1552 fprintf(debug, "PLINCS communication before each iteration: %d\n",
1553 int{li->bCommIter});
1556 /* LINCS can run on any number of threads.
1557 * Currently the number is fixed for the whole simulation,
1558 * but it could be set in set_lincs().
1559 * The current constraint to task assignment code can create independent
1560 * tasks only when not more than two constraints are connected sequentially.
1562 li->ntask = gmx_omp_nthreads_get(emntLINCS);
1563 li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1566 fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n",
1567 li->ntask, li->bTaskDep ? "" : "in");
1575 /* Allocate an extra elements for "task-overlap" constraints */
1576 li->task.resize(li->ntask + 1);
1579 if (bPLINCS || li->ncg_triangle > 0)
1581 please_cite(fplog, "Hess2008a");
1585 please_cite(fplog, "Hess97a");
1590 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1593 fprintf(fplog, "There are inter charge-group constraints,\n"
1594 "will communicate selected coordinates each lincs iteration\n");
1596 if (li->ncg_triangle > 0)
1599 "%d constraints are involved in constraint triangles,\n"
1600 "will apply an additional matrix expansion of order %d for couplings\n"
1601 "between constraints inside triangles\n",
1602 li->ncg_triangle, li->nOrder);
1609 void done_lincs(Lincs *li)
1614 /*! \brief Sets up the work division over the threads. */
1615 static void lincs_thread_setup(Lincs *li, int natoms)
1622 if (natoms > li->atf_nalloc)
1624 li->atf_nalloc = over_alloc_large(natoms);
1625 srenew(li->atf, li->atf_nalloc);
1629 /* Clear the atom flags */
1630 for (a = 0; a < natoms; a++)
1632 bitmask_clear(&atf[a]);
1635 if (li->ntask > BITMASK_SIZE)
1637 gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1640 for (th = 0; th < li->ntask; th++)
1645 li_task = &li->task[th];
1647 /* For each atom set a flag for constraints from each */
1648 for (b = li_task->b0; b < li_task->b1; b++)
1650 bitmask_set_bit(&atf[li->bla[b*2 ]], th);
1651 bitmask_set_bit(&atf[li->bla[b*2 + 1]], th);
1655 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1656 for (th = 0; th < li->ntask; th++)
1664 li_task = &li->task[th];
1666 if (li_task->b1 - li_task->b0 > li_task->ind_nalloc)
1668 li_task->ind_nalloc = over_alloc_large(li_task->b1-li_task->b0);
1669 srenew(li_task->ind, li_task->ind_nalloc);
1670 srenew(li_task->ind_r, li_task->ind_nalloc);
1673 bitmask_init_low_bits(&mask, th);
1676 li_task->nind_r = 0;
1677 for (b = li_task->b0; b < li_task->b1; b++)
1679 /* We let the constraint with the lowest thread index
1680 * operate on atoms with constraints from multiple threads.
1682 if (bitmask_is_disjoint(atf[li->bla[b*2]], mask) &&
1683 bitmask_is_disjoint(atf[li->bla[b*2+1]], mask))
1685 /* Add the constraint to the local atom update index */
1686 li_task->ind[li_task->nind++] = b;
1690 /* Add the constraint to the rest block */
1691 li_task->ind_r[li_task->nind_r++] = b;
1695 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1698 /* We need to copy all constraints which have not be assigned
1699 * to a thread to a separate list which will be handled by one thread.
1701 li_m = &li->task[li->ntask];
1704 for (th = 0; th < li->ntask; th++)
1709 li_task = &li->task[th];
1711 if (li_m->nind + li_task->nind_r > li_m->ind_nalloc)
1713 li_m->ind_nalloc = over_alloc_large(li_m->nind+li_task->nind_r);
1714 srenew(li_m->ind, li_m->ind_nalloc);
1717 for (b = 0; b < li_task->nind_r; b++)
1719 li_m->ind[li_m->nind++] = li_task->ind_r[b];
1724 fprintf(debug, "LINCS thread %d: %d constraints\n",
1731 fprintf(debug, "LINCS thread r: %d constraints\n",
1736 /*! \brief There is no realloc with alignment, so here we make one for reals.
1737 * Note that this function does not preserve the contents of the memory.
1739 static void resize_real_aligned(real **ptr, int nelem)
1741 sfree_aligned(*ptr);
1742 snew_aligned(*ptr, nelem, align_bytes);
1745 //! Assign a constraint.
1746 static void assign_constraint(Lincs *li,
1747 int constraint_index,
1749 real lenA, real lenB,
1750 const t_blocka *at2con)
1756 /* Make an mapping of local topology constraint index to LINCS index */
1757 li->con_index[constraint_index] = con;
1759 li->bllen0[con] = lenA;
1760 li->ddist[con] = lenB - lenA;
1761 /* Set the length to the topology A length */
1762 li->bllen[con] = lenA;
1763 li->bla[2*con] = a1;
1764 li->bla[2*con+1] = a2;
1766 /* Make space in the constraint connection matrix for constraints
1767 * connected to both end of the current constraint.
1770 at2con->index[a1 + 1] - at2con->index[a1] - 1 +
1771 at2con->index[a2 + 1] - at2con->index[a2] - 1;
1773 li->blnr[con + 1] = li->ncc;
1775 /* Increase the constraint count */
1779 /*! \brief Check if constraint with topology index constraint_index is connected
1780 * to other constraints, and if so add those connected constraints to our task. */
1781 static void check_assign_connected(Lincs *li,
1782 const t_iatom *iatom,
1786 const t_blocka *at2con)
1788 /* Currently this function only supports constraint groups
1789 * in which all constraints share at least one atom
1790 * (e.g. H-bond constraints).
1791 * Check both ends of the current constraint for
1792 * connected constraints. We need to assign those
1797 for (end = 0; end < 2; end++)
1801 a = (end == 0 ? a1 : a2);
1803 for (k = at2con->index[a]; k < at2con->index[a + 1]; k++)
1808 /* Check if constraint cc has not yet been assigned */
1809 if (li->con_index[cc] == -1)
1815 lenA = idef.iparams[type].constr.dA;
1816 lenB = idef.iparams[type].constr.dB;
1818 if (bDynamics || lenA != 0 || lenB != 0)
1820 assign_constraint(li, cc, iatom[3*cc + 1], iatom[3*cc + 2], lenA, lenB, at2con);
1827 /*! \brief Check if constraint with topology index constraint_index is involved
1828 * in a constraint triangle, and if so add the other two constraints
1829 * in the triangle to our task. */
1830 static void check_assign_triangle(Lincs *li,
1831 const t_iatom *iatom,
1834 int constraint_index,
1836 const t_blocka *at2con)
1838 int nca, cc[32], ca[32], k;
1839 int c_triangle[2] = { -1, -1 };
1842 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1847 if (c != constraint_index)
1851 aa1 = iatom[c*3 + 1];
1852 aa2 = iatom[c*3 + 2];
1868 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1873 if (c != constraint_index)
1877 aa1 = iatom[c*3 + 1];
1878 aa2 = iatom[c*3 + 2];
1881 for (i = 0; i < nca; i++)
1885 c_triangle[0] = cc[i];
1892 for (i = 0; i < nca; i++)
1896 c_triangle[0] = cc[i];
1904 if (c_triangle[0] >= 0)
1908 for (end = 0; end < 2; end++)
1910 /* Check if constraint c_triangle[end] has not yet been assigned */
1911 if (li->con_index[c_triangle[end]] == -1)
1916 i = c_triangle[end]*3;
1918 lenA = idef.iparams[type].constr.dA;
1919 lenB = idef.iparams[type].constr.dB;
1921 if (bDynamics || lenA != 0 || lenB != 0)
1923 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1930 //! Sets matrix indices.
1931 static void set_matrix_indices(Lincs *li,
1932 const Task *li_task,
1933 const t_blocka *at2con,
1938 for (b = li_task->b0; b < li_task->b1; b++)
1943 a2 = li->bla[b*2 + 1];
1946 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1950 concon = li->con_index[at2con->a[k]];
1953 li->blbnb[i++] = concon;
1956 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1960 concon = li->con_index[at2con->a[k]];
1963 li->blbnb[i++] = concon;
1969 /* Order the blbnb matrix to optimize memory access */
1970 qsort(&(li->blbnb[li->blnr[b]]), li->blnr[b + 1] - li->blnr[b],
1971 sizeof(li->blbnb[0]), int_comp);
1976 void set_lincs(const t_idef &idef,
1977 const t_mdatoms &md,
1979 const t_commrec *cr,
1985 int i, ncc_alloc_old, ncon_tot;
1990 /* Zero the thread index ranges.
1991 * Otherwise without local constraints we could return with old ranges.
1993 for (i = 0; i < li->ntask; i++)
1997 li->task[i].nind = 0;
2001 li->task[li->ntask].nind = 0;
2004 /* This is the local topology, so there are only F_CONSTR constraints */
2005 if (idef.il[F_CONSTR].nr == 0)
2007 /* There are no constraints,
2008 * we do not need to fill any data structures.
2015 fprintf(debug, "Building the LINCS connectivity\n");
2018 if (DOMAINDECOMP(cr))
2020 if (cr->dd->constraints)
2024 dd_get_constraint_range(cr->dd, &start, &natoms);
2028 natoms = dd_numHomeAtoms(*cr->dd);
2036 at2con = make_at2con(natoms, idef.il, idef.iparams,
2037 flexibleConstraintTreatment(bDynamics));
2039 ncon_tot = idef.il[F_CONSTR].nr/3;
2041 /* Ensure we have enough padding for aligned loads for each thread */
2042 if (ncon_tot + li->ntask*simd_width > li->nc_alloc || li->nc_alloc == 0)
2044 li->nc_alloc = over_alloc_dd(ncon_tot + li->ntask*simd_width);
2045 srenew(li->con_index, li->nc_alloc);
2046 resize_real_aligned(&li->bllen0, li->nc_alloc);
2047 resize_real_aligned(&li->ddist, li->nc_alloc);
2048 srenew(li->bla, 2*li->nc_alloc);
2049 resize_real_aligned(&li->blc, li->nc_alloc);
2050 resize_real_aligned(&li->blc1, li->nc_alloc);
2051 srenew(li->blnr, li->nc_alloc + 1);
2052 resize_real_aligned(&li->bllen, li->nc_alloc);
2053 srenew(li->tmpv, li->nc_alloc);
2054 if (DOMAINDECOMP(cr))
2056 srenew(li->nlocat, li->nc_alloc);
2058 resize_real_aligned(&li->tmp1, li->nc_alloc);
2059 resize_real_aligned(&li->tmp2, li->nc_alloc);
2060 resize_real_aligned(&li->tmp3, li->nc_alloc);
2061 resize_real_aligned(&li->tmp4, li->nc_alloc);
2062 resize_real_aligned(&li->mlambda, li->nc_alloc);
2065 iatom = idef.il[F_CONSTR].iatoms;
2067 ncc_alloc_old = li->ncc_alloc;
2068 li->blnr[0] = li->ncc;
2070 /* Assign the constraints for li->ntask LINCS tasks.
2071 * We target a uniform distribution of constraints over the tasks.
2072 * Note that when flexible constraints are present, but are removed here
2073 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
2074 * happen during normal MD, that's ok.
2076 int ncon_assign, ncon_target, con, th;
2078 /* Determine the number of constraints we need to assign here */
2079 ncon_assign = ncon_tot;
2082 /* With energy minimization, flexible constraints are ignored
2083 * (and thus minimized, as they should be).
2085 ncon_assign -= countFlexibleConstraints(idef.il, idef.iparams);
2088 /* Set the target constraint count per task to exactly uniform,
2089 * this might be overridden below.
2091 ncon_target = (ncon_assign + li->ntask - 1)/li->ntask;
2093 /* Mark all constraints as unassigned by setting their index to -1 */
2094 for (con = 0; con < ncon_tot; con++)
2096 li->con_index[con] = -1;
2100 for (th = 0; th < li->ntask; th++)
2104 li_task = &li->task[th];
2106 #if GMX_SIMD_HAVE_REAL
2107 /* With indepedent tasks we likely have H-bond constraints or constraint
2108 * pairs. The connected constraints will be pulled into the task, so the
2109 * constraints per task will often exceed ncon_target.
2110 * Triangle constraints can also increase the count, but there are
2111 * relatively few of those, so we usually expect to get ncon_target.
2115 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
2116 * since otherwise a lot of operations can be wasted.
2117 * There are several ways to round here, we choose the one
2118 * that alternates block sizes, which helps with Intel HT.
2120 ncon_target = ((ncon_assign*(th + 1))/li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1) & ~(GMX_SIMD_REAL_WIDTH - 1);
2122 #endif // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
2124 /* Continue filling the arrays where we left off with the previous task,
2125 * including padding for SIMD.
2127 li_task->b0 = li->nc;
2129 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2131 if (li->con_index[con] == -1)
2136 type = iatom[3*con];
2137 a1 = iatom[3*con + 1];
2138 a2 = iatom[3*con + 2];
2139 lenA = idef.iparams[type].constr.dA;
2140 lenB = idef.iparams[type].constr.dB;
2141 /* Skip the flexible constraints when not doing dynamics */
2142 if (bDynamics || lenA != 0 || lenB != 0)
2144 assign_constraint(li, con, a1, a2, lenA, lenB, &at2con);
2146 if (li->ntask > 1 && !li->bTaskDep)
2148 /* We can generate independent tasks. Check if we
2149 * need to assign connected constraints to our task.
2151 check_assign_connected(li, iatom, idef, bDynamics,
2154 if (li->ntask > 1 && li->ncg_triangle > 0)
2156 /* Ensure constraints in one triangle are assigned
2159 check_assign_triangle(li, iatom, idef, bDynamics,
2160 con, a1, a2, &at2con);
2168 li_task->b1 = li->nc;
2172 /* Copy the last atom pair indices and lengths for constraints
2173 * up to a multiple of simd_width, such that we can do all
2174 * SIMD operations without having to worry about end effects.
2178 li->nc = ((li_task->b1 + simd_width - 1)/simd_width)*simd_width;
2179 last = li_task->b1 - 1;
2180 for (i = li_task->b1; i < li->nc; i++)
2182 li->bla[i*2 ] = li->bla[last*2 ];
2183 li->bla[i*2 + 1] = li->bla[last*2 + 1];
2184 li->bllen0[i] = li->bllen0[last];
2185 li->ddist[i] = li->ddist[last];
2186 li->bllen[i] = li->bllen[last];
2187 li->blnr[i + 1] = li->blnr[last + 1];
2191 /* Keep track of how many constraints we assigned */
2192 li->nc_real += li_task->b1 - li_task->b0;
2196 fprintf(debug, "LINCS task %d constraints %d - %d\n",
2197 th, li_task->b0, li_task->b1);
2201 assert(li->nc_real == ncon_assign);
2205 /* Without DD we order the blbnb matrix to optimize memory access.
2206 * With DD the overhead of sorting is more than the gain during access.
2208 bSortMatrix = !DOMAINDECOMP(cr);
2210 if (li->ncc > li->ncc_alloc)
2212 li->ncc_alloc = over_alloc_small(li->ncc);
2213 srenew(li->blbnb, li->ncc_alloc);
2216 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2217 for (th = 0; th < li->ntask; th++)
2223 li_task = &li->task[th];
2225 if (li->ncg_triangle > 0 &&
2226 li_task->b1 - li_task->b0 > li_task->tri_alloc)
2228 /* This is allocating too much, but it is difficult to improve */
2229 li_task->tri_alloc = over_alloc_dd(li_task->b1 - li_task->b0);
2230 srenew(li_task->triangle, li_task->tri_alloc);
2231 srenew(li_task->tri_bits, li_task->tri_alloc);
2234 set_matrix_indices(li, li_task, &at2con, bSortMatrix);
2236 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2239 done_blocka(&at2con);
2241 if (cr->dd == nullptr)
2243 /* Since the matrix is static, we should free some memory */
2244 li->ncc_alloc = li->ncc;
2245 srenew(li->blbnb, li->ncc_alloc);
2248 if (li->ncc_alloc > ncc_alloc_old)
2250 srenew(li->blmf, li->ncc_alloc);
2251 srenew(li->blmf1, li->ncc_alloc);
2252 srenew(li->tmpncc, li->ncc_alloc);
2255 gmx::ArrayRef<const int> nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2256 if (!nlocat_dd.empty())
2258 /* Convert nlocat from local topology to LINCS constraint indexing */
2259 for (con = 0; con < ncon_tot; con++)
2261 li->nlocat[li->con_index[con]] = nlocat_dd[con];
2266 li->nlocat = nullptr;
2271 fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n",
2272 li->nc_real, li->nc, li->ncc);
2277 lincs_thread_setup(li, md.nr);
2280 set_lincs_matrix(li, md.invmass, md.lambda);
2283 //! Issues a warning when LINCS constraints cannot be satisfied.
2284 static void lincs_warning(gmx_domdec_t *dd, const rvec *x, rvec *xprime, t_pbc *pbc,
2285 int ncons, const int *bla, real *bllen, real wangle,
2286 int maxwarn, int *warncount)
2290 real wfac, d0, d1, cosine;
2292 wfac = std::cos(DEG2RAD*wangle);
2295 "bonds that rotated more than %g degrees:\n"
2296 " atom 1 atom 2 angle previous, current, constraint length\n",
2299 for (b = 0; b < ncons; b++)
2305 pbc_dx_aiuc(pbc, x[i], x[j], v0);
2306 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2310 rvec_sub(x[i], x[j], v0);
2311 rvec_sub(xprime[i], xprime[j], v1);
2315 cosine = ::iprod(v0, v1)/(d0*d1);
2319 " %6d %6d %5.1f %8.4f %8.4f %8.4f\n",
2320 ddglatnr(dd, i), ddglatnr(dd, j),
2321 RAD2DEG*std::acos(cosine), d0, d1, bllen[b]);
2322 if (!std::isfinite(d1))
2324 gmx_fatal(FARGS, "Bond length not finite.");
2330 if (*warncount > maxwarn)
2332 too_many_constraint_warnings(econtLINCS, *warncount);
2336 //! Determine how well the constraints have been satisfied.
2337 static void cconerr(const Lincs *lincsd,
2338 rvec *x, t_pbc *pbc,
2339 real *ncons_loc, real *ssd, real *max, int *imax)
2341 const int *bla, *nlocat;
2344 int count, im, task;
2347 bllen = lincsd->bllen;
2348 nlocat = lincsd->nlocat;
2354 for (task = 0; task < lincsd->ntask; task++)
2358 for (b = lincsd->task[task].b0; b < lincsd->task[task].b1; b++)
2365 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
2369 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
2372 len = r2*gmx::invsqrt(r2);
2373 d = std::abs(len/bllen[b]-1);
2374 if (d > ma && (nlocat == nullptr || nlocat[b]))
2379 if (nlocat == nullptr)
2386 ssd2 += nlocat[b]*d*d;
2392 *ncons_loc = (nlocat ? 0.5 : 1)*count;
2393 *ssd = (nlocat ? 0.5 : 1)*ssd2;
2398 bool constrain_lincs(bool computeRmsd,
2399 const t_inputrec &ir,
2401 Lincs *lincsd, const t_mdatoms &md,
2402 const t_commrec *cr,
2403 const gmx_multisim_t &ms,
2404 const rvec *x, rvec *xprime, rvec *min_proj,
2405 matrix box, t_pbc *pbc,
2406 real lambda, real *dvdlambda,
2407 real invdt, rvec *v,
2408 bool bCalcVir, tensor vir_r_m_dr,
2409 ConstraintVariable econq,
2411 int maxwarn, int *warncount)
2414 char buf2[22], buf3[STRLEN];
2416 real ncons_loc, p_ssd, p_max = 0;
2422 /* This boolean should be set by a flag passed to this routine.
2423 * We can also easily check if any constraint length is changed,
2424 * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2426 bCalcDHDL = (ir.efep != efepNO && dvdlambda != nullptr);
2428 if (lincsd->nc == 0 && cr->dd == nullptr)
2432 lincsd->rmsdData = {{0}};
2438 if (econq == ConstraintVariable::Positions)
2440 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2441 * also with efep!=fepNO.
2443 if (ir.efep != efepNO)
2445 if (md.nMassPerturbed && lincsd->matlam != md.lambda)
2447 set_lincs_matrix(lincsd, md.invmass, md.lambda);
2450 for (i = 0; i < lincsd->nc; i++)
2452 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
2456 if (lincsd->ncg_flex)
2458 /* Set the flexible constraint lengths to the old lengths */
2461 for (i = 0; i < lincsd->nc; i++)
2463 if (lincsd->bllen[i] == 0)
2465 pbc_dx_aiuc(pbc, x[lincsd->bla[2*i]], x[lincsd->bla[2*i+1]], dx);
2466 lincsd->bllen[i] = norm(dx);
2472 for (i = 0; i < lincsd->nc; i++)
2474 if (lincsd->bllen[i] == 0)
2477 std::sqrt(distance2(x[lincsd->bla[2*i]],
2478 x[lincsd->bla[2*i+1]]));
2486 cconerr(lincsd, xprime, pbc,
2487 &ncons_loc, &p_ssd, &p_max, &p_imax);
2490 /* This bWarn var can be updated by multiple threads
2491 * at the same time. But as we only need to detect
2492 * if a warning occurred or not, this is not an issue.
2496 /* The OpenMP parallel region of constrain_lincs for coords */
2497 #pragma omp parallel num_threads(lincsd->ntask)
2501 int th = gmx_omp_get_thread_num();
2503 clear_mat(lincsd->task[th].vir_r_m_dr);
2505 do_lincs(x, xprime, box, pbc, lincsd, th,
2508 ir.LincsWarnAngle, &bWarn,
2510 th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2512 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2515 if (debug && lincsd->nc > 0)
2517 fprintf(debug, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2518 fprintf(debug, " Before LINCS %.6f %.6f %6d %6d\n",
2519 std::sqrt(p_ssd/ncons_loc), p_max,
2520 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2521 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2523 if (computeRmsd || debug)
2525 cconerr(lincsd, xprime, pbc,
2526 &ncons_loc, &p_ssd, &p_max, &p_imax);
2527 lincsd->rmsdData[0] = ncons_loc;
2528 lincsd->rmsdData[1] = p_ssd;
2532 lincsd->rmsdData = {{0}};
2534 if (debug && lincsd->nc > 0)
2537 " After LINCS %.6f %.6f %6d %6d\n\n",
2538 std::sqrt(p_ssd/ncons_loc), p_max,
2539 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2540 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2545 if (maxwarn < INT_MAX)
2547 cconerr(lincsd, xprime, pbc,
2548 &ncons_loc, &p_ssd, &p_max, &p_imax);
2549 if (isMultiSim(&ms))
2551 sprintf(buf3, " in simulation %d", ms.sim);
2558 "\nStep %s, time %g (ps) LINCS WARNING%s\n"
2559 "relative constraint deviation after LINCS:\n"
2560 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2561 gmx_step_str(step, buf2), ir.init_t+step*ir.delta_t,
2563 std::sqrt(p_ssd/ncons_loc), p_max,
2564 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2565 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2567 lincs_warning(cr->dd, x, xprime, pbc,
2568 lincsd->nc, lincsd->bla, lincsd->bllen,
2569 ir.LincsWarnAngle, maxwarn, warncount);
2571 bOK = (p_max < 0.5);
2574 if (lincsd->ncg_flex)
2576 for (i = 0; (i < lincsd->nc); i++)
2578 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2580 lincsd->bllen[i] = 0;
2587 /* The OpenMP parallel region of constrain_lincs for derivatives */
2588 #pragma omp parallel num_threads(lincsd->ntask)
2592 int th = gmx_omp_get_thread_num();
2594 do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
2595 md.invmass, econq, bCalcDHDL,
2596 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2598 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2604 /* Reduce the dH/dlambda contributions over the threads */
2609 for (th = 0; th < lincsd->ntask; th++)
2611 dhdlambda += lincsd->task[th].dhdlambda;
2613 if (econq == ConstraintVariable::Positions)
2615 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2616 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2617 dhdlambda /= ir.delta_t*ir.delta_t;
2619 *dvdlambda += dhdlambda;
2622 if (bCalcVir && lincsd->ntask > 1)
2624 for (i = 1; i < lincsd->ntask; i++)
2626 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2630 /* count assuming nit=1 */
2631 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2632 inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
2633 if (lincsd->ntriangle > 0)
2635 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
2639 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real*2);
2643 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);