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! */
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/legacyheaders/constr.h"
51 #include "gromacs/legacyheaders/copyrite.h"
52 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
53 #include "gromacs/legacyheaders/mdrun.h"
54 #include "gromacs/legacyheaders/nrnb.h"
55 #include "gromacs/legacyheaders/types/commrec.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/pbcutil/pbc-simd.h"
59 #include "gromacs/pbcutil/pbc.h"
60 #include "gromacs/simd/simd.h"
61 #include "gromacs/simd/simd_math.h"
62 #include "gromacs/simd/vector_operations.h"
63 #include "gromacs/topology/block.h"
64 #include "gromacs/topology/mtop_util.h"
65 #include "gromacs/utility/bitmask.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/gmxomp.h"
68 #include "gromacs/utility/smalloc.h"
70 /* MSVC 2010 produces buggy SIMD PBC code, disable SIMD for MSVC <= 2010 */
71 #if defined GMX_SIMD_HAVE_REAL && !(defined _MSC_VER && _MSC_VER < 1700)
76 int b0; /* first constraint for this task */
77 int b1; /* b1-1 is the last constraint for this task */
78 int ntriangle; /* the number of constraints in triangles */
79 int *triangle; /* the list of triangle constraints */
80 int *tri_bits; /* the bits tell if the matrix element should be used */
81 int tri_alloc; /* allocation size of triangle and tri_bits */
82 int nind; /* number of indices */
83 int *ind; /* constraint index for updating atom data */
84 int nind_r; /* number of indices */
85 int *ind_r; /* constraint index for updating atom data */
86 int ind_nalloc; /* allocation size of ind and ind_r */
87 tensor vir_r_m_dr; /* temporary variable for virial calculation */
88 real dhdlambda; /* temporary variable for lambda derivative */
89 real *simd_buf; /* Aligned work array for SIMD */
92 typedef struct gmx_lincsdata {
93 int ncg; /* the global number of constraints */
94 int ncg_flex; /* the global number of flexible constraints */
95 int ncg_triangle; /* the global number of constraints in triangles */
96 int nIter; /* the number of iterations */
97 int nOrder; /* the order of the matrix expansion */
98 int max_connect; /* the maximum number of constrains connected to a single atom */
100 int nc_real; /* the number of real constraints */
101 int nc; /* the number of constraints including padding for SIMD */
102 int nc_alloc; /* the number we allocated memory for */
103 int ncc; /* the number of constraint connections */
104 int ncc_alloc; /* the number we allocated memory for */
105 real matlam; /* the FE lambda value used for filling blc and blmf */
106 int *con_index; /* mapping from topology to LINCS constraints */
107 real *bllen0; /* the reference distance in topology A */
108 real *ddist; /* the reference distance in top B - the r.d. in top A */
109 int *bla; /* the atom pairs involved in the constraints */
110 real *blc; /* 1/sqrt(invmass1 + invmass2) */
111 real *blc1; /* as blc, but with all masses 1 */
112 int *blnr; /* index into blbnb and blmf */
113 int *blbnb; /* list of constraint connections */
114 int ntriangle; /* the local number of constraints in triangles */
115 int ncc_triangle; /* the number of constraint connections in triangles */
116 gmx_bool bCommIter; /* communicate before each LINCS interation */
117 real *blmf; /* matrix of mass factors for constraint connections */
118 real *blmf1; /* as blmf, but with all masses 1 */
119 real *bllen; /* the reference bond length */
120 int *nlocat; /* the local atom count per constraint, can be NULL */
122 int ntask; /* The number of tasks = #threads for LINCS */
123 lincs_task_t *task; /* LINCS thread division */
124 gmx_bitmask_t *atf; /* atom flags for thread parallelization */
125 int atf_nalloc; /* allocation size of atf */
126 gmx_bool bTaskDep; /* are the LINCS tasks interdependent? */
127 /* arrays for temporary storage in the LINCS algorithm */
134 real *mlambda; /* the Lagrange multipliers * -1 */
135 /* storage for the constraint RMS relative deviation output */
139 /* Define simd_width for memory allocation used for SIMD code */
141 static const int simd_width = GMX_SIMD_REAL_WIDTH;
143 static const int simd_width = 1;
145 /* We can't use small memory alignments on many systems, so use min 64 bytes*/
146 static const int align_bytes = std::max<int>(64, simd_width*sizeof(real));
148 real *lincs_rmsd_data(struct gmx_lincsdata *lincsd)
150 return lincsd->rmsd_data;
153 real lincs_rmsd(struct gmx_lincsdata *lincsd, gmx_bool bSD2)
155 if (lincsd->rmsd_data[0] > 0)
157 return sqrt(lincsd->rmsd_data[bSD2 ? 2 : 1]/lincsd->rmsd_data[0]);
165 /* Do a set of nrec LINCS matrix multiplications.
166 * This function will return with up to date thread-local
167 * constraint data, without an OpenMP barrier.
169 static void lincs_matrix_expand(const struct gmx_lincsdata *lincsd,
170 const lincs_task_t *li_task,
172 real *rhs1, real *rhs2, real *sol)
174 int b0, b1, nrec, rec;
175 const int *blnr = lincsd->blnr;
176 const int *blbnb = lincsd->blbnb;
180 nrec = lincsd->nOrder;
182 for (rec = 0; rec < nrec; rec++)
186 if (lincsd->bTaskDep)
190 for (b = b0; b < b1; b++)
196 for (n = blnr[b]; n < blnr[b+1]; n++)
198 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
201 sol[b] = sol[b] + mvb;
209 } /* nrec*(ncons+2*nrtot) flops */
211 if (lincsd->ntriangle > 0)
213 /* Perform an extra nrec recursions for only the constraints
214 * involved in rigid triangles.
215 * In this way their accuracy should come close to those of the other
216 * constraints, since traingles of constraints can produce eigenvalues
217 * around 0.7, while the effective eigenvalue for bond constraints
218 * is around 0.4 (and 0.7*0.7=0.5).
221 if (lincsd->bTaskDep)
223 /* We need a barrier here, since other threads might still be
224 * reading the contents of rhs1 and/o rhs2.
225 * We could avoid this barrier by introducing two extra rhs
226 * arrays for the triangle constraints only.
231 /* Constraints involved in a triangle are ensured to be in the same
232 * LINCS task. This means no barriers are required during the extra
233 * iterations for the triangle constraints.
235 const int *triangle = li_task->triangle;
236 const int *tri_bits = li_task->tri_bits;
238 for (rec = 0; rec < nrec; rec++)
242 for (tb = 0; tb < li_task->ntriangle; tb++)
244 int b, bits, nr0, nr1, n;
252 for (n = nr0; n < nr1; n++)
254 if (bits & (1 << (n - nr0)))
256 mvb = mvb + blcc[n]*rhs1[blbnb[n]];
260 sol[b] = sol[b] + mvb;
268 } /* nrec*(ntriangle + ncc_triangle*2) flops */
272 static void lincs_update_atoms_noind(int ncons, const int *bla,
274 const real *fac, rvec *r,
279 real mvb, im1, im2, tmp0, tmp1, tmp2;
283 for (b = 0; b < ncons; b++)
299 } /* 16 ncons flops */
303 for (b = 0; b < ncons; b++)
321 static void lincs_update_atoms_ind(int ncons, const int *ind, const int *bla,
323 const real *fac, rvec *r,
328 real mvb, im1, im2, tmp0, tmp1, tmp2;
332 for (bi = 0; bi < ncons; bi++)
349 } /* 16 ncons flops */
353 for (bi = 0; bi < ncons; bi++)
368 } /* 16 ncons flops */
372 static void lincs_update_atoms(struct gmx_lincsdata *li, int th,
374 const real *fac, rvec *r,
380 /* Single thread, we simply update for all constraints */
381 lincs_update_atoms_noind(li->nc_real,
382 li->bla, prefac, fac, r, invmass, x);
386 /* Update the atom vector components for our thread local
387 * constraints that only access our local atom range.
388 * This can be done without a barrier.
390 lincs_update_atoms_ind(li->task[th].nind, li->task[th].ind,
391 li->bla, prefac, fac, r, invmass, x);
393 if (li->task[li->ntask].nind > 0)
395 /* Update the constraints that operate on atoms
396 * in multiple thread atom blocks on the master thread.
401 lincs_update_atoms_ind(li->task[li->ntask].nind,
402 li->task[li->ntask].ind,
403 li->bla, prefac, fac, r, invmass, x);
410 /* Calculate distances between indexed atom coordinate pairs
411 * and store the result in 3 SIMD registers through an aligned buffer.
412 * Start from index is and load GMX_SIMD_REAL_WIDTH elements.
413 * Note that pair_index must contain valid indices up to is+GMX_SIMD_REAL_WIDTH.
415 static gmx_inline void gmx_simdcall
416 rvec_dist_to_simd(const rvec *x,
418 const int *pair_index,
426 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
428 /* Store the non PBC corrected distances packed and aligned */
429 for (m = 0; m < DIM; m++)
431 buf[m*GMX_SIMD_REAL_WIDTH + s] =
432 x[pair_index[2*(is+s)]][m] - x[pair_index[2*(is+s) + 1]][m];
435 *dx = gmx_simd_load_r(buf + 0*GMX_SIMD_REAL_WIDTH);
436 *dy = gmx_simd_load_r(buf + 1*GMX_SIMD_REAL_WIDTH);
437 *dz = gmx_simd_load_r(buf + 2*GMX_SIMD_REAL_WIDTH);
440 /* Store a SIMD vector in GMX_SIMD_REAL_WIDTH rvecs through an aligned buffer */
441 static gmx_inline void gmx_simdcall
442 copy_simd_vec_to_rvec(gmx_simd_real_t x,
451 gmx_simd_store_r(buf + 0*GMX_SIMD_REAL_WIDTH, x);
452 gmx_simd_store_r(buf + 1*GMX_SIMD_REAL_WIDTH, y);
453 gmx_simd_store_r(buf + 2*GMX_SIMD_REAL_WIDTH, z);
455 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
457 for (m = 0; m < DIM; m++)
459 v[is + s][m] = buf[m*GMX_SIMD_REAL_WIDTH + s];
464 /* Calculate the constraint distance vectors r to project on from x.
465 * Determine the right-hand side of the matrix equation using quantity f.
466 * This function only differs from calc_dr_x_xp_simd below in that
467 * no constraint length is subtracted and no PBC is used for f.
469 static void gmx_simdcall
470 calc_dr_x_f_simd(int b0,
473 const rvec * gmx_restrict x,
474 const rvec * gmx_restrict f,
475 const real * gmx_restrict blc,
476 const pbc_simd_t * pbc_simd,
477 real * gmx_restrict vbuf1,
478 real * gmx_restrict vbuf2,
479 rvec * gmx_restrict r,
480 real * gmx_restrict rhs,
481 real * gmx_restrict sol)
485 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
487 for (bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
489 gmx_simd_real_t rx_S, ry_S, rz_S, n2_S, il_S;
490 gmx_simd_real_t fx_S, fy_S, fz_S, ip_S, rhs_S;
492 rvec_dist_to_simd(x, bs, bla, vbuf1, &rx_S, &ry_S, &rz_S);
494 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
496 n2_S = gmx_simd_norm2_r(rx_S, ry_S, rz_S);
497 il_S = gmx_simd_invsqrt_r(n2_S);
499 rx_S = gmx_simd_mul_r(rx_S, il_S);
500 ry_S = gmx_simd_mul_r(ry_S, il_S);
501 rz_S = gmx_simd_mul_r(rz_S, il_S);
503 copy_simd_vec_to_rvec(rx_S, ry_S, rz_S, vbuf1, r, bs);
505 rvec_dist_to_simd(f, bs, bla, vbuf2, &fx_S, &fy_S, &fz_S);
507 ip_S = gmx_simd_iprod_r(rx_S, ry_S, rz_S,
510 rhs_S = gmx_simd_mul_r(gmx_simd_load_r(blc + bs), ip_S);
512 gmx_simd_store_r(rhs + bs, rhs_S);
513 gmx_simd_store_r(sol + bs, rhs_S);
516 #endif /* LINCS_SIMD */
518 /* LINCS projection, works on derivatives of the coordinates */
519 static void do_lincsp(rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
520 struct gmx_lincsdata *lincsd, int th,
522 int econq, gmx_bool bCalcDHDL,
523 gmx_bool bCalcVir, tensor rmdf)
526 int *bla, *blnr, *blbnb;
528 real *blc, *blmf, *blcc, *rhs1, *rhs2, *sol;
530 b0 = lincsd->task[th].b0;
531 b1 = lincsd->task[th].b1;
536 blbnb = lincsd->blbnb;
537 if (econq != econqForce)
539 /* Use mass-weighted parameters */
545 /* Use non mass-weighted parameters */
547 blmf = lincsd->blmf1;
549 blcc = lincsd->tmpncc;
556 /* This SIMD code does the same as the plain-C code after the #else.
557 * The only difference is that we always call pbc code, as with SIMD
558 * the overhead of pbc computation (when not needed) is small.
562 /* Convert the pbc struct for SIMD */
563 set_pbc_simd(pbc, &pbc_simd);
565 /* Compute normalized x i-j vectors, store in r.
566 * Compute the inner product of r and xp i-j and store in rhs1.
568 calc_dr_x_f_simd(b0, b1, bla, x, f, blc,
570 lincsd->task[th].simd_buf,
571 lincsd->task[th].simd_buf + GMX_SIMD_REAL_WIDTH*DIM,
574 #else /* LINCS_SIMD */
576 /* Compute normalized i-j vectors */
579 for (b = b0; b < b1; b++)
583 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
589 for (b = b0; b < b1; b++)
593 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
595 } /* 16 ncons flops */
598 for (b = b0; b < b1; b++)
605 mvb = blc[b]*(r[b][0]*(f[i][0] - f[j][0]) +
606 r[b][1]*(f[i][1] - f[j][1]) +
607 r[b][2]*(f[i][2] - f[j][2]));
613 #endif /* LINCS_SIMD */
615 if (lincsd->bTaskDep)
617 /* We need a barrier, since the matrix construction below
618 * can access entries in r of other threads.
623 /* Construct the (sparse) LINCS matrix */
624 for (b = b0; b < b1; b++)
628 for (n = blnr[b]; n < blnr[b+1]; n++)
630 blcc[n] = blmf[n]*iprod(r[b], r[blbnb[n]]);
633 /* Together: 23*ncons + 6*nrtot flops */
635 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
636 /* nrec*(ncons+2*nrtot) flops */
638 if (econq == econqDeriv_FlexCon)
640 /* We only want to constraint the flexible constraints,
641 * so we mask out the normal ones by setting sol to 0.
643 for (b = b0; b < b1; b++)
645 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
652 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
653 for (b = b0; b < b1; b++)
658 /* When constraining forces, we should not use mass weighting,
659 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
661 lincs_update_atoms(lincsd, th, 1.0, sol, r,
662 (econq != econqForce) ? invmass : NULL, fp);
669 for (b = b0; b < b1; b++)
671 dhdlambda -= sol[b]*lincsd->ddist[b];
674 lincsd->task[th].dhdlambda = dhdlambda;
679 /* Constraint virial,
680 * determines sum r_bond x delta f,
681 * where delta f is the constraint correction
682 * of the quantity that is being constrained.
684 for (b = b0; b < b1; b++)
689 mvb = lincsd->bllen[b]*sol[b];
690 for (i = 0; i < DIM; i++)
693 for (j = 0; j < DIM; j++)
695 rmdf[i][j] += tmp1*r[b][j];
698 } /* 23 ncons flops */
703 /* Calculate the constraint distance vectors r to project on from x.
704 * Determine the right-hand side of the matrix equation using coordinates xp.
706 static void gmx_simdcall
707 calc_dr_x_xp_simd(int b0,
710 const rvec * gmx_restrict x,
711 const rvec * gmx_restrict xp,
712 const real * gmx_restrict bllen,
713 const real * gmx_restrict blc,
714 const pbc_simd_t * pbc_simd,
715 real * gmx_restrict vbuf1,
716 real * gmx_restrict vbuf2,
717 rvec * gmx_restrict r,
718 real * gmx_restrict rhs,
719 real * gmx_restrict sol)
723 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
725 for (bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
727 gmx_simd_real_t rx_S, ry_S, rz_S, n2_S, il_S;
728 gmx_simd_real_t rxp_S, ryp_S, rzp_S, ip_S, rhs_S;
730 rvec_dist_to_simd(x, bs, bla, vbuf1, &rx_S, &ry_S, &rz_S);
732 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
734 n2_S = gmx_simd_norm2_r(rx_S, ry_S, rz_S);
735 il_S = gmx_simd_invsqrt_r(n2_S);
737 rx_S = gmx_simd_mul_r(rx_S, il_S);
738 ry_S = gmx_simd_mul_r(ry_S, il_S);
739 rz_S = gmx_simd_mul_r(rz_S, il_S);
741 copy_simd_vec_to_rvec(rx_S, ry_S, rz_S, vbuf1, r, bs);
743 rvec_dist_to_simd(xp, bs, bla, vbuf2, &rxp_S, &ryp_S, &rzp_S);
745 pbc_correct_dx_simd(&rxp_S, &ryp_S, &rzp_S, pbc_simd);
747 ip_S = gmx_simd_iprod_r(rx_S, ry_S, rz_S,
748 rxp_S, ryp_S, rzp_S);
750 rhs_S = gmx_simd_mul_r(gmx_simd_load_r(blc + bs),
751 gmx_simd_sub_r(ip_S, gmx_simd_load_r(bllen + bs)));
753 gmx_simd_store_r(rhs + bs, rhs_S);
754 gmx_simd_store_r(sol + bs, rhs_S);
757 #endif /* LINCS_SIMD */
759 /* Determine the distances and right-hand side for the next iteration */
760 static void calc_dist_iter(int b0,
763 const rvec * gmx_restrict xp,
764 const real * gmx_restrict bllen,
765 const real * gmx_restrict blc,
768 real * gmx_restrict rhs,
769 real * gmx_restrict sol,
774 for (b = b0; b < b1; b++)
776 real len, len2, dlen2, mvb;
782 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
786 rvec_sub(xp[bla[2*b]], xp[bla[2*b+1]], dx);
789 dlen2 = 2*len2 - norm2(dx);
790 if (dlen2 < wfac*len2)
792 /* not race free - see detailed comment in caller */
797 mvb = blc[b]*(len - dlen2*gmx_invsqrt(dlen2));
805 } /* 20*ncons flops */
809 /* As the function above, but using SIMD intrinsics */
810 static void gmx_simdcall
811 calc_dist_iter_simd(int b0,
814 const rvec * gmx_restrict x,
815 const real * gmx_restrict bllen,
816 const real * gmx_restrict blc,
817 const pbc_simd_t * pbc_simd,
819 real * gmx_restrict vbuf,
820 real * gmx_restrict rhs,
821 real * gmx_restrict sol,
824 gmx_simd_real_t min_S = gmx_simd_set1_r(GMX_REAL_MIN);
825 gmx_simd_real_t two_S = gmx_simd_set1_r(2.0);
826 gmx_simd_real_t wfac_S = gmx_simd_set1_r(wfac);
827 gmx_simd_bool_t warn_S;
831 assert(b0 % GMX_SIMD_REAL_WIDTH == 0);
833 /* Initialize all to FALSE */
834 warn_S = gmx_simd_cmplt_r(two_S, gmx_simd_setzero_r());
836 for (bs = b0; bs < b1; bs += GMX_SIMD_REAL_WIDTH)
838 gmx_simd_real_t rx_S, ry_S, rz_S, n2_S;
839 gmx_simd_real_t len_S, len2_S, dlen2_S, lc_S, blc_S;
841 rvec_dist_to_simd(x, bs, bla, vbuf, &rx_S, &ry_S, &rz_S);
843 pbc_correct_dx_simd(&rx_S, &ry_S, &rz_S, pbc_simd);
845 n2_S = gmx_simd_norm2_r(rx_S, ry_S, rz_S);
847 len_S = gmx_simd_load_r(bllen + bs);
848 len2_S = gmx_simd_mul_r(len_S, len_S);
850 dlen2_S = gmx_simd_fmsub_r(two_S, len2_S, n2_S);
852 warn_S = gmx_simd_or_b(warn_S,
853 gmx_simd_cmplt_r(dlen2_S,
854 gmx_simd_mul_r(wfac_S, len2_S)));
856 /* Avoid 1/0 by taking the max with REAL_MIN.
857 * Note: when dlen2 is close to zero (90 degree constraint rotation),
858 * the accuracy of the algorithm is no longer relevant.
860 dlen2_S = gmx_simd_max_r(dlen2_S, min_S);
862 lc_S = gmx_simd_fnmadd_r(dlen2_S, gmx_simd_invsqrt_r(dlen2_S), len_S);
864 blc_S = gmx_simd_load_r(blc + bs);
866 lc_S = gmx_simd_mul_r(blc_S, lc_S);
868 gmx_simd_store_r(rhs + bs, lc_S);
869 gmx_simd_store_r(sol + bs, lc_S);
872 if (gmx_simd_anytrue_b(warn_S))
877 #endif /* LINCS_SIMD */
879 static void do_lincs(rvec *x, rvec *xp, matrix box, t_pbc *pbc,
880 struct gmx_lincsdata *lincsd, int th,
884 real wangle, gmx_bool *bWarn,
885 real invdt, rvec * gmx_restrict v,
886 gmx_bool bCalcVir, tensor vir_r_m_dr)
888 int b0, b1, b, i, j, n, iter;
889 int *bla, *blnr, *blbnb;
891 real *blc, *blmf, *bllen, *blcc, *rhs1, *rhs2, *sol, *blc_sol, *mlambda;
894 b0 = lincsd->task[th].b0;
895 b1 = lincsd->task[th].b1;
900 blbnb = lincsd->blbnb;
903 bllen = lincsd->bllen;
904 blcc = lincsd->tmpncc;
908 blc_sol = lincsd->tmp4;
909 mlambda = lincsd->mlambda;
910 nlocat = lincsd->nlocat;
914 /* This SIMD code does the same as the plain-C code after the #else.
915 * The only difference is that we always call pbc code, as with SIMD
916 * the overhead of pbc computation (when not needed) is small.
920 /* Convert the pbc struct for SIMD */
921 set_pbc_simd(pbc, &pbc_simd);
923 /* Compute normalized x i-j vectors, store in r.
924 * Compute the inner product of r and xp i-j and store in rhs1.
926 calc_dr_x_xp_simd(b0, b1, bla, x, xp, bllen, blc,
928 lincsd->task[th].simd_buf,
929 lincsd->task[th].simd_buf + GMX_SIMD_REAL_WIDTH*DIM,
932 #else /* LINCS_SIMD */
936 /* Compute normalized i-j vectors */
937 for (b = b0; b < b1; b++)
942 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
945 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
946 mvb = blc[b]*(iprod(r[b], dx) - bllen[b]);
953 /* Compute normalized i-j vectors */
954 for (b = b0; b < b1; b++)
956 real tmp0, tmp1, tmp2, rlen, mvb;
960 tmp0 = x[i][0] - x[j][0];
961 tmp1 = x[i][1] - x[j][1];
962 tmp2 = x[i][2] - x[j][2];
963 rlen = gmx_invsqrt(tmp0*tmp0 + tmp1*tmp1 + tmp2*tmp2);
971 mvb = blc[b]*(r[b][0]*(xp[i][0] - xp[j][0]) +
972 r[b][1]*(xp[i][1] - xp[j][1]) +
973 r[b][2]*(xp[i][2] - xp[j][2]) - bllen[b]);
978 /* Together: 26*ncons + 6*nrtot flops */
981 #endif /* LINCS_SIMD */
983 if (lincsd->bTaskDep)
985 /* We need a barrier, since the matrix construction below
986 * can access entries in r of other threads.
991 /* Construct the (sparse) LINCS matrix */
992 for (b = b0; b < b1; b++)
994 for (n = blnr[b]; n < blnr[b+1]; n++)
996 blcc[n] = blmf[n]*iprod(r[b], r[blbnb[n]]);
999 /* Together: 26*ncons + 6*nrtot flops */
1001 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1002 /* nrec*(ncons+2*nrtot) flops */
1005 for (b = b0; b < b1; b++)
1007 mlambda[b] = blc[b]*sol[b];
1010 for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1012 gmx_simd_store_r(mlambda + b,
1013 gmx_simd_mul_r(gmx_simd_load_r(blc + b),
1014 gmx_simd_load_r(sol + b)));
1018 /* Update the coordinates */
1019 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
1022 ******** Correction for centripetal effects ********
1027 wfac = cos(DEG2RAD*wangle);
1030 for (iter = 0; iter < lincsd->nIter; iter++)
1032 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints))
1037 /* Communicate the corrected non-local coordinates */
1038 if (DOMAINDECOMP(cr))
1040 dd_move_x_constraints(cr->dd, box, xp, NULL, FALSE);
1045 else if (lincsd->bTaskDep)
1051 calc_dist_iter_simd(b0, b1, bla, xp, bllen, blc, &pbc_simd, wfac,
1052 lincsd->task[th].simd_buf, rhs1, sol, bWarn);
1054 calc_dist_iter(b0, b1, bla, xp, bllen, blc, pbc, wfac,
1056 /* 20*ncons flops */
1059 lincs_matrix_expand(lincsd, &lincsd->task[th], blcc, rhs1, rhs2, sol);
1060 /* nrec*(ncons+2*nrtot) flops */
1063 for (b = b0; b < b1; b++)
1067 mvb = blc[b]*sol[b];
1072 for (b = b0; b < b1; b += GMX_SIMD_REAL_WIDTH)
1074 gmx_simd_real_t mvb;
1076 mvb = gmx_simd_mul_r(gmx_simd_load_r(blc + b),
1077 gmx_simd_load_r(sol + b));
1078 gmx_simd_store_r(blc_sol + b, mvb);
1079 gmx_simd_store_r(mlambda + b,
1080 gmx_simd_add_r(gmx_simd_load_r(mlambda + b), mvb));
1084 /* Update the coordinates */
1085 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
1087 /* nit*ncons*(37+9*nrec) flops */
1091 /* Update the velocities */
1092 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
1093 /* 16 ncons flops */
1096 if (nlocat != NULL && (bCalcDHDL || bCalcVir))
1098 if (lincsd->bTaskDep)
1100 /* In lincs_update_atoms threads might cross-read mlambda */
1104 /* Only account for local atoms */
1105 for (b = b0; b < b1; b++)
1107 mlambda[b] *= 0.5*nlocat[b];
1116 for (b = b0; b < b1; b++)
1118 /* Note that this this is dhdl*dt^2, the dt^2 factor is corrected
1119 * later after the contributions are reduced over the threads.
1121 dhdl -= lincsd->mlambda[b]*lincsd->ddist[b];
1123 lincsd->task[th].dhdlambda = dhdl;
1128 /* Constraint virial */
1129 for (b = b0; b < b1; b++)
1133 tmp0 = -bllen[b]*mlambda[b];
1134 for (i = 0; i < DIM; i++)
1136 tmp1 = tmp0*r[b][i];
1137 for (j = 0; j < DIM; j++)
1139 vir_r_m_dr[i][j] -= tmp1*r[b][j];
1142 } /* 22 ncons flops */
1146 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
1147 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
1149 * (26+nrec)*ncons + (6+2*nrec)*nrtot
1150 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
1152 * (63+nrec)*ncons + (6+4*nrec)*nrtot
1156 /* Sets the elements in the LINCS matrix for task li_task */
1157 static void set_lincs_matrix_task(struct gmx_lincsdata *li,
1158 lincs_task_t *li_task,
1159 const real *invmass,
1164 /* Construct the coupling coefficient matrix blmf */
1165 li_task->ntriangle = 0;
1167 for (i = li_task->b0; i < li_task->b1; i++)
1172 a2 = li->bla[2*i+1];
1173 for (n = li->blnr[i]; (n < li->blnr[i+1]); n++)
1175 int k, sign, center, end;
1179 /* If we are using multiple, independent tasks for LINCS,
1180 * the calls to check_assign_connected should have
1181 * put all connected constraints in our task.
1183 assert(li->bTaskDep || (k >= li_task->b0 && k < li_task->b1));
1185 if (a1 == li->bla[2*k] || a2 == li->bla[2*k+1])
1193 if (a1 == li->bla[2*k] || a1 == li->bla[2*k+1])
1203 li->blmf[n] = sign*invmass[center]*li->blc[i]*li->blc[k];
1204 li->blmf1[n] = sign*0.5;
1205 if (li->ncg_triangle > 0)
1209 /* Look for constraint triangles */
1210 for (nk = li->blnr[k]; (nk < li->blnr[k+1]); nk++)
1213 if (kk != i && kk != k &&
1214 (li->bla[2*kk] == end || li->bla[2*kk+1] == end))
1216 /* If we are using multiple tasks for LINCS,
1217 * the calls to check_assign_triangle should have
1218 * put all constraints in the triangle in our task.
1220 assert(k >= li_task->b0 && k < li_task->b1);
1221 assert(kk >= li_task->b0 && kk < li_task->b1);
1223 if (li_task->ntriangle == 0 ||
1224 li_task->triangle[li_task->ntriangle - 1] < i)
1226 /* Add this constraint to the triangle list */
1227 li_task->triangle[li_task->ntriangle] = i;
1228 li_task->tri_bits[li_task->ntriangle] = 0;
1229 li_task->ntriangle++;
1230 if (li->blnr[i+1] - li->blnr[i] > static_cast<int>(sizeof(li_task->tri_bits[0])*8 - 1))
1232 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %d allowed for constraints participating in triangles",
1233 li->blnr[i+1] - li->blnr[i],
1234 sizeof(li_task->tri_bits[0])*8-1);
1237 li_task->tri_bits[li_task->ntriangle-1] |= (1 << (n - li->blnr[i]));
1246 /* Sets the elements in the LINCS matrix */
1247 void set_lincs_matrix(struct gmx_lincsdata *li, real *invmass, real lambda)
1250 const real invsqrt2 = 0.7071067811865475244;
1252 for (i = 0; (i < li->nc); i++)
1257 a2 = li->bla[2*i+1];
1258 li->blc[i] = gmx_invsqrt(invmass[a1] + invmass[a2]);
1259 li->blc1[i] = invsqrt2;
1262 /* Construct the coupling coefficient matrix blmf */
1263 int th, ntriangle = 0, ncc_triangle = 0;
1264 #pragma omp parallel for reduction(+: ntriangle, ncc_triangle) num_threads(li->ntask) schedule(static)
1265 for (th = 0; th < li->ntask; th++)
1267 set_lincs_matrix_task(li, &li->task[th], invmass, &ncc_triangle);
1268 ntriangle = li->task[th].ntriangle;
1270 li->ntriangle = ntriangle;
1271 li->ncc_triangle = ncc_triangle;
1275 fprintf(debug, "Of the %d constraints %d participate in triangles\n",
1276 li->nc, li->ntriangle);
1277 fprintf(debug, "There are %d couplings of which %d in triangles\n",
1278 li->ncc, li->ncc_triangle);
1282 * so we know with which lambda value the masses have been set.
1284 li->matlam = lambda;
1287 static int count_triangle_constraints(t_ilist *ilist, t_blocka *at2con)
1289 int ncon1, ncon_tot;
1290 int c0, a00, a01, n1, c1, a10, a11, ac1, n2, c2, a20, a21;
1293 t_iatom *ia1, *ia2, *iap;
1295 ncon1 = ilist[F_CONSTR].nr/3;
1296 ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
1298 ia1 = ilist[F_CONSTR].iatoms;
1299 ia2 = ilist[F_CONSTRNC].iatoms;
1302 for (c0 = 0; c0 < ncon_tot; c0++)
1305 iap = constr_iatomptr(ncon1, ia1, ia2, c0);
1308 for (n1 = at2con->index[a01]; n1 < at2con->index[a01+1]; n1++)
1313 iap = constr_iatomptr(ncon1, ia1, ia2, c1);
1324 for (n2 = at2con->index[ac1]; n2 < at2con->index[ac1+1]; n2++)
1327 if (c2 != c0 && c2 != c1)
1329 iap = constr_iatomptr(ncon1, ia1, ia2, c2);
1332 if (a20 == a00 || a21 == a00)
1346 return ncon_triangle;
1349 static gmx_bool more_than_two_sequential_constraints(const t_ilist *ilist,
1350 const t_blocka *at2con)
1352 t_iatom *ia1, *ia2, *iap;
1353 int ncon1, ncon_tot, c;
1355 gmx_bool bMoreThanTwoSequentialConstraints;
1357 ncon1 = ilist[F_CONSTR].nr/3;
1358 ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
1360 ia1 = ilist[F_CONSTR].iatoms;
1361 ia2 = ilist[F_CONSTRNC].iatoms;
1363 bMoreThanTwoSequentialConstraints = FALSE;
1364 for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
1366 iap = constr_iatomptr(ncon1, ia1, ia2, c);
1369 /* Check if this constraint has constraints connected at both atoms */
1370 if (at2con->index[a1+1] - at2con->index[a1] > 1 &&
1371 at2con->index[a2+1] - at2con->index[a2] > 1)
1373 bMoreThanTwoSequentialConstraints = TRUE;
1377 return bMoreThanTwoSequentialConstraints;
1380 static int int_comp(const void *a, const void *b)
1382 return (*(int *)a) - (*(int *)b);
1385 gmx_lincsdata_t init_lincs(FILE *fplog, gmx_mtop_t *mtop,
1386 int nflexcon_global, t_blocka *at2con,
1387 gmx_bool bPLINCS, int nIter, int nProjOrder)
1389 struct gmx_lincsdata *li;
1391 gmx_moltype_t *molt;
1392 gmx_bool bMoreThanTwoSeq;
1396 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
1397 bPLINCS ? " Parallel" : "");
1403 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1404 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1405 li->ncg_flex = nflexcon_global;
1408 li->nOrder = nProjOrder;
1410 li->max_connect = 0;
1411 for (mt = 0; mt < mtop->nmoltype; mt++)
1415 molt = &mtop->moltype[mt];
1416 for (a = 0; a < molt->atoms.nr; a++)
1418 li->max_connect = std::max(li->max_connect,
1419 at2con[mt].index[a + 1] - at2con[mt].index[a]);
1423 li->ncg_triangle = 0;
1424 bMoreThanTwoSeq = FALSE;
1425 for (mb = 0; mb < mtop->nmolblock; mb++)
1427 molt = &mtop->moltype[mtop->molblock[mb].type];
1430 mtop->molblock[mb].nmol*
1431 count_triangle_constraints(molt->ilist,
1432 &at2con[mtop->molblock[mb].type]);
1434 if (!bMoreThanTwoSeq &&
1435 more_than_two_sequential_constraints(molt->ilist, &at2con[mtop->molblock[mb].type]))
1437 bMoreThanTwoSeq = TRUE;
1441 /* Check if we need to communicate not only before LINCS,
1442 * but also before each iteration.
1443 * The check for only two sequential constraints is only
1444 * useful for the common case of H-bond only constraints.
1445 * With more effort we could also make it useful for small
1446 * molecules with nr. sequential constraints <= nOrder-1.
1448 li->bCommIter = (bPLINCS && (li->nOrder < 1 || bMoreThanTwoSeq));
1450 if (debug && bPLINCS)
1452 fprintf(debug, "PLINCS communication before each iteration: %d\n",
1456 /* LINCS can run on any number of threads.
1457 * Currently the number is fixed for the whole simulation,
1458 * but it could be set in set_lincs().
1459 * The current constraint to task assignment code can create independent
1460 * tasks only when not more than two constraints are connected sequentially.
1462 li->ntask = gmx_omp_nthreads_get(emntLINCS);
1463 li->bTaskDep = (li->ntask > 1 && bMoreThanTwoSeq);
1466 fprintf(debug, "LINCS: using %d threads, tasks are %sdependent\n",
1467 li->ntask, li->bTaskDep ? "" : "in");
1475 /* Allocate an extra elements for "task-overlap" constraints */
1476 snew(li->task, li->ntask + 1);
1479 #pragma omp parallel for num_threads(li->ntask)
1480 for (th = 0; th < li->ntask; th++)
1482 /* Per thread SIMD load buffer for loading 2 simd_width rvecs */
1483 snew_aligned(li->task[th].simd_buf, 2*simd_width*DIM,
1487 if (bPLINCS || li->ncg_triangle > 0)
1489 please_cite(fplog, "Hess2008a");
1493 please_cite(fplog, "Hess97a");
1498 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1501 fprintf(fplog, "There are inter charge-group constraints,\n"
1502 "will communicate selected coordinates each lincs iteration\n");
1504 if (li->ncg_triangle > 0)
1507 "%d constraints are involved in constraint triangles,\n"
1508 "will apply an additional matrix expansion of order %d for couplings\n"
1509 "between constraints inside triangles\n",
1510 li->ncg_triangle, li->nOrder);
1517 /* Sets up the work division over the threads */
1518 static void lincs_thread_setup(struct gmx_lincsdata *li, int natoms)
1525 if (natoms > li->atf_nalloc)
1527 li->atf_nalloc = over_alloc_large(natoms);
1528 srenew(li->atf, li->atf_nalloc);
1532 /* Clear the atom flags */
1533 for (a = 0; a < natoms; a++)
1535 bitmask_clear(&atf[a]);
1538 if (li->ntask > BITMASK_SIZE)
1540 gmx_fatal(FARGS, "More than %d threads is not supported for LINCS.", BITMASK_SIZE);
1543 for (th = 0; th < li->ntask; th++)
1545 lincs_task_t *li_task;
1548 li_task = &li->task[th];
1550 /* For each atom set a flag for constraints from each */
1551 for (b = li_task->b0; b < li_task->b1; b++)
1553 bitmask_set_bit(&atf[li->bla[b*2 ]], th);
1554 bitmask_set_bit(&atf[li->bla[b*2 + 1]], th);
1558 #pragma omp parallel for num_threads(li->ntask) schedule(static)
1559 for (th = 0; th < li->ntask; th++)
1561 lincs_task_t *li_task;
1565 li_task = &li->task[th];
1567 if (li_task->b1 - li_task->b0 > li_task->ind_nalloc)
1569 li_task->ind_nalloc = over_alloc_large(li_task->b1-li_task->b0);
1570 srenew(li_task->ind, li_task->ind_nalloc);
1571 srenew(li_task->ind_r, li_task->ind_nalloc);
1574 bitmask_init_low_bits(&mask, th);
1577 li_task->nind_r = 0;
1578 for (b = li_task->b0; b < li_task->b1; b++)
1580 /* We let the constraint with the lowest thread index
1581 * operate on atoms with constraints from multiple threads.
1583 if (bitmask_is_disjoint(atf[li->bla[b*2]], mask) &&
1584 bitmask_is_disjoint(atf[li->bla[b*2+1]], mask))
1586 /* Add the constraint to the local atom update index */
1587 li_task->ind[li_task->nind++] = b;
1591 /* Add the constraint to the rest block */
1592 li_task->ind_r[li_task->nind_r++] = b;
1597 /* We need to copy all constraints which have not be assigned
1598 * to a thread to a separate list which will be handled by one thread.
1600 li_m = &li->task[li->ntask];
1603 for (th = 0; th < li->ntask; th++)
1605 lincs_task_t *li_task;
1608 li_task = &li->task[th];
1610 if (li_m->nind + li_task->nind_r > li_m->ind_nalloc)
1612 li_m->ind_nalloc = over_alloc_large(li_m->nind+li_task->nind_r);
1613 srenew(li_m->ind, li_m->ind_nalloc);
1616 for (b = 0; b < li_task->nind_r; b++)
1618 li_m->ind[li_m->nind++] = li_task->ind_r[b];
1623 fprintf(debug, "LINCS thread %d: %d constraints\n",
1630 fprintf(debug, "LINCS thread r: %d constraints\n",
1635 /* There is no realloc with alignment, so here we make one for reals.
1636 * Note that this function does not preserve the contents of the memory.
1638 static void resize_real_aligned(real **ptr, int nelem)
1640 sfree_aligned(*ptr);
1641 snew_aligned(*ptr, nelem, align_bytes);
1644 static void assign_constraint(struct gmx_lincsdata *li,
1645 int constraint_index,
1647 real lenA, real lenB,
1648 const t_blocka *at2con)
1654 /* Make an mapping of local topology constraint index to LINCS index */
1655 li->con_index[constraint_index] = con;
1657 li->bllen0[con] = lenA;
1658 li->ddist[con] = lenB - lenA;
1659 /* Set the length to the topology A length */
1660 li->bllen[con] = lenA;
1661 li->bla[2*con] = a1;
1662 li->bla[2*con+1] = a2;
1664 /* Make space in the constraint connection matrix for constraints
1665 * connected to both end of the current constraint.
1668 at2con->index[a1 + 1] - at2con->index[a1] - 1 +
1669 at2con->index[a2 + 1] - at2con->index[a2] - 1;
1671 li->blnr[con + 1] = li->ncc;
1673 /* Increase the constraint count */
1677 /* Check if constraint with topology index constraint_index is connected
1678 * to other constraints, and if so add those connected constraints to our task.
1680 static void check_assign_connected(struct gmx_lincsdata *li,
1681 const t_iatom *iatom,
1685 const t_blocka *at2con)
1687 /* Currently this function only supports constraint groups
1688 * in which all constraints share at least one atom
1689 * (e.g. H-bond constraints).
1690 * Check both ends of the current constraint for
1691 * connected constraints. We need to assign those
1696 for (end = 0; end < 2; end++)
1700 a = (end == 0 ? a1 : a2);
1702 for (k = at2con->index[a]; k < at2con->index[a + 1]; k++)
1707 /* Check if constraint cc has not yet been assigned */
1708 if (li->con_index[cc] == -1)
1714 lenA = idef->iparams[type].constr.dA;
1715 lenB = idef->iparams[type].constr.dB;
1717 if (bDynamics || lenA != 0 || lenB != 0)
1719 assign_constraint(li, cc, iatom[3*cc + 1], iatom[3*cc + 2], lenA, lenB, at2con);
1726 /* Check if constraint with topology index constraint_index is involved
1727 * in a constraint triangle, and if so add the other two constraints
1728 * in the triangle to our task.
1730 static void check_assign_triangle(struct gmx_lincsdata *li,
1731 const t_iatom *iatom,
1734 int constraint_index,
1736 const t_blocka *at2con)
1738 int nca, cc[32], ca[32], k;
1739 int c_triangle[2] = { -1, -1 };
1742 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1747 if (c != constraint_index)
1751 aa1 = iatom[c*3 + 1];
1752 aa2 = iatom[c*3 + 2];
1768 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1773 if (c != constraint_index)
1777 aa1 = iatom[c*3 + 1];
1778 aa2 = iatom[c*3 + 2];
1781 for (i = 0; i < nca; i++)
1785 c_triangle[0] = cc[i];
1792 for (i = 0; i < nca; i++)
1796 c_triangle[0] = cc[i];
1804 if (c_triangle[0] >= 0)
1808 for (end = 0; end < 2; end++)
1810 /* Check if constraint c_triangle[end] has not yet been assigned */
1811 if (li->con_index[c_triangle[end]] == -1)
1816 i = c_triangle[end]*3;
1818 lenA = idef->iparams[type].constr.dA;
1819 lenB = idef->iparams[type].constr.dB;
1821 if (bDynamics || lenA != 0 || lenB != 0)
1823 assign_constraint(li, c_triangle[end], iatom[i + 1], iatom[i + 2], lenA, lenB, at2con);
1830 static void set_matrix_indices(struct gmx_lincsdata *li,
1831 const lincs_task_t *li_task,
1832 const t_blocka *at2con,
1833 gmx_bool bSortMatrix)
1837 for (b = li_task->b0; b < li_task->b1; b++)
1842 a2 = li->bla[b*2 + 1];
1845 for (k = at2con->index[a1]; k < at2con->index[a1 + 1]; k++)
1849 concon = li->con_index[at2con->a[k]];
1852 li->blbnb[i++] = concon;
1855 for (k = at2con->index[a2]; k < at2con->index[a2 + 1]; k++)
1859 concon = li->con_index[at2con->a[k]];
1862 li->blbnb[i++] = concon;
1868 /* Order the blbnb matrix to optimize memory access */
1869 qsort(&(li->blbnb[li->blnr[b]]), li->blnr[b + 1] - li->blnr[b],
1870 sizeof(li->blbnb[0]), int_comp);
1875 void set_lincs(t_idef *idef, t_mdatoms *md,
1876 gmx_bool bDynamics, t_commrec *cr,
1877 struct gmx_lincsdata *li)
1879 int natoms, nflexcon;
1882 int i, ncc_alloc_old, ncon_tot;
1887 /* Zero the thread index ranges.
1888 * Otherwise without local constraints we could return with old ranges.
1890 for (i = 0; i < li->ntask; i++)
1894 li->task[i].nind = 0;
1898 li->task[li->ntask].nind = 0;
1901 /* This is the local topology, so there are only F_CONSTR constraints */
1902 if (idef->il[F_CONSTR].nr == 0)
1904 /* There are no constraints,
1905 * we do not need to fill any data structures.
1912 fprintf(debug, "Building the LINCS connectivity\n");
1915 if (DOMAINDECOMP(cr))
1917 if (cr->dd->constraints)
1921 dd_get_constraint_range(cr->dd, &start, &natoms);
1925 natoms = cr->dd->nat_home;
1930 natoms = md->homenr;
1932 at2con = make_at2con(0, natoms, idef->il, idef->iparams, bDynamics,
1935 ncon_tot = idef->il[F_CONSTR].nr/3;
1937 /* Ensure we have enough padding for aligned loads for each thread */
1938 if (ncon_tot + li->ntask*simd_width > li->nc_alloc || li->nc_alloc == 0)
1940 li->nc_alloc = over_alloc_dd(ncon_tot + li->ntask*simd_width);
1941 srenew(li->con_index, li->nc_alloc);
1942 resize_real_aligned(&li->bllen0, li->nc_alloc);
1943 resize_real_aligned(&li->ddist, li->nc_alloc);
1944 srenew(li->bla, 2*li->nc_alloc);
1945 resize_real_aligned(&li->blc, li->nc_alloc);
1946 resize_real_aligned(&li->blc1, li->nc_alloc);
1947 srenew(li->blnr, li->nc_alloc + 1);
1948 resize_real_aligned(&li->bllen, li->nc_alloc);
1949 srenew(li->tmpv, li->nc_alloc);
1950 if (DOMAINDECOMP(cr))
1952 srenew(li->nlocat, li->nc_alloc);
1954 resize_real_aligned(&li->tmp1, li->nc_alloc);
1955 resize_real_aligned(&li->tmp2, li->nc_alloc);
1956 resize_real_aligned(&li->tmp3, li->nc_alloc);
1957 resize_real_aligned(&li->tmp4, li->nc_alloc);
1958 resize_real_aligned(&li->mlambda, li->nc_alloc);
1961 iatom = idef->il[F_CONSTR].iatoms;
1963 ncc_alloc_old = li->ncc_alloc;
1964 li->blnr[0] = li->ncc;
1966 /* Assign the constraints for li->ntask LINCS tasks.
1967 * We target a uniform distribution of constraints over the tasks.
1968 * Note that when flexible constraints are present, but are removed here
1969 * (e.g. because we are doing EM) we get imbalance, but since that doesn't
1970 * happen during normal MD, that's ok.
1972 int ncon_assign, ncon_target, con, th;
1974 /* Determine the number of constraints we need to assign here */
1975 ncon_assign = ncon_tot;
1978 /* With energy minimization, flexible constraints are ignored
1979 * (and thus minimized, as they should be).
1981 ncon_assign -= nflexcon;
1984 /* Set the target constraint count per task to exactly uniform,
1985 * this might be overridden below.
1987 ncon_target = (ncon_assign + li->ntask - 1)/li->ntask;
1989 /* Mark all constraints as unassigned by setting their index to -1 */
1990 for (con = 0; con < ncon_tot; con++)
1992 li->con_index[con] = -1;
1996 for (th = 0; th < li->ntask; th++)
1998 lincs_task_t *li_task;
2000 li_task = &li->task[th];
2003 /* With indepedent tasks we likely have H-bond constraints or constraint
2004 * pairs. The connected constraints will be pulled into the task, so the
2005 * constraints per task will often exceed ncon_target.
2006 * Triangle constraints can also increase the count, but there are
2007 * relatively few of those, so we usually expect to get ncon_target.
2011 /* We round ncon_target to a multiple of GMX_SIMD_WIDTH,
2012 * since otherwise a lot of operations can be wasted.
2013 * There are several ways to round here, we choose the one
2014 * that alternates block sizes, which helps with Intel HT.
2016 ncon_target = ((ncon_assign*(th + 1))/li->ntask - li->nc_real + GMX_SIMD_REAL_WIDTH - 1) & ~(GMX_SIMD_REAL_WIDTH - 1);
2020 /* Continue filling the arrays where we left off with the previous task,
2021 * including padding for SIMD.
2023 li_task->b0 = li->nc;
2025 while (con < ncon_tot && li->nc - li_task->b0 < ncon_target)
2027 if (li->con_index[con] == -1)
2032 type = iatom[3*con];
2033 a1 = iatom[3*con + 1];
2034 a2 = iatom[3*con + 2];
2035 lenA = idef->iparams[type].constr.dA;
2036 lenB = idef->iparams[type].constr.dB;
2037 /* Skip the flexible constraints when not doing dynamics */
2038 if (bDynamics || lenA != 0 || lenB != 0)
2040 assign_constraint(li, con, a1, a2, lenA, lenB, &at2con);
2042 if (li->ntask > 1 && !li->bTaskDep)
2044 /* We can generate independent tasks. Check if we
2045 * need to assign connected constraints to our task.
2047 check_assign_connected(li, iatom, idef, bDynamics,
2050 if (li->ntask > 1 && li->ncg_triangle > 0)
2052 /* Ensure constraints in one triangle are assigned
2055 check_assign_triangle(li, iatom, idef, bDynamics,
2056 con, a1, a2, &at2con);
2064 li_task->b1 = li->nc;
2068 /* Copy the last atom pair indices and lengths for constraints
2069 * up to a multiple of simd_width, such that we can do all
2070 * SIMD operations without having to worry about end effects.
2074 li->nc = ((li_task->b1 + simd_width - 1)/simd_width)*simd_width;
2075 last = li_task->b1 - 1;
2076 for (i = li_task->b1; i < li->nc; i++)
2078 li->bla[i*2 ] = li->bla[last*2 ];
2079 li->bla[i*2 + 1] = li->bla[last*2 + 1];
2080 li->bllen0[i] = li->bllen0[last];
2081 li->ddist[i] = li->ddist[last];
2082 li->bllen[i] = li->bllen[last];
2083 li->blnr[i + 1] = li->blnr[last + 1];
2087 /* Keep track of how many constraints we assigned */
2088 li->nc_real += li_task->b1 - li_task->b0;
2092 fprintf(debug, "LINCS task %d constraints %d - %d\n",
2093 th, li_task->b0, li_task->b1);
2097 assert(li->nc_real == ncon_assign);
2099 gmx_bool bSortMatrix;
2101 /* Without DD we order the blbnb matrix to optimize memory access.
2102 * With DD the overhead of sorting is more than the gain during access.
2104 bSortMatrix = !DOMAINDECOMP(cr);
2106 if (li->ncc > li->ncc_alloc)
2108 li->ncc_alloc = over_alloc_small(li->ncc);
2109 srenew(li->blbnb, li->ncc_alloc);
2112 #pragma omp parallel for num_threads(li->ntask) schedule(static)
2113 for (th = 0; th < li->ntask; th++)
2115 lincs_task_t *li_task;
2117 li_task = &li->task[th];
2119 if (li->ncg_triangle > 0 &&
2120 li_task->b1 - li_task->b0 > li_task->tri_alloc)
2122 /* This is allocating too much, but it is difficult to improve */
2123 li_task->tri_alloc = over_alloc_dd(li_task->b1 - li_task->b0);
2124 srenew(li_task->triangle, li_task->tri_alloc);
2125 srenew(li_task->tri_bits, li_task->tri_alloc);
2128 set_matrix_indices(li, li_task, &at2con, bSortMatrix);
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 (!gmx_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 = fabs(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 if (ir->eI == eiSD2 && v == NULL)
2343 lincsd->rmsd_data[i] = 0;
2349 if (econq == econqCoord)
2351 /* We can't use bCalcDHDL here, since NULL can be passed for dvdlambda
2352 * also with efep!=fepNO.
2354 if (ir->efep != efepNO)
2356 if (md->nMassPerturbed && lincsd->matlam != md->lambda)
2358 set_lincs_matrix(lincsd, md->invmass, md->lambda);
2361 for (i = 0; i < lincsd->nc; i++)
2363 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
2367 if (lincsd->ncg_flex)
2369 /* Set the flexible constraint lengths to the old lengths */
2372 for (i = 0; i < lincsd->nc; i++)
2374 if (lincsd->bllen[i] == 0)
2376 pbc_dx_aiuc(pbc, x[lincsd->bla[2*i]], x[lincsd->bla[2*i+1]], dx);
2377 lincsd->bllen[i] = norm(dx);
2383 for (i = 0; i < lincsd->nc; i++)
2385 if (lincsd->bllen[i] == 0)
2388 sqrt(distance2(x[lincsd->bla[2*i]],
2389 x[lincsd->bla[2*i+1]]));
2397 cconerr(lincsd, xprime, pbc,
2398 &ncons_loc, &p_ssd, &p_max, &p_imax);
2401 /* This bWarn var can be updated by multiple threads
2402 * at the same time. But as we only need to detect
2403 * if a warning occured or not, this is not an issue.
2407 /* The OpenMP parallel region of constrain_lincs for coords */
2408 #pragma omp parallel num_threads(lincsd->ntask)
2410 int th = gmx_omp_get_thread_num();
2412 clear_mat(lincsd->task[th].vir_r_m_dr);
2414 do_lincs(x, xprime, box, pbc, lincsd, th,
2417 ir->LincsWarnAngle, &bWarn,
2419 th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2422 if (bLog && fplog && lincsd->nc > 0)
2424 fprintf(fplog, " Rel. Constraint Deviation: RMS MAX between atoms\n");
2425 fprintf(fplog, " Before LINCS %.6f %.6f %6d %6d\n",
2426 sqrt(p_ssd/ncons_loc), p_max,
2427 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2428 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2432 cconerr(lincsd, xprime, pbc,
2433 &ncons_loc, &p_ssd, &p_max, &p_imax);
2434 /* Check if we are doing the second part of SD */
2435 if (ir->eI == eiSD2 && v == NULL)
2443 lincsd->rmsd_data[0] = ncons_loc;
2444 lincsd->rmsd_data[i] = p_ssd;
2448 lincsd->rmsd_data[0] = 0;
2449 lincsd->rmsd_data[1] = 0;
2450 lincsd->rmsd_data[2] = 0;
2452 if (bLog && fplog && lincsd->nc > 0)
2455 " After LINCS %.6f %.6f %6d %6d\n\n",
2456 sqrt(p_ssd/ncons_loc), p_max,
2457 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2458 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2465 cconerr(lincsd, xprime, pbc,
2466 &ncons_loc, &p_ssd, &p_max, &p_imax);
2469 sprintf(buf3, " in simulation %d", cr->ms->sim);
2475 sprintf(buf, "\nStep %s, time %g (ps) LINCS WARNING%s\n"
2476 "relative constraint deviation after LINCS:\n"
2477 "rms %.6f, max %.6f (between atoms %d and %d)\n",
2478 gmx_step_str(step, buf2), ir->init_t+step*ir->delta_t,
2480 sqrt(p_ssd/ncons_loc), p_max,
2481 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
2482 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
2485 fprintf(fplog, "%s", buf);
2487 fprintf(stderr, "%s", buf);
2488 lincs_warning(fplog, cr->dd, x, xprime, pbc,
2489 lincsd->nc, lincsd->bla, lincsd->bllen,
2490 ir->LincsWarnAngle, maxwarn, warncount);
2492 bOK = (p_max < 0.5);
2495 if (lincsd->ncg_flex)
2497 for (i = 0; (i < lincsd->nc); i++)
2499 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
2501 lincsd->bllen[i] = 0;
2508 /* The OpenMP parallel region of constrain_lincs for derivatives */
2509 #pragma omp parallel num_threads(lincsd->ntask)
2511 int th = gmx_omp_get_thread_num();
2513 do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
2514 md->invmass, econq, bCalcDHDL,
2515 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->task[th].vir_r_m_dr);
2521 /* Reduce the dH/dlambda contributions over the threads */
2526 for (th = 0; th < lincsd->ntask; th++)
2528 dhdlambda += lincsd->task[th].dhdlambda;
2532 /* dhdlambda contains dH/dlambda*dt^2, correct for this */
2533 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
2534 dhdlambda /= ir->delta_t*ir->delta_t;
2536 *dvdlambda += dhdlambda;
2539 if (bCalcVir && lincsd->ntask > 1)
2541 for (i = 1; i < lincsd->ntask; i++)
2543 m_add(vir_r_m_dr, lincsd->task[i].vir_r_m_dr, vir_r_m_dr);
2547 /* count assuming nit=1 */
2548 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc_real);
2549 inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
2550 if (lincsd->ntriangle > 0)
2552 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
2556 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc_real*2);
2560 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc_real);