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, 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.
37 /* This file is completely threadsafe - keep it that way! */
49 #include "gromacs/domdec/domdec.h"
50 #include "gromacs/domdec/domdec_struct.h"
51 #include "gromacs/gmxlib/nrnb.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdlib/constr.h"
56 #include "gromacs/mdlib/gmx_omp_nthreads.h"
57 #include "gromacs/mdlib/mdrun.h"
58 #include "gromacs/mdtypes/commrec.h"
59 #include "gromacs/mdtypes/inputrec.h"
60 #include "gromacs/mdtypes/md_enums.h"
61 #include "gromacs/pbcutil/pbc.h"
62 #include "gromacs/pbcutil/pbc-simd.h"
63 #include "gromacs/simd/simd.h"
64 #include "gromacs/simd/simd_math.h"
65 #include "gromacs/simd/vector_operations.h"
66 #include "gromacs/topology/block.h"
67 #include "gromacs/topology/mtop_util.h"
68 #include "gromacs/utility/basedefinitions.h"
69 #include "gromacs/utility/bitmask.h"
70 #include "gromacs/utility/cstringutil.h"
71 #include "gromacs/utility/exceptions.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/gmxomp.h"
74 #include "gromacs/utility/pleasecite.h"
75 #include "gromacs/utility/smalloc.h"
77 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
80 int b0; /* first constraint for this task */
81 int b1; /* b1-1 is the last constraint for this task */
82 int ntriangle; /* the number of constraints in triangles */
83 int *triangle; /* the list of triangle constraints */
84 int *tri_bits; /* the bits tell if the matrix element should be used */
85 int tri_alloc; /* allocation size of triangle and tri_bits */
86 int nind; /* number of indices */
87 int *ind; /* constraint index for updating atom data */
88 int nind_r; /* number of indices */
89 int *ind_r; /* constraint index for updating atom data */
90 int ind_nalloc; /* allocation size of ind and ind_r */
91 tensor vir_r_m_dr; /* temporary variable for virial calculation */
92 real dhdlambda; /* temporary variable for lambda derivative */
95 typedef struct gmx_lincsdata {
96 int ncg; /* the global number of constraints */
97 int ncg_flex; /* the global number of flexible constraints */
98 int ncg_triangle; /* the global number of constraints in triangles */
99 int nIter; /* the number of iterations */
100 int nOrder; /* the order of the matrix expansion */
101 int max_connect; /* the maximum number of constrains connected to a single atom */
103 int nc_real; /* the number of real constraints */
104 int nc; /* the number of constraints including padding for SIMD */
105 int nc_alloc; /* the number we allocated memory for */
106 int ncc; /* the number of constraint connections */
107 int ncc_alloc; /* the number we allocated memory for */
108 real matlam; /* the FE lambda value used for filling blc and blmf */
109 int *con_index; /* mapping from topology to LINCS constraints */
110 real *bllen0; /* the reference distance in topology A */
111 real *ddist; /* the reference distance in top B - the r.d. in top A */
112 int *bla; /* the atom pairs involved in the constraints */
113 real *blc; /* 1/sqrt(invmass1 + invmass2) */
114 real *blc1; /* as blc, but with all masses 1 */
115 int *blnr; /* index into blbnb and blmf */
116 int *blbnb; /* list of constraint connections */
117 int ntriangle; /* the local number of constraints in triangles */
118 int ncc_triangle; /* the number of constraint connections in triangles */
119 gmx_bool bCommIter; /* communicate before each LINCS interation */
120 real *blmf; /* matrix of mass factors for constraint connections */
121 real *blmf1; /* as blmf, but with all masses 1 */
122 real *bllen; /* the reference bond length */
123 int *nlocat; /* the local atom count per constraint, can be NULL */
125 int ntask; /* The number of tasks = #threads for LINCS */
126 lincs_task_t *task; /* LINCS thread division */
127 gmx_bitmask_t *atf; /* atom flags for thread parallelization */
128 int atf_nalloc; /* allocation size of atf */
129 gmx_bool bTaskDep; /* are the LINCS tasks interdependent? */
130 /* arrays for temporary storage in the LINCS algorithm */
137 real *mlambda; /* the Lagrange multipliers * -1 */
138 /* storage for the constraint RMS relative deviation output */
142 /* Define simd_width for memory allocation used for SIMD code */
143 #if GMX_SIMD_HAVE_REAL
144 static const int simd_width = GMX_SIMD_REAL_WIDTH;
146 static const int simd_width = 1;
149 /* Align to 128 bytes, consistent with the current implementation of
150 AlignedAllocator, which currently forces 128 byte alignment. */
151 static const int align_bytes = 128;
153 real *lincs_rmsd_data(struct gmx_lincsdata *lincsd)
155 return lincsd->rmsd_data;
158 real lincs_rmsd(struct gmx_lincsdata *lincsd)
160 if (lincsd->rmsd_data[0] > 0)
162 return std::sqrt(lincsd->rmsd_data[1]/lincsd->rmsd_data[0]);
170 /* Do a set of nrec LINCS matrix multiplications.
171 * This function will return with up to date thread-local
172 * constraint data, without an OpenMP barrier.
174 static void lincs_matrix_expand(const struct gmx_lincsdata *lincsd,
175 const lincs_task_t *li_task,
177 real *rhs1, real *rhs2, real *sol)
179 int b0, b1, nrec, rec;
180 const int *blnr = lincsd->blnr;
181 const int *blbnb = lincsd->blbnb;
185 nrec = lincsd->nOrder;
187 for (rec = 0; rec < nrec; rec++)
191 if (lincsd->bTaskDep)
195 for (b = b0; b < b1; b++)
201 for (n = blnr[b]; n < blnr[b+1]; n++)
203 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
206 sol[b] = sol[b] + mvb;
214 } /* nrec*(ncons+2*nrtot) flops */
216 if (lincsd->ntriangle > 0)
218 /* Perform an extra nrec recursions for only the constraints
219 * involved in rigid triangles.
220 * In this way their accuracy should come close to those of the other
221 * constraints, since traingles of constraints can produce eigenvalues
222 * around 0.7, while the effective eigenvalue for bond constraints
223 * is around 0.4 (and 0.7*0.7=0.5).
226 if (lincsd->bTaskDep)
228 /* We need a barrier here, since other threads might still be
229 * reading the contents of rhs1 and/o rhs2.
230 * We could avoid this barrier by introducing two extra rhs
231 * arrays for the triangle constraints only.
236 /* Constraints involved in a triangle are ensured to be in the same
237 * LINCS task. This means no barriers are required during the extra
238 * iterations for the triangle constraints.
240 const int *triangle = li_task->triangle;
241 const int *tri_bits = li_task->tri_bits;
243 for (rec = 0; rec < nrec; rec++)
247 for (tb = 0; tb < li_task->ntriangle; tb++)
249 int b, bits, nr0, nr1, n;
257 for (n = nr0; n < nr1; n++)
259 if (bits & (1 << (n - nr0)))
261 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
265 sol[b] = sol[b] + mvb;
273 } /* nrec*(ntriangle + ncc_triangle*2) flops */
277 static void lincs_update_atoms_noind(int ncons, const int *bla,
279 const real *fac, rvec *r,
284 real mvb, im1, im2, tmp0, tmp1, tmp2;
288 for (b = 0; b < ncons; b++)
304 } /* 16 ncons flops */
308 for (b = 0; b < ncons; b++)
326 static void lincs_update_atoms_ind(int ncons, const int *ind, const int *bla,
328 const real *fac, rvec *r,
333 real mvb, im1, im2, tmp0, tmp1, tmp2;
337 for (bi = 0; bi < ncons; bi++)
354 } /* 16 ncons flops */
358 for (bi = 0; bi < ncons; bi++)
373 } /* 16 ncons flops */
377 static void lincs_update_atoms(struct gmx_lincsdata *li, int th,
379 const real *fac, rvec *r,
385 /* Single thread, we simply update for all constraints */
386 lincs_update_atoms_noind(li->nc_real,
387 li->bla, prefac, fac, r, invmass, x);
391 /* Update the atom vector components for our thread local
392 * constraints that only access our local atom range.
393 * This can be done without a barrier.
395 lincs_update_atoms_ind(li->task[th].nind, li->task[th].ind,
396 li->bla, prefac, fac, r, invmass, x);
398 if (li->task[li->ntask].nind > 0)
400 /* Update the constraints that operate on atoms
401 * in multiple thread atom blocks on the master thread.
406 lincs_update_atoms_ind(li->task[li->ntask].nind,
407 li->task[li->ntask].ind,
408 li->bla, prefac, fac, r, invmass, x);
414 #if GMX_SIMD_HAVE_REAL
415 /* Calculate the constraint distance vectors r to project on from x.
416 * Determine the right-hand side of the matrix equation using quantity f.
417 * This function only differs from calc_dr_x_xp_simd below in that
418 * no constraint length is subtracted and no PBC is used for f.
420 static void gmx_simdcall
421 calc_dr_x_f_simd(int b0,
424 const rvec * gmx_restrict x,
425 const rvec * gmx_restrict f,
426 const real * gmx_restrict blc,
427 const real * pbc_simd,
428 rvec * gmx_restrict r,
429 real * gmx_restrict rhs,
430 real * gmx_restrict sol)
432 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
434 GMX_ALIGNED(int, GMX_SIMD_REAL_WIDTH) offset2[GMX_SIMD_REAL_WIDTH];
436 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
441 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
443 SimdReal x0_S, y0_S, z0_S;
444 SimdReal x1_S, y1_S, z1_S;
445 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
446 SimdReal fx_S, fy_S, fz_S, ip_S, rhs_S;
447 GMX_ALIGNED(int, GMX_SIMD_REAL_WIDTH) offset0[GMX_SIMD_REAL_WIDTH];
448 GMX_ALIGNED(int, GMX_SIMD_REAL_WIDTH) offset1[GMX_SIMD_REAL_WIDTH];
450 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
452 offset0[i] = bla[bs*2 + i*2];
453 offset1[i] = bla[bs*2 + i*2 + 1];
456 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
457 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
462 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
464 n2_S = norm2(rx_S, ry_S, rz_S);
465 il_S = invsqrt(n2_S);
471 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
473 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(f), offset0, &x0_S, &y0_S, &z0_S);
474 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(f), offset1, &x1_S, &y1_S, &z1_S);
479 ip_S = iprod(rx_S, ry_S, rz_S, fx_S, fy_S, fz_S);
481 rhs_S = load(blc + bs) * ip_S;
483 store(rhs + bs, rhs_S);
484 store(sol + bs, rhs_S);
487 #endif // GMX_SIMD_HAVE_REAL
489 /* LINCS projection, works on derivatives of the coordinates */
490 static void do_lincsp(rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
491 struct gmx_lincsdata *lincsd, int th,
493 int econq, gmx_bool bCalcDHDL,
494 gmx_bool bCalcVir, tensor rmdf)
497 int *bla, *blnr, *blbnb;
499 real *blc, *blmf, *blcc, *rhs1, *rhs2, *sol;
501 b0 = lincsd->task[th].b0;
502 b1 = lincsd->task[th].b1;
507 blbnb = lincsd->blbnb;
508 if (econq != econqForce)
510 /* Use mass-weighted parameters */
516 /* Use non mass-weighted parameters */
518 blmf = lincsd->blmf1;
520 blcc = lincsd->tmpncc;
525 #if GMX_SIMD_HAVE_REAL
526 /* This SIMD code does the same as the plain-C code after the #else.
527 * The only difference is that we always call pbc code, as with SIMD
528 * the overhead of pbc computation (when not needed) is small.
530 GMX_ALIGNED(real, GMX_SIMD_REAL_WIDTH) pbc_simd[9*GMX_SIMD_REAL_WIDTH];
532 /* Convert the pbc struct for SIMD */
533 set_pbc_simd(pbc, pbc_simd);
535 /* Compute normalized x i-j vectors, store in r.
536 * Compute the inner product of r and xp i-j and store in rhs1.
538 calc_dr_x_f_simd(b0, b1, bla, x, f, blc,
542 #else // GMX_SIMD_HAVE_REAL
544 /* Compute normalized i-j vectors */
547 for (b = b0; b < b1; b++)
551 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
557 for (b = b0; b < b1; b++)
561 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
563 } /* 16 ncons flops */
566 for (b = b0; b < b1; b++)
573 mvb = blc[b]*(r[b][0]*(f[i][0] - f[j][0]) +
574 r[b][1]*(f[i][1] - f[j][1]) +
575 r[b][2]*(f[i][2] - f[j][2]));
581 #endif // GMX_SIMD_HAVE_REAL
583 if (lincsd->bTaskDep)
585 /* We need a barrier, since the matrix construction below
586 * can access entries in r of other threads.
591 /* Construct the (sparse) LINCS matrix */
592 for (b = b0; b < b1; b++)
596 for (n = blnr[b]; n < blnr[b+1]; n++)
598 blcc[n] = blmf[n]*iprod(r[b], r[blbnb[n]]);
601 /* Together: 23*ncons + 6*nrtot flops */
603 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
604 /* nrec*(ncons+2*nrtot) flops */
606 if (econq == econqDeriv_FlexCon)
608 /* We only want to constraint the flexible constraints,
609 * so we mask out the normal ones by setting sol to 0.
611 for (b = b0; b < b1; b++)
613 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
620 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
621 for (b = b0; b < b1; b++)
626 /* When constraining forces, we should not use mass weighting,
627 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
629 lincs_update_atoms(lincsd, th, 1.0, sol, r,
630 (econq != econqForce) ? invmass : NULL, fp);
637 for (b = b0; b < b1; b++)
639 dhdlambda -= sol[b]*lincsd->ddist[b];
642 lincsd->task[th].dhdlambda = dhdlambda;
647 /* Constraint virial,
648 * determines sum r_bond x delta f,
649 * where delta f is the constraint correction
650 * of the quantity that is being constrained.
652 for (b = b0; b < b1; b++)
657 mvb = lincsd->bllen[b]*sol[b];
658 for (i = 0; i < DIM; i++)
661 for (j = 0; j < DIM; j++)
663 rmdf[i][j] += tmp1*r[b][j];
666 } /* 23 ncons flops */
670 #if GMX_SIMD_HAVE_REAL
671 /* Calculate the constraint distance vectors r to project on from x.
672 * Determine the right-hand side of the matrix equation using coordinates xp.
674 static void gmx_simdcall
675 calc_dr_x_xp_simd(int b0,
678 const rvec * gmx_restrict x,
679 const rvec * gmx_restrict xp,
680 const real * gmx_restrict bllen,
681 const real * gmx_restrict blc,
682 const real * pbc_simd,
683 rvec * gmx_restrict r,
684 real * gmx_restrict rhs,
685 real * gmx_restrict sol)
687 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
688 GMX_ALIGNED(int, GMX_SIMD_REAL_WIDTH) offset2[GMX_SIMD_REAL_WIDTH];
690 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
695 for (int bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
697 SimdReal x0_S, y0_S, z0_S;
698 SimdReal x1_S, y1_S, z1_S;
699 SimdReal rx_S, ry_S, rz_S, n2_S, il_S;
700 SimdReal rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
701 GMX_ALIGNED(int, GMX_SIMD_REAL_WIDTH) offset0[GMX_SIMD_REAL_WIDTH];
702 GMX_ALIGNED(int, GMX_SIMD_REAL_WIDTH) offset1[GMX_SIMD_REAL_WIDTH];
704 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
706 offset0[i] = bla[bs*2 + i*2];
707 offset1[i] = bla[bs*2 + i*2 + 1];
710 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
711 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
716 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
718 n2_S = norm2(rx_S, ry_S, rz_S);
719 il_S = invsqrt(n2_S);
725 transposeScatterStoreU<3>(reinterpret_cast<real *>(r + bs), offset2, rx_S, ry_S, rz_S);
727 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(xp), offset0, &x0_S, &y0_S, &z0_S);
728 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(xp), offset1, &x1_S, &y1_S, &z1_S);
733 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
735 ip_S = iprod(rx_S, ry_S, rz_S, rxp_S, ryp_S, rzp_S);
737 rhs_S = load(blc + bs) * (ip_S - load(bllen + bs));
739 store(rhs + bs, rhs_S);
740 store(sol + bs, rhs_S);
743 #endif // GMX_SIMD_HAVE_REAL
745 /* Determine the distances and right-hand side for the next iteration */
746 static void calc_dist_iter(int b0,
749 const rvec * gmx_restrict xp,
750 const real * gmx_restrict bllen,
751 const real * gmx_restrict blc,
754 real * gmx_restrict rhs,
755 real * gmx_restrict sol,
760 for (b = b0; b < b1; b++)
762 real len, len2, dlen2, mvb;
768 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
772 rvec_sub(xp[bla[2*b]], xp[bla[2*b+1]], dx);
775 dlen2 = 2*len2 - norm2(dx);
776 if (dlen2 < wfac*len2)
778 /* not race free - see detailed comment in caller */
783 mvb = blc[b]*(len - dlen2*gmx::invsqrt(dlen2));
791 } /* 20*ncons flops */
794 #if GMX_SIMD_HAVE_REAL
795 /* As the function above, but using SIMD intrinsics */
796 static void gmx_simdcall
797 calc_dist_iter_simd(int b0,
800 const rvec * gmx_restrict x,
801 const real * gmx_restrict bllen,
802 const real * gmx_restrict blc,
803 const real * pbc_simd,
805 real * gmx_restrict rhs,
806 real * gmx_restrict sol,
809 SimdReal min_S(GMX_REAL_MIN);
811 SimdReal wfac_S(wfac);
816 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
818 /* Initialize all to FALSE */
819 warn_S = (two_S < setZero());
821 for (bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
823 SimdReal x0_S, y0_S, z0_S;
824 SimdReal x1_S, y1_S, z1_S;
825 SimdReal rx_S, ry_S, rz_S, n2_S;
826 SimdReal len_S, len2_S, dlen2_S, lc_S, blc_S;
827 GMX_ALIGNED(int, GMX_SIMD_REAL_WIDTH) offset0[GMX_SIMD_REAL_WIDTH];
828 GMX_ALIGNED(int, GMX_SIMD_REAL_WIDTH) offset1[GMX_SIMD_REAL_WIDTH];
830 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
832 offset0[i] = bla[bs*2 + i*2];
833 offset1[i] = bla[bs*2 + i*2 + 1];
836 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset0, &x0_S, &y0_S, &z0_S);
837 gatherLoadUTranspose<3>(reinterpret_cast<const real *>(x), offset1, &x1_S, &y1_S, &z1_S);
842 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
844 n2_S = norm2(rx_S, ry_S, rz_S);
846 len_S = load(bllen + bs);
847 len2_S = len_S * len_S;
849 dlen2_S = fms(two_S, len2_S, n2_S);
851 warn_S = warn_S || (dlen2_S < (wfac_S * len2_S));
853 /* Avoid 1/0 by taking the max with REAL_MIN.
854 * Note: when dlen2 is close to zero (90 degree constraint rotation),
855 * the accuracy of the algorithm is no longer relevant.
857 dlen2_S = max(dlen2_S, min_S);
859 lc_S = fnma(dlen2_S, invsqrt(dlen2_S), len_S);
861 blc_S = load(blc + bs);
865 store(rhs + bs, lc_S);
866 store(sol + bs, lc_S);
874 #endif // GMX_SIMD_HAVE_REAL
876 static void do_lincs(rvec *x, rvec *xp, matrix box, t_pbc *pbc,
877 struct gmx_lincsdata *lincsd, int th,
881 real wangle, gmx_bool *bWarn,
882 real invdt, rvec * gmx_restrict v,
883 gmx_bool bCalcVir, tensor vir_r_m_dr)
885 int b0, b1, b, i, j, n, iter;
886 int *bla, *blnr, *blbnb;
888 real *blc, *blmf, *bllen, *blcc, *rhs1, *rhs2, *sol, *blc_sol, *mlambda;
891 b0 = lincsd->task[th].b0;
892 b1 = lincsd->task[th].b1;
897 blbnb = lincsd->blbnb;
900 bllen = lincsd->bllen;
901 blcc = lincsd->tmpncc;
905 blc_sol = lincsd->tmp4;
906 mlambda = lincsd->mlambda;
907 nlocat = lincsd->nlocat;
909 #if GMX_SIMD_HAVE_REAL
911 /* This SIMD code does the same as the plain-C code after the #else.
912 * The only difference is that we always call pbc code, as with SIMD
913 * the overhead of pbc computation (when not needed) is small.
915 GMX_ALIGNED(real, GMX_SIMD_REAL_WIDTH) pbc_simd[9*GMX_SIMD_REAL_WIDTH];
917 /* Convert the pbc struct for SIMD */
918 set_pbc_simd(pbc, pbc_simd);
920 /* Compute normalized x i-j vectors, store in r.
921 * Compute the inner product of r and xp i-j and store in rhs1.
923 calc_dr_x_xp_simd(b0, b1, bla, x, xp, bllen, blc,
927 #else // GMX_SIMD_HAVE_REAL
931 /* Compute normalized i-j vectors */
932 for (b = b0; b < b1; b++)
937 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
940 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
941 mvb = blc[b]*(iprod(r[b], dx) - bllen[b]);
948 /* Compute normalized i-j vectors */
949 for (b = b0; b < b1; b++)
951 real tmp0, tmp1, tmp2, rlen, mvb;
955 tmp0 = x[i][0] - x[j][0];
956 tmp1 = x[i][1] - x[j][1];
957 tmp2 = x[i][2] - x[j][2];
958 rlen = gmx::invsqrt(tmp0*tmp0 + tmp1*tmp1 + tmp2*tmp2);
966 mvb = blc[b]*(r[b][0]*(xp[i][0] - xp[j][0]) +
967 r[b][1]*(xp[i][1] - xp[j][1]) +
968 r[b][2]*(xp[i][2] - xp[j][2]) - bllen[b]);
973 /* Together: 26*ncons + 6*nrtot flops */
976 #endif // GMX_SIMD_HAVE_REAL
978 if (lincsd->bTaskDep)
980 /* We need a barrier, since the matrix construction below
981 * can access entries in r of other threads.
986 /* Construct the (sparse) LINCS matrix */
987 for (b = b0; b < b1; b++)
989 for (n = blnr[b]; n < blnr[b+1]; n++)
991 blcc[n] = blmf[n]*iprod(r[b], r[blbnb[n]]);
994 /* Together: 26*ncons + 6*nrtot flops */
996 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
997 /* nrec*(ncons+2*nrtot) flops */
999 #if GMX_SIMD_HAVE_REAL
1000 for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1002 SimdReal t1 = load(blc + b);
1003 SimdReal t2 = load(sol + b);
1004 store(mlambda + b, t1 * t2);
1007 for (b = b0; b < b1; b++)
1009 mlambda[b] = blc[b]*sol[b];
1011 #endif // GMX_SIMD_HAVE_REAL
1013 /* Update the coordinates */
1014 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1017 ******** Correction for centripetal effects ********
1022 wfac = cos(DEG2RAD*wangle);
1025 for (iter = 0; iter < lincsd->nIter; iter++)
1027 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
1032 /* Communicate the corrected non-local coordinates */
1033 if (DOMAINDECOMP(cr))
1035 dd_move_x_constraints(cr->dd, box, xp, NULL, FALSE);
1040 else if (lincsd->bTaskDep)
1045 #if GMX_SIMD_HAVE_REAL
1046 calc_dist_iter_simd(b0, b1, bla, xp, bllen, blc, pbc_simd, wfac,
1049 calc_dist_iter(b0, b1, bla, xp, bllen, blc, pbc, wfac,
1051 /* 20*ncons flops */
1052 #endif // GMX_SIMD_HAVE_REAL
1054 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1055 /* nrec*(ncons+2*nrtot) flops */
1057 #if GMX_SIMD_HAVE_REAL
1058 for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1060 SimdReal t1 = load(blc + b);
1061 SimdReal t2 = load(sol + b);
1062 SimdReal mvb = t1 * t2;
1063 store(blc_sol + b, mvb);
1064 store(mlambda + b, load(mlambda + b) + mvb);
1067 for (b = b0; b < b1; b++)
1071 mvb = blc[b]*sol[b];
1075 #endif // GMX_SIMD_HAVE_REAL
1077 /* Update the coordinates */
1078 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1080 /* nit*ncons*(37+9*nrec) flops */
1084 /* Update the velocities */
1085 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1086 /* 16 ncons flops */
1089 if (nlocat != NULL && (bCalcDHDL || bCalcVir))
1091 if (lincsd->bTaskDep)
1093 /* In lincs_update_atoms threads might cross-read mlambda */
1097 /* Only account for local atoms */
1098 for (b = b0; b < b1; b++)
1100 mlambda[b] *= 0.5*nlocat[b];
1109 for (b = b0; b < b1; b++)
1111 /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1112 * later after the contributions are reduced over the threads.
1114 dhdl -= lincsd->mlambda[b]*lincsd->ddist[b];
1116 lincsd->task[th].dhdlambda = dhdl;
1121 /* Constraint virial */
1122 for (b = b0; b < b1; b++)
1126 tmp0 = -bllen[b]*mlambda[b];
1127 for (i = 0; i < DIM; i++)
1129 tmp1 = tmp0*r[b][i];
1130 for (j = 0; j < DIM; j++)
1132 vir_r_m_dr[i][j] -= tmp1*r[b][j];
1135 } /* 22 ncons flops */
1139 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1140 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1142 * (26+nrec)*ncons + (6+2*nrec)*nrtot
1143 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1145 * (63+nrec)*ncons + (6+4*nrec)*nrtot
1149 /* Sets the elements in the LINCS matrix for task li_task */
1150 static void set_lincs_matrix_task(struct gmx_lincsdata *li,
1151 lincs_task_t *li_task,
1152 const real *invmass,
1157 /* Construct the coupling coefficient matrix blmf */
1158 li_task->ntriangle = 0;
1160 for (i = li_task->b0; i < li_task->b1; i++)
1165 a2 = li->bla[2*i+1];
1166 for (n = li->blnr[i]; (n < li->blnr[i+1]); n++)
1168 int k, sign, center, end;
1172 /* If we are using multiple, independent tasks for LINCS,
1173 * the calls to check_assign_connected should have
1174 * put all connected constraints in our task.
1176 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1178 if (a1 == li->bla[2*k] || a2 == li->bla[2*k+1])
1186 if (a1 == li->bla[2*k] || a1 == li->bla[2*k+1])
1196 li->blmf[n] = sign*invmass[center]*li->blc[i]*li->blc[k];
1197 li->blmf1[n] = sign*0.5;
1198 if (li->ncg_triangle > 0)
1202 /* Look for constraint triangles */
1203 for (nk = li->blnr[k]; (nk < li->blnr[k+1]); nk++)
1206 if (kk != i && kk != k &&
1207 (li->bla[2*kk] == end || li->bla[2*kk+1] == end))
1209 /* If we are using multiple tasks for LINCS,
1210 * the calls to check_assign_triangle should have
1211 * put all constraints in the triangle in our task.
1213 assert(k >= li_task->b0 && k < li_task->b1);
1214 assert(kk >= li_task->b0 && kk < li_task->b1);
1216 if (li_task->ntriangle == 0 ||
1217 li_task->triangle[li_task->ntriangle - 1] < i)
1219 /* Add this constraint to the triangle list */
1220 li_task->triangle[li_task->ntriangle] = i;
1221 li_task->tri_bits[li_task->ntriangle] = 0;
1222 li_task->ntriangle++;
1223 if (li->blnr[i+1] - li->blnr[i] > static_cast<int>(sizeof(li_task->tri_bits[0])*8 - 1))
1225 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %d allowed for constraints participating in triangles",
1226 li->blnr[i+1] - li->blnr[i],
1227 sizeof(li_task->tri_bits[0])*8-1);
1230 li_task->tri_bits[li_task->ntriangle-1] |= (1 << (n - li->blnr[i]));
1239 /* Sets the elements in the LINCS matrix */
1240 void set_lincs_matrix(struct gmx_lincsdata *li, real *invmass, real lambda)
1243 const real invsqrt2 = 0.7071067811865475244;
1245 for (i = 0; (i < li->nc); i++)
1250 a2 = li->bla[2*i+1];
1251 li->blc[i] = gmx::invsqrt(invmass[a1] + invmass[a2]);
1252 li->blc1[i] = invsqrt2;
1255 /* Construct the coupling coefficient matrix blmf */
1256 int th, ntriangle = 0, ncc_triangle = 0;
1257 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle) num_threads(li->ntask) schedule(static)
1258 for (th = 0; th < li->ntask; th++)
1262 set_lincs_matrix_task(li, &li->task[th], invmass, &ncc_triangle);
1263 ntriangle = li->task[th].ntriangle;
1265 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1267 li->ntriangle = ntriangle;
1268 li->ncc_triangle = ncc_triangle;
1272 fprintf(debug, "Of the %d constraints %d participate in triangles\n",
1273 li->nc, li->ntriangle);
1274 fprintf(debug, "There are %d couplings of which %d in triangles\n",
1275 li->ncc, li->ncc_triangle);
1279 * so we know with which lambda value the masses have been set.
1281 li->matlam = lambda;
1284 static int count_triangle_constraints(const t_ilist *ilist,
1285 const t_blocka *at2con)
1287 int ncon1, ncon_tot;
1288 int c0, a00, a01, n1, c1, a10, a11, ac1, n2, c2, a20, a21;
1291 t_iatom *ia1, *ia2, *iap;
1293 ncon1 = ilist[F_CONSTR].nr/3;
1294 ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
1296 ia1 = ilist[F_CONSTR].iatoms;
1297 ia2 = ilist[F_CONSTRNC].iatoms;
1300 for (c0 = 0; c0 < ncon_tot; c0++)
1303 iap = constr_iatomptr(ncon1, ia1, ia2, c0);
1306 for (n1 = at2con->index[a01]; n1 < at2con->index[a01+1]; n1++)
1311 iap = constr_iatomptr(ncon1, ia1, ia2, c1);
1322 for (n2 = at2con->index[ac1]; n2 < at2con->index[ac1+1]; n2++)
1325 if (c2 != c0 && c2 != c1)
1327 iap = constr_iatomptr(ncon1, ia1, ia2, c2);
1330 if (a20 == a00 || a21 == a00)
1344 return ncon_triangle;
1347 static gmx_bool more_than_two_sequential_constraints(const t_ilist *ilist,
1348 const t_blocka *at2con)
1350 t_iatom *ia1, *ia2, *iap;
1351 int ncon1, ncon_tot, c;
1353 gmx_bool bMoreThanTwoSequentialConstraints;
1355 ncon1 = ilist[F_CONSTR].nr/3;
1356 ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
1358 ia1 = ilist[F_CONSTR].iatoms;
1359 ia2 = ilist[F_CONSTRNC].iatoms;
1361 bMoreThanTwoSequentialConstraints = FALSE;
1362 for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
1364 iap = constr_iatomptr(ncon1, ia1, ia2, c);
1367 /* Check if this constraint has constraints connected at both atoms */
1368 if (at2con->index[a1+1] - at2con->index[a1] > 1 &&
1369 at2con->index[a2+1] - at2con->index[a2] > 1)
1371 bMoreThanTwoSequentialConstraints = TRUE;
1375 return bMoreThanTwoSequentialConstraints;
1378 static int int_comp(const void *a, const void *b)
1380 return (*(int *)a) - (*(int *)b);
1383 gmx_lincsdata_t init_lincs(FILE *fplog, const gmx_mtop_t *mtop,
1384 int nflexcon_global, const t_blocka *at2con,
1385 gmx_bool bPLINCS, int nIter, int nProjOrder)
1387 struct gmx_lincsdata *li;
1389 gmx_moltype_t *molt;
1390 gmx_bool bMoreThanTwoSeq;
1394 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
1395 bPLINCS ? " Parallel" : "");
1401 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1402 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1403 li->ncg_flex = nflexcon_global;
1406 li->nOrder = nProjOrder;
1408 li->max_connect = 0;
1409 for (mt = 0; mt < mtop->nmoltype; mt++)
1413 molt = &mtop->moltype[mt];
1414 for (a = 0; a < molt->atoms.nr; a++)
1416 li->max_connect = std::max(li->max_connect,
1417 at2con[mt].index[a + 1] - at2con[mt].index[a]);
1421 li->ncg_triangle = 0;
1422 bMoreThanTwoSeq = FALSE;
1423 for (mb = 0; mb < mtop->nmolblock; mb++)
1425 molt = &mtop->moltype[mtop->molblock[mb].type];
1428 mtop->molblock[mb].nmol*
1429 count_triangle_constraints(molt->ilist,
1430 &at2con[mtop->molblock[mb].type]);
1432 if (!bMoreThanTwoSeq &&
1433 more_than_two_sequential_constraints(molt->ilist, &at2con[mtop->molblock[mb].type]))
1435 bMoreThanTwoSeq = TRUE;
1439 /* Check if we need to communicate not only before LINCS,
1440 * but also before each iteration.
1441 * The check for only two sequential constraints is only
1442 * useful for the common case of H-bond only constraints.
1443 * With more effort we could also make it useful for small
1444 * molecules with nr. sequential constraints <= nOrder-1.
1446 li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1448 if (debug && bPLINCS)
1450 fprintf(debug, "PLINCS communication before each iteration: %d\n",
1454 /* LINCS can run on any number of threads.
1455 * Currently the number is fixed for the whole simulation,
1456 * but it could be set in set_lincs().
1457 * The current constraint to task assignment code can create independent
1458 * tasks only when not more than two constraints are connected sequentially.
1460 li->ntask = gmx_omp_nthreads_get(emntLINCS);
1461 li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1464 fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n",
1465 li->ntask, li->bTaskDep ? "" : "in");
1473 /* Allocate an extra elements for "task-overlap" constraints */
1474 snew(li->task, li->ntask + 1);
1477 if (bPLINCS || li->ncg_triangle > 0)
1479 please_cite(fplog, "Hess2008a");
1483 please_cite(fplog, "Hess97a");
1488 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1491 fprintf(fplog, "There are inter charge-group constraints,\n"
1492 "will communicate selected coordinates each lincs iteration\n");
1494 if (li->ncg_triangle > 0)
1497 "%d constraints are involved in constraint triangles,\n"
1498 "will apply an additional matrix expansion of order %d for couplings\n"
1499 "between constraints inside triangles\n",
1500 li->ncg_triangle, li->nOrder);
1507 /* Sets up the work division over the threads */
1508 static void lincs_thread_setup(struct gmx_lincsdata *li, int natoms)
1515 if (natoms > li->atf_nalloc)
1517 li->atf_nalloc = over_alloc_large(natoms);
1518 srenew(li->atf, li->atf_nalloc);
1522 /* Clear the atom flags */
1523 for (a = 0; a < natoms; a++)
1525 bitmask_clear(&atf[a]);
1528 if (li->ntask > BITMASK_SIZE)
1530 gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1533 for (th = 0; th < li->ntask; th++)
1535 lincs_task_t *li_task;
1538 li_task = &li->task[th];
1540 /* For each atom set a flag for constraints from each */
1541 for (b = li_task->b0; b < li_task->b1; b++)
1543 bitmask_set_bit(&atf[li->bla[b*2 ]], th);
1544 bitmask_set_bit(&atf[li->bla[b*2 + 1]], th);
1548 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1549 for (th = 0; th < li->ntask; th++)
1553 lincs_task_t *li_task;
1557 li_task = &li->task[th];
1559 if (li_task->b1 - li_task->b0 > li_task->ind_nalloc)
1561 li_task->ind_nalloc = over_alloc_large(li_task->b1-li_task->b0);
1562 srenew(li_task->ind, li_task->ind_nalloc);
1563 srenew(li_task->ind_r, li_task->ind_nalloc);
1566 bitmask_init_low_bits(&mask, th);
1569 li_task->nind_r = 0;
1570 for (b = li_task->b0; b < li_task->b1; b++)
1572 /* We let the constraint with the lowest thread index
1573 * operate on atoms with constraints from multiple threads.
1575 if (bitmask_is_disjoint(atf[li->bla[b*2]], mask) &&
1576 bitmask_is_disjoint(atf[li->bla[b*2+1]], mask))
1578 /* Add the constraint to the local atom update index */
1579 li_task->ind[li_task->nind++] = b;
1583 /* Add the constraint to the rest block */
1584 li_task->ind_r[li_task->nind_r++] = b;
1588 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1591 /* We need to copy all constraints which have not be assigned
1592 * to a thread to a separate list which will be handled by one thread.
1594 li_m = &li->task[li->ntask];
1597 for (th = 0; th < li->ntask; th++)
1599 lincs_task_t *li_task;
1602 li_task = &li->task[th];
1604 if (li_m->nind + li_task->nind_r > li_m->ind_nalloc)
1606 li_m->ind_nalloc = over_alloc_large(li_m->nind+li_task->nind_r);
1607 srenew(li_m->ind, li_m->ind_nalloc);
1610 for (b = 0; b < li_task->nind_r; b++)
1612 li_m->ind[li_m->nind++] = li_task->ind_r[b];
1617 fprintf(debug, "LINCS thread %d: %d constraints\n",
1624 fprintf(debug, "LINCS thread r: %d constraints\n",
1629 /* There is no realloc with alignment, so here we make one for reals.
1630 * Note that this function does not preserve the contents of the memory.
1632 static void resize_real_aligned(real **ptr, int nelem)
1634 sfree_aligned(*ptr);
1635 snew_aligned(*ptr, nelem, align_bytes);
1638 static void assign_constraint(struct gmx_lincsdata *li,
1639 int constraint_index,
1641 real lenA, real lenB,
1642 const t_blocka *at2con)
1648 /* Make an mapping of local topology constraint index to LINCS index */
1649 li->con_index[constraint_index] = con;
1651 li->bllen0[con] = lenA;
1652 li->ddist[con] = lenB - lenA;
1653 /* Set the length to the topology A length */
1654 li->bllen[con] = lenA;
1655 li->bla[2*con] = a1;
1656 li->bla[2*con+1] = a2;
1658 /* Make space in the constraint connection matrix for constraints
1659 * connected to both end of the current constraint.
1662 at2con->index[a1 + 1] - at2con->index[a1] - 1 +
1663 at2con->index[a2 + 1] - at2con->index[a2] - 1;
1665 li->blnr[con + 1] = li->ncc;
1667 /* Increase the constraint count */
1671 /* Check if constraint with topology index constraint_index is connected
1672 * to other constraints, and if so add those connected constraints to our task.
1674 static void check_assign_connected(struct gmx_lincsdata *li,
1675 const t_iatom *iatom,
1679 const t_blocka *at2con)
1681 /* Currently this function only supports constraint groups
1682 * in which all constraints share at least one atom
1683 * (e.g. H-bond constraints).
1684 * Check both ends of the current constraint for
1685 * connected constraints. We need to assign those
1690 for (end = 0; end < 2; end++)
1694 a = (end == 0 ? a1 : a2);
1696 for (k = at2con->index[a]; k < at2con->index[a + 1]; k++)
1701 /* Check if constraint cc has not yet been assigned */
1702 if (li->con_index[cc] == -1)
1708 lenA = idef->iparams[type].constr.dA;
1709 lenB = idef->iparams[type].constr.dB;
1711 if (bDynamics || lenA != 0 || lenB != 0)
1713 assign_constraint(li, cc, iatom[3*cc + 1], iatom[3*cc + 2], lenA, lenB, at2con);
1720 /* Check if constraint with topology index constraint_index is involved
1721 * in a constraint triangle, and if so add the other two constraints
1722 * in the triangle to our task.
1724 static void check_assign_triangle(struct gmx_lincsdata *li,
1725 const t_iatom *iatom,
1728 int constraint_index,
1730 const t_blocka *at2con)
1732 int nca, cc[32], ca[32], k;
1733 int c_triangle[2] = { -1, -1 };
1736 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1741 if (c != constraint_index)
1745 aa1 = iatom[c*3 + 1];
1746 aa2 = iatom[c*3 + 2];
1762 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1767 if (c != constraint_index)
1771 aa1 = iatom[c*3 + 1];
1772 aa2 = iatom[c*3 + 2];
1775 for (i = 0; i < nca; i++)
1779 c_triangle[0] = cc[i];
1786 for (i = 0; i < nca; i++)
1790 c_triangle[0] = cc[i];
1798 if (c_triangle[0] >= 0)
1802 for (end = 0; end < 2; end++)
1804 /* Check if constraint c_triangle[end] has not yet been assigned */
1805 if (li->con_index[c_triangle[end]] == -1)
1810 i = c_triangle[end]*3;
1812 lenA = idef->iparams[type].constr.dA;
1813 lenB = idef->iparams[type].constr.dB;
1815 if (bDynamics || lenA != 0 || lenB != 0)
1817 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1824 static void set_matrix_indices(struct gmx_lincsdata *li,
1825 const lincs_task_t *li_task,
1826 const t_blocka *at2con,
1827 gmx_bool bSortMatrix)
1831 for (b = li_task->b0; b < li_task->b1; b++)
1836 a2 = li->bla[b*2 + 1];
1839 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1843 concon = li->con_index[at2con->a[k]];
1846 li->blbnb[i++] = concon;
1849 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1853 concon = li->con_index[at2con->a[k]];
1856 li->blbnb[i++] = concon;
1862 /* Order the blbnb matrix to optimize memory access */
1863 qsort(&(li->blbnb[li->blnr[b]]), li->blnr[b + 1] - li->blnr[b],
1864 sizeof(li->blbnb[0]), int_comp);
1869 void set_lincs(const t_idef *idef,
1870 const t_mdatoms *md,
1873 struct gmx_lincsdata *li)
1875 int natoms, nflexcon;
1878 int i, ncc_alloc_old, ncon_tot;
1883 /* Zero the thread index ranges.
1884 * Otherwise without local constraints we could return with old ranges.
1886 for (i = 0; i < li->ntask; i++)
1890 li->task[i].nind = 0;
1894 li->task[li->ntask].nind = 0;
1897 /* This is the local topology, so there are only F_CONSTR constraints */
1898 if (idef->il[F_CONSTR].nr == 0)
1900 /* There are no constraints,
1901 * we do not need to fill any data structures.
1908 fprintf(debug, "Building the LINCS connectivity\n");
1911 if (DOMAINDECOMP(cr))
1913 if (cr->dd->constraints)
1917 dd_get_constraint_range(cr->dd, &start, &natoms);
1921 natoms = cr->dd->nat_home;
1926 natoms = md->homenr;
1928 at2con = make_at2con(0, natoms, idef->il, idef->iparams, bDynamics,
1931 ncon_tot = idef->il[F_CONSTR].nr/3;
1933 /* Ensure we have enough padding for aligned loads for each thread */
1934 if (ncon_tot + li->ntask*simd_width > li->nc_alloc || li->nc_alloc == 0)
1936 li->nc_alloc = over_alloc_dd(ncon_tot + li->ntask*simd_width);
1937 srenew(li->con_index, li->nc_alloc);
1938 resize_real_aligned(&li->bllen0, li->nc_alloc);
1939 resize_real_aligned(&li->ddist, li->nc_alloc);
1940 srenew(li->bla, 2*li->nc_alloc);
1941 resize_real_aligned(&li->blc, li->nc_alloc);
1942 resize_real_aligned(&li->blc1, li->nc_alloc);
1943 srenew(li->blnr, li->nc_alloc + 1);
1944 resize_real_aligned(&li->bllen, li->nc_alloc);
1945 srenew(li->tmpv, li->nc_alloc);
1946 if (DOMAINDECOMP(cr))
1948 srenew(li->nlocat, li->nc_alloc);
1950 resize_real_aligned(&li->tmp1, li->nc_alloc);
1951 resize_real_aligned(&li->tmp2, li->nc_alloc);
1952 resize_real_aligned(&li->tmp3, li->nc_alloc);
1953 resize_real_aligned(&li->tmp4, li->nc_alloc);
1954 resize_real_aligned(&li->mlambda, li->nc_alloc);
1957 iatom = idef->il[F_CONSTR].iatoms;
1959 ncc_alloc_old = li->ncc_alloc;
1960 li->blnr[0] = li->ncc;
1962 /* Assign the constraints for li->ntask LINCS tasks.
1963 * We target a uniform distribution of constraints over the tasks.
1964 * Note that when flexible constraints are present, but are removed here
1965 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
1966 * happen during normal MD, that's ok.
1968 int ncon_assign, ncon_target, con, th;
1970 /* Determine the number of constraints we need to assign here */
1971 ncon_assign = ncon_tot;
1974 /* With energy minimization, flexible constraints are ignored
1975 * (and thus minimized, as they should be).
1977 ncon_assign -= nflexcon;
1980 /* Set the target constraint count per task to exactly uniform,
1981 * this might be overridden below.
1983 ncon_target = (ncon_assign + li->ntask - 1)/li->ntask;
1985 /* Mark all constraints as unassigned by setting their index to -1 */
1986 for (con = 0; con < ncon_tot; con++)
1988 li->con_index[con] = -1;
1992 for (th = 0; th < li->ntask; th++)
1994 lincs_task_t *li_task;
1996 li_task = &li->task[th];
1998 #if GMX_SIMD_HAVE_REAL
1999 /* With indepedent tasks we likely have H-bond constraints or constraint
2000 * pairs. The connected constraints will be pulled into the task, so the
2001 * constraints per task will often exceed ncon_target.
2002 * Triangle constraints can also increase the count, but there are
2003 * relatively few of those, so we usually expect to get ncon_target.
2007 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
2008 * since otherwise a lot of operations can be wasted.
2009 * There are several ways to round here, we choose the one
2010 * that alternates block sizes, which helps with Intel HT.
2012 ncon_target = ((ncon_assign*(th + 1))/li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1) & ~(GMX_SIMD_REAL_WIDTH - 1);
2014 #endif // GMX_SIMD==2 && GMX_SIMD_HAVE_REAL
2016 /* Continue filling the arrays where we left off with the previous task,
2017 * including padding for SIMD.
2019 li_task->b0 = li->nc;
2021 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2023 if (li->con_index[con] == -1)
2028 type = iatom[3*con];
2029 a1 = iatom[3*con + 1];
2030 a2 = iatom[3*con + 2];
2031 lenA = idef->iparams[type].constr.dA;
2032 lenB = idef->iparams[type].constr.dB;
2033 /* Skip the flexible constraints when not doing dynamics */
2034 if (bDynamics || lenA != 0 || lenB != 0)
2036 assign_constraint(li, con, a1, a2, lenA, lenB, &at2con);
2038 if (li->ntask > 1 && !li->bTaskDep)
2040 /* We can generate independent tasks. Check if we
2041 * need to assign connected constraints to our task.
2043 check_assign_connected(li, iatom, idef, bDynamics,
2046 if (li->ntask > 1 && li->ncg_triangle > 0)
2048 /* Ensure constraints in one triangle are assigned
2051 check_assign_triangle(li, iatom, idef, bDynamics,
2052 con, a1, a2, &at2con);
2060 li_task->b1 = li->nc;
2064 /* Copy the last atom pair indices and lengths for constraints
2065 * up to a multiple of simd_width, such that we can do all
2066 * SIMD operations without having to worry about end effects.
2070 li->nc = ((li_task->b1 + simd_width - 1)/simd_width)*simd_width;
2071 last = li_task->b1 - 1;
2072 for (i = li_task->b1; i < li->nc; i++)
2074 li->bla[i*2 ] = li->bla[last*2 ];
2075 li->bla[i*2 + 1] = li->bla[last*2 + 1];
2076 li->bllen0[i] = li->bllen0[last];
2077 li->ddist[i] = li->ddist[last];
2078 li->bllen[i] = li->bllen[last];
2079 li->blnr[i + 1] = li->blnr[last + 1];
2083 /* Keep track of how many constraints we assigned */
2084 li->nc_real += li_task->b1 - li_task->b0;
2088 fprintf(debug, "LINCS task %d constraints %d - %d\n",
2089 th, li_task->b0, li_task->b1);
2093 assert(li->nc_real == ncon_assign);
2095 gmx_bool bSortMatrix;
2097 /* Without DD we order the blbnb matrix to optimize memory access.
2098 * With DD the overhead of sorting is more than the gain during access.
2100 bSortMatrix = !DOMAINDECOMP(cr);
2102 if (li->ncc > li->ncc_alloc)
2104 li->ncc_alloc = over_alloc_small(li->ncc);
2105 srenew(li->blbnb, li->ncc_alloc);
2108 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2109 for (th = 0; th < li->ntask; th++)
2113 lincs_task_t *li_task;
2115 li_task = &li->task[th];
2117 if (li->ncg_triangle > 0 &&
2118 li_task->b1 - li_task->b0 > li_task->tri_alloc)
2120 /* This is allocating too much, but it is difficult to improve */
2121 li_task->tri_alloc = over_alloc_dd(li_task->b1 - li_task->b0);
2122 srenew(li_task->triangle, li_task->tri_alloc);
2123 srenew(li_task->tri_bits, li_task->tri_alloc);
2126 set_matrix_indices(li, li_task, &at2con, bSortMatrix);
2128 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2131 done_blocka(&at2con);
2135 /* Since the matrix is static, we should free some memory */
2136 li->ncc_alloc = li->ncc;
2137 srenew(li->blbnb, li->ncc_alloc);
2140 if (li->ncc_alloc > ncc_alloc_old)
2142 srenew(li->blmf, li->ncc_alloc);
2143 srenew(li->blmf1, li->ncc_alloc);
2144 srenew(li->tmpncc, li->ncc_alloc);
2147 if (DOMAINDECOMP(cr) && dd_constraints_nlocalatoms(cr->dd) != NULL)
2151 nlocat_dd = dd_constraints_nlocalatoms(cr->dd);
2153 /* Convert nlocat from local topology to LINCS constraint indexing */
2154 for (con = 0; con < ncon_tot; con++)
2156 li->nlocat[li->con_index[con]] = nlocat_dd[con];
2166 fprintf(debug, "Number of constraints is %d, padded %d, couplings %d\n",
2167 li->nc_real, li->nc, li->ncc);
2172 lincs_thread_setup(li, md->nr);
2175 set_lincs_matrix(li, md->invmass, md->lambda);
2178 static void lincs_warning(FILE *fplog,
2179 gmx_domdec_t *dd, rvec *x, rvec *xprime, t_pbc *pbc,
2180 int ncons, int *bla, real *bllen, real wangle,
2181 int maxwarn, int *warncount)
2185 real wfac, d0, d1, cosine;
2188 wfac = cos(DEG2RAD*wangle);
2190 sprintf(buf, "bonds that rotated more than %g degrees:\n"
2191 " atom 1 atom 2 angle previous, current, constraint length\n",
2193 fprintf(stderr, "%s", buf);
2196 fprintf(fplog, "%s", buf);
2199 for (b = 0; b < ncons; b++)
2205 pbc_dx_aiuc(pbc, x[i], x[j], v0);
2206 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
2210 rvec_sub(x[i], x[j], v0);
2211 rvec_sub(xprime[i], xprime[j], v1);
2215 cosine = iprod(v0, v1)/(d0*d1);
2218 sprintf(buf, " %6d %6d %5.1f %8.4f %8.4f %8.4f\n",
2219 ddglatnr(dd, i), ddglatnr(dd, j),
2220 RAD2DEG*acos(cosine), d0, d1, bllen[b]);
2221 fprintf(stderr, "%s", buf);
2224 fprintf(fplog, "%s", buf);
2226 if (!std::isfinite(d1))
2228 gmx_fatal(FARGS, "Bond length not finite.");
2234 if (*warncount > maxwarn)
2236 too_many_constraint_warnings(econtLINCS, *warncount);
2240 static void cconerr(const struct gmx_lincsdata *lincsd,
2241 rvec *x, t_pbc *pbc,
2242 real *ncons_loc, real *ssd, real *max, int *imax)
2244 const int *bla, *nlocat;
2247 int count, im, task;
2250 bllen = lincsd->bllen;
2251 nlocat = lincsd->nlocat;
2257 for (task = 0; task < lincsd->ntask; task++)
2261 for (b = lincsd->task[task].b0; b < lincsd->task[task].b1; b++)
2268 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
2272 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
2275 len = r2*gmx::invsqrt(r2);
2276 d = std::abs(len/bllen[b]-1);
2277 if (d > ma && (nlocat == NULL || nlocat[b]))
2289 ssd2 += nlocat[b]*d*d;
2295 *ncons_loc = (nlocat ? 0.5 : 1)*count;
2296 *ssd = (nlocat ? 0.5 : 1)*ssd2;
2301 gmx_bool constrain_lincs(FILE *fplog, gmx_bool bLog, gmx_bool bEner,
2304 struct gmx_lincsdata *lincsd, t_mdatoms *md,
2306 rvec *x, rvec *xprime, rvec *min_proj,
2307 matrix box, t_pbc *pbc,
2308 real lambda, real *dvdlambda,
2309 real invdt, rvec *v,
2310 gmx_bool bCalcVir, tensor vir_r_m_dr,
2313 int maxwarn, int *warncount)
2316 char buf[STRLEN], buf2[22], buf3[STRLEN];
2318 real ncons_loc, p_ssd, p_max = 0;
2320 gmx_bool bOK, bWarn;
2324 /* This boolean should be set by a flag passed to this routine.
2325 * We can also easily check if any constraint length is changed,
2326 * if not dH/dlambda=0 and we can also set the boolean to FALSE.
2328 bCalcDHDL = (ir->efep != efepNO && dvdlambda != NULL);
2330 if (lincsd->nc == 0 && cr->dd == NULL)
2334 lincsd->rmsd_data[0] = 0;
2335 lincsd->rmsd_data[1] = 0;
2341 if (econq == econqCoord)
2343 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2344 * also with efep!=fepNO.
2346 if (ir->efep != efepNO)
2348 if (md->nMassPerturbed && lincsd->matlam != md->lambda)
2350 set_lincs_matrix(lincsd, md->invmass, md->lambda);
2353 for (i = 0; i < lincsd->nc; i++)
2355 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
2359 if (lincsd->ncg_flex)
2361 /* Set the flexible constraint lengths to the old lengths */
2364 for (i = 0; i < lincsd->nc; i++)
2366 if (lincsd->bllen[i] == 0)
2368 pbc_dx_aiuc(pbc, x[lincsd->bla[2*i]], x[lincsd->bla[2*i+1]], dx);
2369 lincsd->bllen[i] = norm(dx);
2375 for (i = 0; i < lincsd->nc; i++)
2377 if (lincsd->bllen[i] == 0)
2380 std::sqrt(distance2(x[lincsd->bla[2*i]],
2381 x[lincsd->bla[2*i+1]]));
2389 cconerr(lincsd, xprime, pbc,
2390 &ncons_loc, &p_ssd, &p_max, &p_imax);
2393 /* This bWarn var can be updated by multiple threads
2394 * at the same time. But as we only need to detect
2395 * if a warning occured or not, this is not an issue.
2399 /* The OpenMP parallel region of constrain_lincs for coords */
2400 #pragma omp parallel num_threads(lincsd->ntask)
2404 int th = gmx_omp_get_thread_num();
2406 clear_mat(lincsd->task[th].vir_r_m_dr);
2408 do_lincs(x, xprime, box, pbc, lincsd, th,
2411 ir->LincsWarnAngle, &bWarn,
2413 th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2415 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2418 if (bLog && fplog && lincsd->nc > 0)
2420 fprintf(fplog, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2421 fprintf(fplog, " Before LINCS %.6f %.6f %6d %6d\n",
2422 std::sqrt(p_ssd/ncons_loc), p_max,
2423 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2424 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2428 cconerr(lincsd, xprime, pbc,
2429 &ncons_loc, &p_ssd, &p_max, &p_imax);
2430 lincsd->rmsd_data[0] = ncons_loc;
2431 lincsd->rmsd_data[1] = p_ssd;
2435 lincsd->rmsd_data[0] = 0;
2436 lincsd->rmsd_data[1] = 0;
2437 lincsd->rmsd_data[2] = 0;
2439 if (bLog && fplog && lincsd->nc > 0)
2442 " After LINCS %.6f %.6f %6d %6d\n\n",
2443 std::sqrt(p_ssd/ncons_loc), p_max,
2444 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2445 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2452 cconerr(lincsd, xprime, pbc,
2453 &ncons_loc, &p_ssd, &p_max, &p_imax);
2456 sprintf(buf3, " in simulation %d", cr->ms->sim);
2462 sprintf(buf, "\nStep %s, time %g (ps) LINCS WARNING%s\n"
2463 "relative constraint deviation after LINCS:\n"
2464 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2465 gmx_step_str(step, buf2), ir->init_t+step*ir->delta_t,
2467 std::sqrt(p_ssd/ncons_loc), p_max,
2468 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2469 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2472 fprintf(fplog, "%s", buf);
2474 fprintf(stderr, "%s", buf);
2475 lincs_warning(fplog, cr->dd, x, xprime, pbc,
2476 lincsd->nc, lincsd->bla, lincsd->bllen,
2477 ir->LincsWarnAngle, maxwarn, warncount);
2479 bOK = (p_max < 0.5);
2482 if (lincsd->ncg_flex)
2484 for (i = 0; (i < lincsd->nc); i++)
2486 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2488 lincsd->bllen[i] = 0;
2495 /* The OpenMP parallel region of constrain_lincs for derivatives */
2496 #pragma omp parallel num_threads(lincsd->ntask)
2500 int th = gmx_omp_get_thread_num();
2502 do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
2503 md->invmass, econq, bCalcDHDL,
2504 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2506 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2512 /* Reduce the dH/dlambda contributions over the threads */
2517 for (th = 0; th < lincsd->ntask; th++)
2519 dhdlambda += lincsd->task[th].dhdlambda;
2521 if (econq == econqCoord)
2523 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2524 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2525 dhdlambda /= ir->delta_t*ir->delta_t;
2527 *dvdlambda += dhdlambda;
2530 if (bCalcVir && lincsd->ntask > 1)
2532 for (i = 1; i < lincsd->ntask; i++)
2534 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2538 /* count assuming nit=1 */
2539 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2540 inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
2541 if (lincsd->ntriangle > 0)
2543 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
2547 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real*2);
2551 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);