Update instructions in containers.rst
[gromacs.git] / src / gromacs / mdlib / shake.cpp
blobee5dee5f4f82ea38d7e3e045ba3582c2d1bb270e
1 /*
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,2017,2018 by the GROMACS development team.
7 * Copyright (c) 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.
38 /*! \internal \file
39 * \brief Defines SHAKE code.
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \author Berk Hess <hess@kth.se>
43 * \author Mark Abraham <mark.j.abraham@gmail.com>
44 * \ingroup module_mdlib
46 #include "gmxpre.h"
48 #include "shake.h"
50 #include <cmath>
51 #include <cstdlib>
53 #include <algorithm>
55 #include "gromacs/gmxlib/nrnb.h"
56 #include "gromacs/math/functions.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdlib/constr.h"
59 #include "gromacs/mdlib/splitter.h"
60 #include "gromacs/mdtypes/inputrec.h"
61 #include "gromacs/mdtypes/md_enums.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/topology/invblock.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/smalloc.h"
67 namespace gmx
70 typedef struct
72 int iatom[3];
73 int blocknr;
74 } t_sortblock;
76 //! Compares sort blocks.
77 static int pcomp(const void* p1, const void* p2)
79 int db;
80 int min1, min2, max1, max2;
81 const t_sortblock* a1 = reinterpret_cast<const t_sortblock*>(p1);
82 const t_sortblock* a2 = reinterpret_cast<const t_sortblock*>(p2);
84 db = a1->blocknr - a2->blocknr;
86 if (db != 0)
88 return db;
91 min1 = std::min(a1->iatom[1], a1->iatom[2]);
92 max1 = std::max(a1->iatom[1], a1->iatom[2]);
93 min2 = std::min(a2->iatom[1], a2->iatom[2]);
94 max2 = std::max(a2->iatom[1], a2->iatom[2]);
96 if (min1 == min2)
98 return max1 - max2;
100 else
102 return min1 - min2;
106 //! Prints sortblocks
107 static void pr_sortblock(FILE* fp, const char* title, int nsb, t_sortblock sb[])
109 int i;
111 fprintf(fp, "%s\n", title);
112 for (i = 0; (i < nsb); i++)
114 fprintf(fp, "i: %5d, iatom: (%5d %5d %5d), blocknr: %5d\n", i, sb[i].iatom[0],
115 sb[i].iatom[1], sb[i].iatom[2], sb[i].blocknr);
119 //! Reallocates a vector.
120 static void resizeLagrangianData(shakedata* shaked, int ncons)
122 shaked->scaled_lagrange_multiplier.resize(ncons);
125 void make_shake_sblock_serial(shakedata* shaked, InteractionDefinitions* idef, const int numAtoms)
127 int i, m, ncons;
128 int bstart, bnr;
129 t_blocka sblocks;
130 t_sortblock* sb;
131 int* inv_sblock;
133 /* Since we are processing the local topology,
134 * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
136 ncons = idef->il[F_CONSTR].size() / 3;
138 init_blocka(&sblocks);
139 sfree(sblocks.index); // To solve memory leak
140 gen_sblocks(nullptr, numAtoms, *idef, &sblocks, FALSE);
143 bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
144 nblocks=blocks->multinr[idef->nodeid] - bstart;
146 bstart = 0;
147 if (debug)
149 fprintf(debug, "ncons: %d, bstart: %d, nblocks: %d\n", ncons, bstart, sblocks.nr);
152 /* Calculate block number for each atom */
153 inv_sblock = make_invblocka(&sblocks, numAtoms);
155 done_blocka(&sblocks);
157 /* Store the block number in temp array and
158 * sort the constraints in order of the sblock number
159 * and the atom numbers, really sorting a segment of the array!
161 int* iatom = idef->il[F_CONSTR].iatoms.data();
162 snew(sb, ncons);
163 for (i = 0; (i < ncons); i++, iatom += 3)
165 for (m = 0; (m < 3); m++)
167 sb[i].iatom[m] = iatom[m];
169 sb[i].blocknr = inv_sblock[iatom[1]];
172 /* Now sort the blocks */
173 if (debug)
175 pr_sortblock(debug, "Before sorting", ncons, sb);
176 fprintf(debug, "Going to sort constraints\n");
179 std::qsort(sb, ncons, sizeof(*sb), pcomp);
181 if (debug)
183 pr_sortblock(debug, "After sorting", ncons, sb);
186 iatom = idef->il[F_CONSTR].iatoms.data();
187 for (i = 0; (i < ncons); i++, iatom += 3)
189 for (m = 0; (m < 3); m++)
191 iatom[m] = sb[i].iatom[m];
195 shaked->sblock.clear();
196 bnr = -2;
197 for (i = 0; (i < ncons); i++)
199 if (sb[i].blocknr != bnr)
201 bnr = sb[i].blocknr;
202 shaked->sblock.push_back(3 * i);
205 /* Last block... */
206 shaked->sblock.push_back(3 * ncons);
208 sfree(sb);
209 sfree(inv_sblock);
210 resizeLagrangianData(shaked, ncons);
213 void make_shake_sblock_dd(shakedata* shaked, const InteractionList& ilcon)
215 int ncons, c, cg;
217 ncons = ilcon.size() / 3;
218 const int* iatom = ilcon.iatoms.data();
219 shaked->sblock.clear();
220 cg = 0;
221 for (c = 0; c < ncons; c++)
223 if (c == 0 || iatom[1] >= cg + 1)
225 shaked->sblock.push_back(3 * c);
226 while (iatom[1] >= cg + 1)
228 cg++;
231 iatom += 3;
233 shaked->sblock.push_back(3 * ncons);
234 resizeLagrangianData(shaked, ncons);
237 /*! \brief Inner kernel for SHAKE constraints
239 * Original implementation from R.C. van Schaik and W.F. van Gunsteren
240 * (ETH Zuerich, June 1992), adapted for GROMACS by David van der
241 * Spoel November 1992.
243 * The algorithm here is based section five of Ryckaert, Ciccotti and
244 * Berendsen, J Comp Phys, 23, 327, 1977.
246 * \param[in] iatom Mini-topology of triplets of constraint type (unused
247 * in this function) and indices of two atoms involved
248 * \param[in] ncon Number of constraints
249 * \param[out] nnit Number of iterations performed
250 * \param[in] maxnit Maximum number of iterations permitted
251 * \param[in] constraint_distance_squared The objective value for each constraint
252 * \param[inout] positions The initial (and final) values of the positions
253 * of all atoms
254 * \param[in] pbc PBC information
255 * \param[in] initial_displacements The initial displacements of each constraint
256 * \param[in] half_of_reduced_mass Half of the reduced mass for each constraint
257 * \param[in] omega SHAKE over-relaxation factor (set non-1.0 by
258 * using shake-sor=yes in the .mdp,
259 * but there is no documentation anywhere)
260 * \param[in] invmass Inverse mass of each atom
261 * \param[in] distance_squared_tolerance Multiplicative tolerance on the difference in the
262 * square of the constrained distance (see code)
263 * \param[out] scaled_lagrange_multiplier Scaled Lagrange multiplier for each constraint
264 * (-2 * eta from p. 336 of the paper, divided by
265 * the constraint distance)
266 * \param[out] nerror Zero upon success, returns one more than the index of
267 * the problematic constraint if the input was malformed
269 * \todo Make SHAKE use better data structures, in particular for iatom. */
270 void cshake(const int iatom[],
271 int ncon,
272 int* nnit,
273 int maxnit,
274 ArrayRef<const real> constraint_distance_squared,
275 ArrayRef<RVec> positions,
276 const t_pbc* pbc,
277 ArrayRef<const RVec> initial_displacements,
278 ArrayRef<const real> half_of_reduced_mass,
279 real omega,
280 const real invmass[],
281 ArrayRef<const real> distance_squared_tolerance,
282 ArrayRef<real> scaled_lagrange_multiplier,
283 int* nerror)
285 /* default should be increased! MRS 8/4/2009 */
286 const real mytol = 1e-10;
288 // TODO nconv is used solely as a boolean, so we should write the
289 // code like that
290 int error = 0;
291 int nconv = 1;
292 int nit;
293 for (nit = 0; (nit < maxnit) && (nconv != 0) && (error == 0); nit++)
295 nconv = 0;
296 for (int ll = 0; (ll < ncon) && (error == 0); ll++)
298 const int l3 = 3 * ll;
299 const real rijx = initial_displacements[ll][XX];
300 const real rijy = initial_displacements[ll][YY];
301 const real rijz = initial_displacements[ll][ZZ];
302 const int i = iatom[l3 + 1];
303 const int j = iatom[l3 + 2];
305 /* Compute r prime between atoms i and j, which is the
306 displacement *before* this update stage */
307 rvec r_prime;
308 if (pbc)
310 pbc_dx(pbc, positions[i], positions[j], r_prime);
312 else
314 rvec_sub(positions[i], positions[j], r_prime);
316 const real r_prime_squared = norm2(r_prime);
317 const real constraint_distance_squared_ll = constraint_distance_squared[ll];
318 const real diff = constraint_distance_squared_ll - r_prime_squared;
320 /* iconvf is less than 1 when the error is smaller than a bound */
321 const real iconvf = std::abs(diff) * distance_squared_tolerance[ll];
323 if (iconvf > 1.0_real)
325 nconv = static_cast<int>(iconvf);
326 const real r_dot_r_prime =
327 (rijx * r_prime[XX] + rijy * r_prime[YY] + rijz * r_prime[ZZ]);
329 if (r_dot_r_prime < constraint_distance_squared_ll * mytol)
331 error = ll + 1;
333 else
335 /* The next line solves equation 5.6 (neglecting
336 the term in g^2), for g */
337 real scaled_lagrange_multiplier_ll =
338 omega * diff * half_of_reduced_mass[ll] / r_dot_r_prime;
339 scaled_lagrange_multiplier[ll] += scaled_lagrange_multiplier_ll;
340 const real xh = rijx * scaled_lagrange_multiplier_ll;
341 const real yh = rijy * scaled_lagrange_multiplier_ll;
342 const real zh = rijz * scaled_lagrange_multiplier_ll;
343 const real im = invmass[i];
344 const real jm = invmass[j];
345 positions[i][XX] += xh * im;
346 positions[i][YY] += yh * im;
347 positions[i][ZZ] += zh * im;
348 positions[j][XX] -= xh * jm;
349 positions[j][YY] -= yh * jm;
350 positions[j][ZZ] -= zh * jm;
355 *nnit = nit;
356 *nerror = error;
359 //! Implements RATTLE (ie. SHAKE for velocity verlet integrators)
360 static void crattle(const int iatom[],
361 int ncon,
362 int* nnit,
363 int maxnit,
364 ArrayRef<const real> constraint_distance_squared,
365 ArrayRef<RVec> vp,
366 ArrayRef<const RVec> rij,
367 ArrayRef<const real> m2,
368 real omega,
369 const real invmass[],
370 ArrayRef<const real> distance_squared_tolerance,
371 ArrayRef<real> scaled_lagrange_multiplier,
372 int* nerror,
373 real invdt)
376 * r.c. van schaik and w.f. van gunsteren
377 * eth zuerich
378 * june 1992
379 * Adapted for use with Gromacs by David van der Spoel november 92 and later.
380 * rattle added by M.R. Shirts, April 2004, from code written by Jay Ponder in TINKER
381 * second part of rattle algorithm
384 // TODO nconv is used solely as a boolean, so we should write the
385 // code like that
386 int error = 0;
387 int nconv = 1;
388 int nit;
389 for (nit = 0; (nit < maxnit) && (nconv != 0) && (error == 0); nit++)
391 nconv = 0;
392 for (int ll = 0; (ll < ncon) && (error == 0); ll++)
394 const int l3 = 3 * ll;
395 const real rijx = rij[ll][XX];
396 const real rijy = rij[ll][YY];
397 const real rijz = rij[ll][ZZ];
398 const int i = iatom[l3 + 1];
399 const int j = iatom[l3 + 2];
400 rvec v;
401 rvec_sub(vp[i], vp[j], v);
403 const real vpijd = v[XX] * rijx + v[YY] * rijy + v[ZZ] * rijz;
404 const real constraint_distance_squared_ll = constraint_distance_squared[ll];
406 /* iconv is zero when the error is smaller than a bound */
407 const real iconvf = fabs(vpijd) * (distance_squared_tolerance[ll] / invdt);
409 if (iconvf > 1.0_real)
411 nconv = static_cast<int>(iconvf);
412 const real fac = omega * 2.0_real * m2[ll] / constraint_distance_squared_ll;
413 const real acor = -fac * vpijd;
414 scaled_lagrange_multiplier[ll] += acor;
415 const real xh = rijx * acor;
416 const real yh = rijy * acor;
417 const real zh = rijz * acor;
419 const real im = invmass[i];
420 const real jm = invmass[j];
422 vp[i][XX] += xh * im;
423 vp[i][YY] += yh * im;
424 vp[i][ZZ] += zh * im;
425 vp[j][XX] -= xh * jm;
426 vp[j][YY] -= yh * jm;
427 vp[j][ZZ] -= zh * jm;
431 *nnit = nit;
432 *nerror = error;
435 //! Applies SHAKE
436 static int vec_shakef(FILE* fplog,
437 shakedata* shaked,
438 const real invmass[],
439 int ncon,
440 ArrayRef<const t_iparams> ip,
441 const int* iatom,
442 real tol,
443 ArrayRef<const RVec> x,
444 ArrayRef<RVec> prime,
445 const t_pbc* pbc,
446 real omega,
447 bool bFEP,
448 real lambda,
449 ArrayRef<real> scaled_lagrange_multiplier,
450 real invdt,
451 ArrayRef<RVec> v,
452 bool bCalcVir,
453 tensor vir_r_m_dr,
454 ConstraintVariable econq)
456 int maxnit = 1000;
457 int nit = 0, ll, i, j, d, d2, type;
458 real L1;
459 int error = 0;
460 real constraint_distance;
462 shaked->rij.resize(ncon);
463 shaked->half_of_reduced_mass.resize(ncon);
464 shaked->distance_squared_tolerance.resize(ncon);
465 shaked->constraint_distance_squared.resize(ncon);
467 ArrayRef<RVec> rij = shaked->rij;
468 ArrayRef<real> half_of_reduced_mass = shaked->half_of_reduced_mass;
469 ArrayRef<real> distance_squared_tolerance = shaked->distance_squared_tolerance;
470 ArrayRef<real> constraint_distance_squared = shaked->constraint_distance_squared;
472 L1 = 1.0_real - lambda;
473 const int* ia = iatom;
474 for (ll = 0; (ll < ncon); ll++, ia += 3)
476 type = ia[0];
477 i = ia[1];
478 j = ia[2];
480 if (pbc)
482 pbc_dx(pbc, x[i], x[j], rij[ll]);
484 else
486 rvec_sub(x[i], x[j], rij[ll]);
488 const real mm = 2.0_real * (invmass[i] + invmass[j]);
489 half_of_reduced_mass[ll] = 1.0_real / mm;
490 if (bFEP)
492 constraint_distance = L1 * ip[type].constr.dA + lambda * ip[type].constr.dB;
494 else
496 constraint_distance = ip[type].constr.dA;
498 constraint_distance_squared[ll] = gmx::square(constraint_distance);
499 distance_squared_tolerance[ll] = 0.5 / (constraint_distance_squared[ll] * tol);
502 switch (econq)
504 case ConstraintVariable::Positions:
505 cshake(iatom, ncon, &nit, maxnit, constraint_distance_squared, prime, pbc, rij,
506 half_of_reduced_mass, omega, invmass, distance_squared_tolerance,
507 scaled_lagrange_multiplier, &error);
508 break;
509 case ConstraintVariable::Velocities:
510 crattle(iatom, ncon, &nit, maxnit, constraint_distance_squared, prime, rij,
511 half_of_reduced_mass, omega, invmass, distance_squared_tolerance,
512 scaled_lagrange_multiplier, &error, invdt);
513 break;
514 default: gmx_incons("Unknown constraint quantity for SHAKE");
517 if (nit >= maxnit)
519 if (fplog)
521 fprintf(fplog, "Shake did not converge in %d steps\n", maxnit);
523 fprintf(stderr, "Shake did not converge in %d steps\n", maxnit);
524 nit = 0;
526 else if (error != 0)
528 if (fplog)
530 fprintf(fplog,
531 "Inner product between old and new vector <= 0.0!\n"
532 "constraint #%d atoms %d and %d\n",
533 error - 1, iatom[3 * (error - 1) + 1] + 1, iatom[3 * (error - 1) + 2] + 1);
535 fprintf(stderr,
536 "Inner product between old and new vector <= 0.0!\n"
537 "constraint #%d atoms %d and %d\n",
538 error - 1, iatom[3 * (error - 1) + 1] + 1, iatom[3 * (error - 1) + 2] + 1);
539 nit = 0;
542 /* Constraint virial and correct the Lagrange multipliers for the length */
544 ia = iatom;
546 for (ll = 0; (ll < ncon); ll++, ia += 3)
548 type = ia[0];
549 i = ia[1];
550 j = ia[2];
552 if ((econq == ConstraintVariable::Positions) && !v.empty())
554 /* Correct the velocities */
555 real mm = scaled_lagrange_multiplier[ll] * invmass[i] * invdt;
556 for (d = 0; d < DIM; d++)
558 v[ia[1]][d] += mm * rij[ll][d];
560 mm = scaled_lagrange_multiplier[ll] * invmass[j] * invdt;
561 for (d = 0; d < DIM; d++)
563 v[ia[2]][d] -= mm * rij[ll][d];
565 /* 16 flops */
568 /* constraint virial */
569 if (bCalcVir)
571 const real mm = scaled_lagrange_multiplier[ll];
572 for (d = 0; d < DIM; d++)
574 const real tmp = mm * rij[ll][d];
575 for (d2 = 0; d2 < DIM; d2++)
577 vir_r_m_dr[d][d2] -= tmp * rij[ll][d2];
580 /* 21 flops */
583 /* cshake and crattle produce Lagrange multipliers scaled by
584 the reciprocal of the constraint length, so fix that */
585 if (bFEP)
587 constraint_distance = L1 * ip[type].constr.dA + lambda * ip[type].constr.dB;
589 else
591 constraint_distance = ip[type].constr.dA;
593 scaled_lagrange_multiplier[ll] *= constraint_distance;
596 return nit;
599 //! Check that constraints are satisfied.
600 static void check_cons(FILE* log,
601 int nc,
602 ArrayRef<const RVec> x,
603 ArrayRef<const RVec> prime,
604 ArrayRef<const RVec> v,
605 const t_pbc* pbc,
606 ArrayRef<const t_iparams> ip,
607 const int* iatom,
608 const real invmass[],
609 ConstraintVariable econq)
611 int ai, aj;
612 int i;
613 real d, dp;
614 rvec dx, dv;
616 GMX_ASSERT(!v.empty(), "Input has to be non-null");
617 fprintf(log, " i mi j mj before after should be\n");
618 const int* ia = iatom;
619 for (i = 0; (i < nc); i++, ia += 3)
621 ai = ia[1];
622 aj = ia[2];
623 rvec_sub(x[ai], x[aj], dx);
624 d = norm(dx);
626 switch (econq)
628 case ConstraintVariable::Positions:
629 if (pbc)
631 pbc_dx(pbc, prime[ai], prime[aj], dx);
633 else
635 rvec_sub(prime[ai], prime[aj], dx);
637 dp = norm(dx);
638 fprintf(log, "%5d %5.2f %5d %5.2f %10.5f %10.5f %10.5f\n", ai + 1,
639 1.0 / invmass[ai], aj + 1, 1.0 / invmass[aj], d, dp, ip[ia[0]].constr.dA);
640 break;
641 case ConstraintVariable::Velocities:
642 rvec_sub(v[ai], v[aj], dv);
643 d = iprod(dx, dv);
644 rvec_sub(prime[ai], prime[aj], dv);
645 dp = iprod(dx, dv);
646 fprintf(log, "%5d %5.2f %5d %5.2f %10.5f %10.5f %10.5f\n", ai + 1,
647 1.0 / invmass[ai], aj + 1, 1.0 / invmass[aj], d, dp, 0.);
648 break;
649 default: gmx_incons("Unknown constraint quantity for SHAKE");
654 //! Applies SHAKE.
655 static bool bshakef(FILE* log,
656 shakedata* shaked,
657 const real invmass[],
658 const InteractionDefinitions& idef,
659 const t_inputrec& ir,
660 ArrayRef<const RVec> x_s,
661 ArrayRef<RVec> prime,
662 const t_pbc* pbc,
663 t_nrnb* nrnb,
664 real lambda,
665 real* dvdlambda,
666 real invdt,
667 ArrayRef<RVec> v,
668 bool bCalcVir,
669 tensor vir_r_m_dr,
670 bool bDumpOnError,
671 ConstraintVariable econq)
673 real dt_2, dvdl;
674 int i, n0, ncon, blen, type, ll;
675 int tnit = 0, trij = 0;
677 ncon = idef.il[F_CONSTR].size() / 3;
679 for (ll = 0; ll < ncon; ll++)
681 shaked->scaled_lagrange_multiplier[ll] = 0;
684 // TODO Rewrite this block so that it is obvious that i, iatoms
685 // and lam are all iteration variables. Is this easier if the
686 // sblock data structure is organized differently?
687 const int* iatoms = &(idef.il[F_CONSTR].iatoms[shaked->sblock[0]]);
688 ArrayRef<real> lam = shaked->scaled_lagrange_multiplier;
689 for (i = 0; (i < shaked->numShakeBlocks());)
691 blen = (shaked->sblock[i + 1] - shaked->sblock[i]);
692 blen /= 3;
693 n0 = vec_shakef(log, shaked, invmass, blen, idef.iparams, iatoms, ir.shake_tol, x_s, prime,
694 pbc, shaked->omega, ir.efep != efepNO, lambda, lam, invdt, v, bCalcVir,
695 vir_r_m_dr, econq);
697 if (n0 == 0)
699 if (bDumpOnError && log)
702 check_cons(log, blen, x_s, prime, v, pbc, idef.iparams, iatoms, invmass, econq);
705 return FALSE;
707 tnit += n0 * blen;
708 trij += blen;
709 iatoms += 3 * blen; /* Increment pointer! */
710 lam = lam.subArray(blen, lam.ssize() - blen);
711 i++;
713 /* only for position part? */
714 if (econq == ConstraintVariable::Positions)
716 if (ir.efep != efepNO)
718 ArrayRef<const t_iparams> iparams = idef.iparams;
720 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
721 dt_2 = 1 / gmx::square(ir.delta_t);
722 dvdl = 0;
723 for (ll = 0; ll < ncon; ll++)
725 type = idef.il[F_CONSTR].iatoms[3 * ll];
727 /* Per equations in the manual, dv/dl = -2 \sum_ll lagrangian_ll * r_ll * (d_B - d_A) */
728 /* The vector scaled_lagrange_multiplier[ll] contains the value -2 r_ll eta_ll
729 (eta_ll is the estimate of the Langrangian, definition on page 336 of Ryckaert et
730 al 1977), so the pre-factors are already present. */
731 const real bondA = iparams[type].constr.dA;
732 const real bondB = iparams[type].constr.dB;
733 dvdl += shaked->scaled_lagrange_multiplier[ll] * dt_2 * (bondB - bondA);
735 *dvdlambda += dvdl;
738 if (ir.bShakeSOR)
740 if (tnit > shaked->gamma)
742 shaked->delta *= -0.5;
744 shaked->omega += shaked->delta;
745 shaked->gamma = tnit;
747 inc_nrnb(nrnb, eNR_SHAKE, tnit);
748 inc_nrnb(nrnb, eNR_SHAKE_RIJ, trij);
749 if (!v.empty())
751 inc_nrnb(nrnb, eNR_CONSTR_V, trij * 2);
753 if (bCalcVir)
755 inc_nrnb(nrnb, eNR_CONSTR_VIR, trij);
758 return TRUE;
761 bool constrain_shake(FILE* log,
762 shakedata* shaked,
763 const real invmass[],
764 const InteractionDefinitions& idef,
765 const t_inputrec& ir,
766 ArrayRef<const RVec> x_s,
767 ArrayRef<RVec> xprime,
768 ArrayRef<RVec> vprime,
769 const t_pbc* pbc,
770 t_nrnb* nrnb,
771 real lambda,
772 real* dvdlambda,
773 real invdt,
774 ArrayRef<RVec> v,
775 bool bCalcVir,
776 tensor vir_r_m_dr,
777 bool bDumpOnError,
778 ConstraintVariable econq)
780 if (shaked->numShakeBlocks() == 0)
782 return true;
784 bool bOK;
785 switch (econq)
787 case (ConstraintVariable::Positions):
788 bOK = bshakef(log, shaked, invmass, idef, ir, x_s, xprime, pbc, nrnb, lambda, dvdlambda,
789 invdt, v, bCalcVir, vir_r_m_dr, bDumpOnError, econq);
790 break;
791 case (ConstraintVariable::Velocities):
792 bOK = bshakef(log, shaked, invmass, idef, ir, x_s, vprime, pbc, nrnb, lambda, dvdlambda,
793 invdt, {}, bCalcVir, vir_r_m_dr, bDumpOnError, econq);
794 break;
795 default:
796 gmx_fatal(FARGS,
797 "Internal error, SHAKE called for constraining something else than "
798 "coordinates");
800 return bOK;
803 } // namespace gmx