2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Defines SETTLE code.
41 * \author Berk Hess <hess@kth.se>
42 * \author Mark Abraham <mark.j.abraham@gmail.com>
43 * \ingroup module_mdlib
55 #include "gromacs/math/functions.h"
56 #include "gromacs/math/invertmatrix.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdlib/constr.h"
59 #include "gromacs/mdtypes/mdatom.h"
60 #include "gromacs/pbcutil/ishift.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/topology/idef.h"
66 #include "gromacs/topology/ifunc.h"
67 #include "gromacs/topology/mtop_util.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/gmxassert.h"
70 #include "gromacs/utility/smalloc.h"
96 settleparam_t massw; /* Parameters for SETTLE for coordinates */
97 settleparam_t mass1; /* Parameters with all masses 1, for forces */
99 int nsettle; /* The number of settles on our rank */
100 int* ow1; /* Index to OW1 atoms, size nsettle + SIMD padding */
101 int* hw2; /* Index to HW2 atoms, size nsettle + SIMD padding */
102 int* hw3; /* Index to HW3 atoms, size nsettle + SIMD padding */
103 real* virfac; /* Virial factor 0 or 1, size nsettle + SIMD pad. */
104 int nalloc; /* Allocation size of ow1, hw2, hw3, virfac */
106 bool bUseSimd; /* Use SIMD intrinsics code, if possible */
110 //! Initializes a projection matrix.
111 static void init_proj_matrix(real invmO, real invmH, real dOH, real dHH, matrix inverseCouplingMatrix)
113 /* We normalize the inverse masses with invmO for the matrix inversion.
114 * so we can keep using masses of almost zero for frozen particles,
115 * without running out of the float range in invertMatrix.
117 double invmORelative = 1.0;
118 double invmHRelative = invmH / static_cast<double>(invmO);
119 double distanceRatio = dHH / static_cast<double>(dOH);
121 /* Construct the constraint coupling matrix */
123 mat[0][0] = invmORelative + invmHRelative;
124 mat[0][1] = invmORelative * (1.0 - 0.5 * gmx::square(distanceRatio));
125 mat[0][2] = invmHRelative * 0.5 * distanceRatio;
126 mat[1][1] = mat[0][0];
127 mat[1][2] = mat[0][2];
128 mat[2][2] = invmHRelative + invmHRelative;
129 mat[1][0] = mat[0][1];
130 mat[2][0] = mat[0][2];
131 mat[2][1] = mat[1][2];
133 invertMatrix(mat, inverseCouplingMatrix);
135 msmul(inverseCouplingMatrix, 1 / invmO, inverseCouplingMatrix);
138 //! Initializes settle parameters.
139 static void settleparam_init(settleparam_t* p, real mO, real mH, real invmO, real invmH, real dOH, real dHH)
141 /* We calculate parameters in double precision to minimize errors.
142 * The velocity correction applied during SETTLE coordinate constraining
143 * introduces a systematic error of approximately 1 bit per atom,
144 * depending on what the compiler does with the code.
150 wohh = mO + 2.0 * mH;
154 double rc = dHH / 2.0;
155 double ra = 2.0 * mH * std::sqrt(dOH * dOH - rc * rc) / wohh;
156 p->rb = std::sqrt(dOH * dOH - rc * rc) - ra;
161 /* For projection: inverse masses and coupling matrix inversion */
165 p->invdOH = 1.0 / dOH;
166 p->invdHH = 1.0 / dHH;
168 init_proj_matrix(invmO, invmH, dOH, dHH, p->invmat);
172 fprintf(debug, "wh =%g, rc = %g, ra = %g\n", p->wh, p->rc, p->ra);
173 fprintf(debug, "rb = %g, irc2 = %g, dHH = %g, dOH = %g\n", p->rb, p->irc2, p->dHH, p->dOH);
177 settledata* settle_init(const gmx_mtop_t& mtop)
179 /* Check that we have only one settle type */
180 int settle_type = -1;
181 gmx_mtop_ilistloop_t iloop = gmx_mtop_ilistloop_init(mtop);
183 const int nral1 = 1 + NRAL(F_SETTLE);
184 while (const InteractionLists* ilists = gmx_mtop_ilistloop_next(iloop, &nmol))
186 const InteractionList& ilist = (*ilists)[F_SETTLE];
187 for (int i = 0; i < ilist.size(); i += nral1)
189 if (settle_type == -1)
191 settle_type = ilist.iatoms[i];
193 else if (ilist.iatoms[i] != settle_type)
196 "The [molecules] section of your topology specifies more than one block "
198 "a [moleculetype] with a [settles] block. Only one such is allowed.\n"
199 "If you are trying to partition your solvent into different *groups*\n"
200 "(e.g. for freezing, T-coupling, etc.), you are using the wrong "
202 "files specify groups. Otherwise, you may wish to change the least-used\n"
203 "block of molecules with SETTLE constraints into 3 normal constraints.");
207 GMX_RELEASE_ASSERT(settle_type >= 0, "settle_init called without settles");
213 /* We will not initialize the normal SETTLE parameters here yet,
214 * since the atom (inv)masses can depend on the integrator and
215 * free-energy perturbation. We set mO=-1 to trigger later initialization.
217 settled->massw.mO = -1;
219 real dOH = mtop.ffparams.iparams[settle_type].settle.doh;
220 real dHH = mtop.ffparams.iparams[settle_type].settle.dhh;
221 settleparam_init(&settled->mass1, 1.0, 1.0, 1.0, 1.0, dOH, dHH);
223 settled->ow1 = nullptr;
224 settled->hw2 = nullptr;
225 settled->hw3 = nullptr;
226 settled->virfac = nullptr;
229 /* Without SIMD configured, this bool is not used */
230 settled->bUseSimd = (getenv("GMX_DISABLE_SIMD_KERNELS") == nullptr);
235 void settle_free(settledata* settled)
237 sfree_aligned(settled->ow1);
238 sfree_aligned(settled->hw2);
239 sfree_aligned(settled->hw3);
240 sfree_aligned(settled->virfac);
244 void settle_set_constraints(settledata* settled, const t_ilist* il_settle, const t_mdatoms& mdatoms)
246 #if GMX_SIMD_HAVE_REAL
247 const int pack_size = GMX_SIMD_REAL_WIDTH;
249 const int pack_size = 1;
252 const int nral1 = 1 + NRAL(F_SETTLE);
253 int nsettle = il_settle->nr / nral1;
254 settled->nsettle = nsettle;
258 const t_iatom* iatoms = il_settle->iatoms;
260 /* Here we initialize the normal SETTLE parameters */
261 if (settled->massw.mO < 0)
263 int firstO = iatoms[1];
264 int firstH = iatoms[2];
265 settleparam_init(&settled->massw, mdatoms.massT[firstO], mdatoms.massT[firstH],
266 mdatoms.invmass[firstO], mdatoms.invmass[firstH], settled->mass1.dOH,
270 if (nsettle + pack_size > settled->nalloc)
272 settled->nalloc = over_alloc_dd(nsettle + pack_size);
273 sfree_aligned(settled->ow1);
274 sfree_aligned(settled->hw2);
275 sfree_aligned(settled->hw3);
276 sfree_aligned(settled->virfac);
277 snew_aligned(settled->ow1, settled->nalloc, 64);
278 snew_aligned(settled->hw2, settled->nalloc, 64);
279 snew_aligned(settled->hw3, settled->nalloc, 64);
280 snew_aligned(settled->virfac, settled->nalloc, 64);
283 for (int i = 0; i < nsettle; i++)
285 settled->ow1[i] = iatoms[i * nral1 + 1];
286 settled->hw2[i] = iatoms[i * nral1 + 2];
287 settled->hw3[i] = iatoms[i * nral1 + 3];
288 /* We should avoid double counting of virial contributions for
289 * SETTLEs that appear in multiple DD domains, so we only count
290 * the contribution on the home range of the oxygen atom.
292 settled->virfac[i] = (iatoms[i * nral1 + 1] < mdatoms.homenr ? 1 : 0);
295 /* Pack the index array to the full SIMD width with copies from
296 * the last normal entry, but with no virial contribution.
298 int end_packed = ((nsettle + pack_size - 1) / pack_size) * pack_size;
299 for (int i = nsettle; i < end_packed; i++)
301 settled->ow1[i] = settled->ow1[nsettle - 1];
302 settled->hw2[i] = settled->hw2[nsettle - 1];
303 settled->hw3[i] = settled->hw3[nsettle - 1];
304 settled->virfac[i] = 0;
309 void settle_proj(settledata* settled,
310 ConstraintVariable econq,
312 const t_iatom iatoms[],
317 int calcvir_atom_end,
320 /* Settle for projection out constraint components
321 * of derivatives of the coordinates.
322 * Berk Hess 2008-1-10
326 real imO, imH, dOH, dHH, invdOH, invdHH;
328 int i, m, m2, ow1, hw2, hw3;
329 rvec roh2, roh3, rhh, dc, fc;
331 calcvir_atom_end *= DIM;
333 if (econq == ConstraintVariable::Force)
343 copy_mat(p->invmat, invmat);
349 const int nral1 = 1 + NRAL(F_SETTLE);
351 for (i = 0; i < nsettle; i++)
353 ow1 = iatoms[i * nral1 + 1];
354 hw2 = iatoms[i * nral1 + 2];
355 hw3 = iatoms[i * nral1 + 3];
359 rvec_sub(x[ow1], x[hw2], roh2);
360 rvec_sub(x[ow1], x[hw3], roh3);
361 rvec_sub(x[hw2], x[hw3], rhh);
365 pbc_dx_aiuc(pbc, x[ow1], x[hw2], roh2);
366 pbc_dx_aiuc(pbc, x[ow1], x[hw3], roh3);
367 pbc_dx_aiuc(pbc, x[hw2], x[hw3], rhh);
369 svmul(invdOH, roh2, roh2);
370 svmul(invdOH, roh3, roh3);
371 svmul(invdHH, rhh, rhh);
374 /* Determine the projections of der on the bonds */
376 for (m = 0; m < DIM; m++)
378 dc[0] += (der[ow1][m] - der[hw2][m]) * roh2[m];
379 dc[1] += (der[ow1][m] - der[hw3][m]) * roh3[m];
380 dc[2] += (der[hw2][m] - der[hw3][m]) * rhh[m];
384 /* Determine the correction for the three bonds */
385 mvmul(invmat, dc, fc);
388 /* Subtract the corrections from derp */
389 for (m = 0; m < DIM; m++)
391 derp[ow1][m] -= imO * (fc[0] * roh2[m] + fc[1] * roh3[m]);
392 derp[hw2][m] -= imH * (-fc[0] * roh2[m] + fc[2] * rhh[m]);
393 derp[hw3][m] -= imH * (-fc[1] * roh3[m] - fc[2] * rhh[m]);
398 if (ow1 < calcvir_atom_end)
400 /* Determining r \dot m der is easy,
401 * since fc contains the mass weighted corrections for der.
404 for (m = 0; m < DIM; m++)
406 for (m2 = 0; m2 < DIM; m2++)
408 vir_r_m_dder[m][m2] += dOH * roh2[m] * roh2[m2] * fc[0]
409 + dOH * roh3[m] * roh3[m2] * fc[1]
410 + dHH * rhh[m] * rhh[m2] * fc[2];
418 /*! \brief The actual settle code, templated for real/SimdReal and for optimization */
419 template<typename T, typename TypeBool, int packSize, typename TypePbc, bool bCorrectVelocity, bool bCalcVirial>
420 static void settleTemplate(const settledata* settled,
427 real* gmx_restrict v,
429 bool* bErrorHasOccurred)
431 /* ******************************************************************* */
433 /* Original code by Shuichi Miyamoto, last update Oct. 1, 1992 ** */
435 /* Algorithm changes by Berk Hess: ** */
436 /* 2004-07-15 Convert COM to double precision to avoid drift ** */
437 /* 2006-10-16 Changed velocity update to use differences ** */
438 /* 2012-09-24 Use oxygen as reference instead of COM ** */
439 /* 2016-02 Complete rewrite of the code for SIMD ** */
441 /* Reference for the SETTLE algorithm ** */
442 /* S. Miyamoto et al., J. Comp. Chem., 13, 952 (1992). ** */
444 /* ******************************************************************* */
446 assert(settleStart % packSize == 0);
447 assert(settleEnd % packSize == 0);
449 TypeBool bError = TypeBool(false);
451 const settleparam_t* p = &settled->massw;
460 T almost_zero = T(1e-12);
462 T sum_r_m_dr[DIM][DIM];
466 for (int d2 = 0; d2 < DIM; d2++)
468 for (int d = 0; d < DIM; d++)
470 sum_r_m_dr[d2][d] = T(0);
475 for (int i = settleStart; i < settleEnd; i += packSize)
477 /* Here we pad up to packSize with copies from the last valid entry.
478 * This gives correct results, since we store (not increment) all
479 * output, so we store the same output multiple times.
481 const int* ow1 = settled->ow1 + i;
482 const int* hw2 = settled->hw2 + i;
483 const int* hw3 = settled->hw3 + i;
485 T x_ow1[DIM], x_hw2[DIM], x_hw3[DIM];
487 gatherLoadUTranspose<3>(x, ow1, &x_ow1[XX], &x_ow1[YY], &x_ow1[ZZ]);
488 gatherLoadUTranspose<3>(x, hw2, &x_hw2[XX], &x_hw2[YY], &x_hw2[ZZ]);
489 gatherLoadUTranspose<3>(x, hw3, &x_hw3[XX], &x_hw3[YY], &x_hw3[ZZ]);
491 T xprime_ow1[DIM], xprime_hw2[DIM], xprime_hw3[DIM];
493 gatherLoadUTranspose<3>(xprime, ow1, &xprime_ow1[XX], &xprime_ow1[YY], &xprime_ow1[ZZ]);
494 gatherLoadUTranspose<3>(xprime, hw2, &xprime_hw2[XX], &xprime_hw2[YY], &xprime_hw2[ZZ]);
495 gatherLoadUTranspose<3>(xprime, hw3, &xprime_hw3[XX], &xprime_hw3[YY], &xprime_hw3[ZZ]);
497 T dist21[DIM], dist31[DIM];
498 T doh2[DIM], doh3[DIM];
499 T sh_hw2[DIM], sh_hw3[DIM];
501 pbc_dx_aiuc(pbc, x_hw2, x_ow1, dist21);
503 pbc_dx_aiuc(pbc, x_hw3, x_ow1, dist31);
505 /* Tedious way of doing pbc */
506 pbc_dx_aiuc(pbc, xprime_hw2, xprime_ow1, doh2);
507 for (int d = 0; d < DIM; d++)
509 sh_hw2[d] = xprime_hw2[d] - (xprime_ow1[d] + doh2[d]);
510 xprime_hw2[d] = xprime_hw2[d] - sh_hw2[d];
512 pbc_dx_aiuc(pbc, xprime_hw3, xprime_ow1, doh3);
513 for (int d = 0; d < DIM; d++)
515 sh_hw3[d] = xprime_hw3[d] - (xprime_ow1[d] + doh3[d]);
516 xprime_hw3[d] = xprime_hw3[d] - sh_hw3[d];
519 /* Not calculating the center of mass using the oxygen position
520 * and the O-H distances, as done below, will make SETTLE
521 * the largest source of energy drift for simulations of water,
522 * as then the oxygen coordinate is multiplied by 0.89 at every step,
523 * which can then transfer a systematic rounding to the oxygen velocity.
526 for (int d = 0; d < DIM; d++)
528 a1[d] = -(doh2[d] + doh3[d]) * wh;
529 com[d] = xprime_ow1[d] - a1[d];
532 for (int d = 0; d < DIM; d++)
534 b1[d] = xprime_hw2[d] - com[d];
537 for (int d = 0; d < DIM; d++)
539 c1[d] = xprime_hw3[d] - com[d];
543 T xakszd = dist21[YY] * dist31[ZZ] - dist21[ZZ] * dist31[YY];
544 T yakszd = dist21[ZZ] * dist31[XX] - dist21[XX] * dist31[ZZ];
545 T zakszd = dist21[XX] * dist31[YY] - dist21[YY] * dist31[XX];
546 T xaksxd = a1[YY] * zakszd - a1[ZZ] * yakszd;
547 T yaksxd = a1[ZZ] * xakszd - a1[XX] * zakszd;
548 T zaksxd = a1[XX] * yakszd - a1[YY] * xakszd;
549 T xaksyd = yakszd * zaksxd - zakszd * yaksxd;
550 T yaksyd = zakszd * xaksxd - xakszd * zaksxd;
551 T zaksyd = xakszd * yaksxd - yakszd * xaksxd;
554 T axlng = gmx::invsqrt(xaksxd * xaksxd + yaksxd * yaksxd + zaksxd * zaksxd);
555 T aylng = gmx::invsqrt(xaksyd * xaksyd + yaksyd * yaksyd + zaksyd * zaksyd);
556 T azlng = gmx::invsqrt(xakszd * xakszd + yakszd * yakszd + zakszd * zakszd);
558 T trns1[DIM], trns2[DIM], trns3[DIM];
560 trns1[XX] = xaksxd * axlng;
561 trns2[XX] = yaksxd * axlng;
562 trns3[XX] = zaksxd * axlng;
563 trns1[YY] = xaksyd * aylng;
564 trns2[YY] = yaksyd * aylng;
565 trns3[YY] = zaksyd * aylng;
566 trns1[ZZ] = xakszd * azlng;
567 trns2[ZZ] = yakszd * azlng;
568 trns3[ZZ] = zakszd * azlng;
573 for (int d = 0; d < 2; d++)
575 b0d[d] = trns1[d] * dist21[XX] + trns2[d] * dist21[YY] + trns3[d] * dist21[ZZ];
576 c0d[d] = trns1[d] * dist31[XX] + trns2[d] * dist31[YY] + trns3[d] * dist31[ZZ];
579 T a1d_z, b1d[DIM], c1d[DIM];
581 a1d_z = trns1[ZZ] * a1[XX] + trns2[ZZ] * a1[YY] + trns3[ZZ] * a1[ZZ];
582 for (int d = 0; d < DIM; d++)
584 b1d[d] = trns1[d] * b1[XX] + trns2[d] * b1[YY] + trns3[d] * b1[ZZ];
585 c1d[d] = trns1[d] * c1[XX] + trns2[d] * c1[YY] + trns3[d] * c1[ZZ];
591 T sinphi = a1d_z * gmx::invsqrt(ra * ra);
592 tmp2 = 1.0 - sinphi * sinphi;
594 /* If tmp2 gets close to or beyond zero we have severly distorted
595 * water molecules and we should terminate the simulation.
596 * Below we take the max with almost_zero to continue the loop.
598 bError = bError || (tmp2 <= almost_zero);
600 tmp2 = max(tmp2, almost_zero);
601 tmp = gmx::invsqrt(tmp2);
602 T cosphi = tmp2 * tmp;
603 T sinpsi = (b1d[ZZ] - c1d[ZZ]) * irc2 * tmp;
604 tmp2 = 1.0 - sinpsi * sinpsi;
606 T cospsi = tmp2 * gmx::invsqrt(tmp2);
609 T a2d_y = ra * cosphi;
610 T b2d_x = -rc * cospsi;
612 T t2 = rc * sinpsi * sinphi;
617 /* --- Step3 al,be,ga --- */
618 T alpha = b2d_x * (b0d[XX] - c0d[XX]) + b0d[YY] * b2d_y + c0d[YY] * c2d_y;
619 T beta = b2d_x * (c0d[YY] - b0d[YY]) + b0d[XX] * b2d_y + c0d[XX] * c2d_y;
620 T gamma = b0d[XX] * b1d[YY] - b1d[XX] * b0d[YY] + c0d[XX] * c1d[YY] - c1d[XX] * c0d[YY];
621 T al2be2 = alpha * alpha + beta * beta;
622 tmp2 = (al2be2 - gamma * gamma);
623 T sinthe = (alpha * gamma - beta * tmp2 * gmx::invsqrt(tmp2)) * gmx::invsqrt(al2be2 * al2be2);
626 /* --- Step4 A3' --- */
627 tmp2 = 1.0 - sinthe * sinthe;
628 T costhe = tmp2 * gmx::invsqrt(tmp2);
630 T a3d[DIM], b3d[DIM], c3d[DIM];
632 a3d[XX] = -a2d_y * sinthe;
633 a3d[YY] = a2d_y * costhe;
635 b3d[XX] = b2d_x * costhe - b2d_y * sinthe;
636 b3d[YY] = b2d_x * sinthe + b2d_y * costhe;
638 c3d[XX] = -b2d_x * costhe - c2d_y * sinthe;
639 c3d[YY] = -b2d_x * sinthe + c2d_y * costhe;
643 /* --- Step5 A3 --- */
644 T a3[DIM], b3[DIM], c3[DIM];
646 a3[XX] = trns1[XX] * a3d[XX] + trns1[YY] * a3d[YY] + trns1[ZZ] * a3d[ZZ];
647 a3[YY] = trns2[XX] * a3d[XX] + trns2[YY] * a3d[YY] + trns2[ZZ] * a3d[ZZ];
648 a3[ZZ] = trns3[XX] * a3d[XX] + trns3[YY] * a3d[YY] + trns3[ZZ] * a3d[ZZ];
649 b3[XX] = trns1[XX] * b3d[XX] + trns1[YY] * b3d[YY] + trns1[ZZ] * b3d[ZZ];
650 b3[YY] = trns2[XX] * b3d[XX] + trns2[YY] * b3d[YY] + trns2[ZZ] * b3d[ZZ];
651 b3[ZZ] = trns3[XX] * b3d[XX] + trns3[YY] * b3d[YY] + trns3[ZZ] * b3d[ZZ];
652 c3[XX] = trns1[XX] * c3d[XX] + trns1[YY] * c3d[YY] + trns1[ZZ] * c3d[ZZ];
653 c3[YY] = trns2[XX] * c3d[XX] + trns2[YY] * c3d[YY] + trns2[ZZ] * c3d[ZZ];
654 c3[ZZ] = trns3[XX] * c3d[XX] + trns3[YY] * c3d[YY] + trns3[ZZ] * c3d[ZZ];
657 /* Compute and store the corrected new coordinate */
658 for (int d = 0; d < DIM; d++)
660 xprime_ow1[d] = com[d] + a3[d];
662 for (int d = 0; d < DIM; d++)
664 xprime_hw2[d] = com[d] + b3[d] + sh_hw2[d];
666 for (int d = 0; d < DIM; d++)
668 xprime_hw3[d] = com[d] + c3[d] + sh_hw3[d];
670 /* 9 flops + 6 pbc flops */
672 transposeScatterStoreU<3>(xprime, ow1, xprime_ow1[XX], xprime_ow1[YY], xprime_ow1[ZZ]);
673 transposeScatterStoreU<3>(xprime, hw2, xprime_hw2[XX], xprime_hw2[YY], xprime_hw2[ZZ]);
674 transposeScatterStoreU<3>(xprime, hw3, xprime_hw3[XX], xprime_hw3[YY], xprime_hw3[ZZ]);
676 if (bCorrectVelocity || bCalcVirial)
678 T da[DIM], db[DIM], dc[DIM];
679 for (int d = 0; d < DIM; d++)
681 da[d] = a3[d] - a1[d];
683 for (int d = 0; d < DIM; d++)
685 db[d] = b3[d] - b1[d];
687 for (int d = 0; d < DIM; d++)
689 dc[d] = c3[d] - c1[d];
693 if (bCorrectVelocity)
695 T v_ow1[DIM], v_hw2[DIM], v_hw3[DIM];
697 gatherLoadUTranspose<3>(v, ow1, &v_ow1[XX], &v_ow1[YY], &v_ow1[ZZ]);
698 gatherLoadUTranspose<3>(v, hw2, &v_hw2[XX], &v_hw2[YY], &v_hw2[ZZ]);
699 gatherLoadUTranspose<3>(v, hw3, &v_hw3[XX], &v_hw3[YY], &v_hw3[ZZ]);
701 /* Add the position correction divided by dt to the velocity */
702 for (int d = 0; d < DIM; d++)
704 v_ow1[d] = gmx::fma(da[d], invdt, v_ow1[d]);
706 for (int d = 0; d < DIM; d++)
708 v_hw2[d] = gmx::fma(db[d], invdt, v_hw2[d]);
710 for (int d = 0; d < DIM; d++)
712 v_hw3[d] = gmx::fma(dc[d], invdt, v_hw3[d]);
716 transposeScatterStoreU<3>(v, ow1, v_ow1[XX], v_ow1[YY], v_ow1[ZZ]);
717 transposeScatterStoreU<3>(v, hw2, v_hw2[XX], v_hw2[YY], v_hw2[ZZ]);
718 transposeScatterStoreU<3>(v, hw3, v_hw3[XX], v_hw3[YY], v_hw3[ZZ]);
723 /* Filter out the non-local settles */
724 T filter = load<T>(settled->virfac + i);
728 T mdo[DIM], mdb[DIM], mdc[DIM];
730 for (int d = 0; d < DIM; d++)
732 mdb[d] = mHf * db[d];
733 mdc[d] = mHf * dc[d];
734 mdo[d] = mOf * da[d] + mdb[d] + mdc[d];
737 for (int d2 = 0; d2 < DIM; d2++)
739 for (int d = 0; d < DIM; d++)
743 - (x_ow1[d2] * mdo[d] + dist21[d2] * mdb[d] + dist31[d2] * mdc[d]);
753 for (int d2 = 0; d2 < DIM; d2++)
755 for (int d = 0; d < DIM; d++)
757 vir_r_m_dr[d2][d] += reduce(sum_r_m_dr[d2][d]);
762 *bErrorHasOccurred = anyTrue(bError);
765 /*! \brief Wrapper template function that divides the settles over threads
766 * and instantiates the core template with instantiated booleans.
768 template<typename T, typename TypeBool, int packSize, typename TypePbc>
769 static void settleTemplateWrapper(settledata* settled,
779 bool* bErrorHasOccurred)
781 /* We need to assign settles to threads in groups of pack_size */
782 int numSettlePacks = (settled->nsettle + packSize - 1) / packSize;
783 /* Round the end value up to give thread 0 more work */
784 int settleStart = ((numSettlePacks * thread + nthread - 1) / nthread) * packSize;
785 int settleEnd = ((numSettlePacks * (thread + 1) + nthread - 1) / nthread) * packSize;
791 settleTemplate<T, TypeBool, packSize, TypePbc, true, false>(
792 settled, settleStart, settleEnd, pbc, x, xprime, invdt, v, nullptr, bErrorHasOccurred);
796 settleTemplate<T, TypeBool, packSize, TypePbc, true, true>(
797 settled, settleStart, settleEnd, pbc, x, xprime, invdt, v, vir_r_m_dr, bErrorHasOccurred);
804 settleTemplate<T, TypeBool, packSize, TypePbc, false, false>(
805 settled, settleStart, settleEnd, pbc, x, xprime, invdt, v, nullptr, bErrorHasOccurred);
809 settleTemplate<T, TypeBool, packSize, TypePbc, false, true>(
810 settled, settleStart, settleEnd, pbc, x, xprime, invdt, v, vir_r_m_dr, bErrorHasOccurred);
815 void csettle(settledata* settled,
825 bool* bErrorHasOccurred)
827 #if GMX_SIMD_HAVE_REAL
828 if (settled->bUseSimd)
830 /* Convert the pbc struct for SIMD */
831 alignas(GMX_SIMD_ALIGNMENT) real pbcSimd[9 * GMX_SIMD_REAL_WIDTH];
832 set_pbc_simd(pbc, pbcSimd);
834 settleTemplateWrapper<SimdReal, SimdBool, GMX_SIMD_REAL_WIDTH, const real*>(
835 settled, nthread, thread, pbcSimd, x, xprime, invdt, v, bCalcVirial, vir_r_m_dr,
841 /* This construct is needed because pbc_dx_aiuc doesn't accept pbc=NULL */
843 const t_pbc* pbcNonNull;
851 set_pbc(&pbcNo, PbcType::No, nullptr);
855 settleTemplateWrapper<real, bool, 1, const t_pbc*>(settled, nthread, thread, pbcNonNull, x,
856 xprime, invdt, v, bCalcVirial,
857 vir_r_m_dr, bErrorHasOccurred);