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, 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.
41 #include "gromacs/legacyheaders/constr.h"
42 #include "gromacs/legacyheaders/nrnb.h"
43 #include "gromacs/legacyheaders/txtdump.h"
44 #include "gromacs/legacyheaders/typedefs.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/utility/smalloc.h"
48 typedef struct gmx_shakedata
61 gmx_shakedata_t shake_init()
73 /* SOR initialization */
81 static void pv(FILE *log, char *s, rvec x)
85 fprintf(log, "%5s:", s);
86 for (m = 0; (m < DIM); m++)
88 fprintf(log, " %10.3f", x[m]);
94 void cshake(atom_id iatom[], int ncon, int *nnit, int maxnit,
95 real dist2[], real xp[], real rij[], real m2[], real omega,
96 real invmass[], real tt[], real lagr[], int *nerror)
99 * r.c. van schaik and w.f. van gunsteren
102 * Adapted for use with Gromacs by David van der Spoel november 92 and later.
104 /* default should be increased! MRS 8/4/2009 */
105 const real mytol = 1e-10;
107 int ll, i, j, i3, j3, l3;
108 int ix, iy, iz, jx, jy, jz;
109 real toler, rpij2, rrpr, tx, ty, tz, diff, acor, im, jm;
110 real xh, yh, zh, rijx, rijy, rijz;
113 int nit, error, nconv;
118 for (nit = 0; (nit < maxnit) && (nconv != 0) && (error == 0); nit++)
121 for (ll = 0; (ll < ncon) && (error == 0); ll++)
141 rpij2 = tx*tx+ty*ty+tz*tz;
145 /* iconvf is less than 1 when the error is smaller than a bound */
146 /* But if tt is too big, then it will result in looping in iconv */
148 iconvf = fabs(diff)*tt[ll];
153 rrpr = rijx*tx+rijy*ty+rijz*tz;
155 if (rrpr < toler*mytol)
161 acor = omega*diff*m2[ll]/rrpr;
182 int vec_shakef(FILE *fplog, gmx_shakedata_t shaked,
183 real invmass[], int ncon,
184 t_iparams ip[], t_iatom *iatom,
185 real tol, rvec x[], rvec prime[], real omega,
186 gmx_bool bFEP, real lambda, real lagr[],
188 gmx_bool bCalcVir, tensor vir_r_m_dr, int econq,
192 real *M2, *tt, *dist2;
194 int nit = 0, ll, i, j, type;
196 real L1, tol2, toler;
199 real g, vscale, rscale, rvscale;
201 if (ncon > shaked->nalloc)
203 shaked->nalloc = over_alloc_dd(ncon);
204 srenew(shaked->rij, shaked->nalloc);
205 srenew(shaked->M2, shaked->nalloc);
206 srenew(shaked->tt, shaked->nalloc);
207 srenew(shaked->dist2, shaked->nalloc);
212 dist2 = shaked->dist2;
217 for (ll = 0; (ll < ncon); ll++, ia += 3)
223 mm = 2*(invmass[i]+invmass[j]);
224 rij[ll][XX] = x[i][XX]-x[j][XX];
225 rij[ll][YY] = x[i][YY]-x[j][YY];
226 rij[ll][ZZ] = x[i][ZZ]-x[j][ZZ];
230 toler = sqr(L1*ip[type].constr.dA + lambda*ip[type].constr.dB);
234 toler = sqr(ip[type].constr.dA);
237 tt[ll] = 1.0/(toler*tol2);
243 cshake(iatom, ncon, &nit, maxnit, dist2, prime[0], rij[0], M2, omega, invmass, tt, lagr, &error);
246 crattle(iatom, ncon, &nit, maxnit, dist2, prime[0], rij[0], M2, omega, invmass, tt, lagr, &error, invdt, vetavar);
254 fprintf(fplog, "Shake did not converge in %d steps\n", maxnit);
256 fprintf(stderr, "Shake did not converge in %d steps\n", maxnit);
263 fprintf(fplog, "Inner product between old and new vector <= 0.0!\n"
264 "constraint #%d atoms %u and %u\n",
265 error-1, iatom[3*(error-1)+1]+1, iatom[3*(error-1)+2]+1);
267 fprintf(stderr, "Inner product between old and new vector <= 0.0!\n"
268 "constraint #%d atoms %u and %u\n",
269 error-1, iatom[3*(error-1)+1]+1, iatom[3*(error-1)+2]+1);
273 /* Constraint virial and correct the lagrange multipliers for the length */
277 for (ll = 0; (ll < ncon); ll++, ia += 3)
280 if ((econq == econqCoord) && v != NULL)
282 /* Correct the velocities */
283 mm = lagr[ll]*invmass[ia[1]]*invdt/vetavar->rscale;
284 for (i = 0; i < DIM; i++)
286 v[ia[1]][i] += mm*rij[ll][i];
288 mm = lagr[ll]*invmass[ia[2]]*invdt/vetavar->rscale;
289 for (i = 0; i < DIM; i++)
291 v[ia[2]][i] -= mm*rij[ll][i];
296 /* constraint virial */
299 if (econq == econqCoord)
301 mm = lagr[ll]/vetavar->rvscale;
303 if (econq == econqVeloc)
305 mm = lagr[ll]/(vetavar->vscale*vetavar->vscale_nhc[0]);
307 for (i = 0; i < DIM; i++)
310 for (j = 0; j < DIM; j++)
312 vir_r_m_dr[i][j] -= tmp*rij[ll][j];
318 /* Correct the lagrange multipliers for the length */
319 /* (more details would be useful here . . . )*/
324 toler = L1*ip[type].constr.dA + lambda*ip[type].constr.dB;
328 toler = ip[type].constr.dA;
336 static void check_cons(FILE *log, int nc, rvec x[], rvec prime[], rvec v[],
337 t_iparams ip[], t_iatom *iatom,
338 real invmass[], int econq)
347 " i mi j mj before after should be\n");
349 for (i = 0; (i < nc); i++, ia += 3)
353 rvec_sub(x[ai], x[aj], dx);
359 rvec_sub(prime[ai], prime[aj], dx);
361 fprintf(log, "%5d %5.2f %5d %5.2f %10.5f %10.5f %10.5f\n",
362 ai+1, 1.0/invmass[ai],
363 aj+1, 1.0/invmass[aj], d, dp, ip[ia[0]].constr.dA);
366 rvec_sub(v[ai], v[aj], dv);
368 rvec_sub(prime[ai], prime[aj], dv);
370 fprintf(log, "%5d %5.2f %5d %5.2f %10.5f %10.5f %10.5f\n",
371 ai+1, 1.0/invmass[ai],
372 aj+1, 1.0/invmass[aj], d, dp, 0.);
378 gmx_bool bshakef(FILE *log, gmx_shakedata_t shaked,
379 real invmass[], int nblocks, int sblock[],
380 t_idef *idef, t_inputrec *ir, rvec x_s[], rvec prime[],
381 t_nrnb *nrnb, real *lagr, real lambda, real *dvdlambda,
382 real invdt, rvec *v, gmx_bool bCalcVir, tensor vir_r_m_dr,
383 gmx_bool bDumpOnError, int econq, t_vetavars *vetavar)
386 real *lam, dt_2, dvdl;
387 int i, n0, ncons, blen, type;
388 int tnit = 0, trij = 0;
391 fprintf(log, "nblocks=%d, sblock[0]=%d\n", nblocks, sblock[0]);
394 ncons = idef->il[F_CONSTR].nr/3;
396 for (i = 0; i < ncons; i++)
401 iatoms = &(idef->il[F_CONSTR].iatoms[sblock[0]]);
403 for (i = 0; (i < nblocks); )
405 blen = (sblock[i+1]-sblock[i]);
407 n0 = vec_shakef(log, shaked, invmass, blen, idef->iparams,
408 iatoms, ir->shake_tol, x_s, prime, shaked->omega,
409 ir->efep != efepNO, lambda, lam, invdt, v, bCalcVir, vir_r_m_dr,
413 check_cons(log, blen, x_s, prime, v, idef->iparams, iatoms, invmass, econq);
418 if (bDumpOnError && log)
421 check_cons(log, blen, x_s, prime, v, idef->iparams, iatoms, invmass, econq);
428 iatoms += 3*blen; /* Increment pointer! */
432 /* only for position part? */
433 if (econq == econqCoord)
435 if (ir->efep != efepNO)
438 dt_2 = 1/sqr(ir->delta_t);
440 for (i = 0; i < ncons; i++)
442 type = idef->il[F_CONSTR].iatoms[3*i];
444 /* dh/dl contribution from constraint force is dh/dr (constraint force) dot dr/dl */
445 /* constraint force is -\sum_i lagr_i* d(constraint)/dr, with constrant = r^2-d^2 */
446 /* constraint force is -\sum_i lagr_i* 2 r */
447 /* so dh/dl = -\sum_i lagr_i* 2 r * dr/dl */
448 /* However, by comparison with lincs and with
449 comparison with a full thermodynamics cycle (see
450 redmine issue #1255), this is off by a factor of
451 two -- the 2r should apparently just be r. Further
452 investigation should be done at some point to
453 understand why and see if there is something deeper
456 bondA = idef->iparams[type].constr.dA;
457 bondB = idef->iparams[type].constr.dB;
458 dvdl += lagr[i] * dt_2 * ((1.0-lambda)*bondA + lambda*bondB) * (bondB-bondA);
464 fprintf(log, "tnit: %5d omega: %10.5f\n", tnit, omega);
468 if (tnit > shaked->gamma)
470 shaked->delta *= -0.5;
472 shaked->omega += shaked->delta;
473 shaked->gamma = tnit;
475 inc_nrnb(nrnb, eNR_SHAKE, tnit);
476 inc_nrnb(nrnb, eNR_SHAKE_RIJ, trij);
479 inc_nrnb(nrnb, eNR_CONSTR_V, trij*2);
483 inc_nrnb(nrnb, eNR_CONSTR_VIR, trij);
489 void crattle(atom_id iatom[], int ncon, int *nnit, int maxnit,
490 real dist2[], real vp[], real rij[], real m2[], real omega,
491 real invmass[], real tt[], real lagr[], int *nerror, real invdt, t_vetavars *vetavar)
494 * r.c. van schaik and w.f. van gunsteren
497 * Adapted for use with Gromacs by David van der Spoel november 92 and later.
498 * rattle added by M.R. Shirts, April 2004, from code written by Jay Ponder in TINKER
499 * second part of rattle algorithm
502 const real mytol = 1e-10;
504 int ll, i, j, i3, j3, l3, ii;
505 int ix, iy, iz, jx, jy, jz;
506 real toler, rijd, vpijd, vx, vy, vz, diff, acor, xdotd, fac, im, jm, imdt, jmdt;
507 real xh, yh, zh, rijx, rijy, rijz;
510 int nit, error, nconv;
511 real veta, vscale_nhc, iconvf;
513 veta = vetavar->veta;
514 vscale_nhc = vetavar->vscale_nhc[0]; /* for now, just use the first state */
518 for (nit = 0; (nit < maxnit) && (nconv != 0) && (error == 0); nit++)
521 for (ll = 0; (ll < ncon) && (error == 0); ll++)
541 vpijd = vx*rijx+vy*rijy+vz*rijz;
543 /* this is r(t+dt) \dotproduct \dot{r}(t+dt) */
544 xdotd = vpijd*vscale_nhc + veta*toler;
546 /* iconv is zero when the error is smaller than a bound */
547 iconvf = fabs(xdotd)*(tt[ll]/invdt);
552 fac = omega*2.0*m2[ll]/toler;
560 im = invmass[i]/vscale_nhc;
561 jm = invmass[j]/vscale_nhc;