Split txtdump.*, part 1
[gromacs.git] / src / gromacs / mdlib / shakef.cpp
blob5dad928503d4fbb45464b1dfc973f8bb6575c9d6
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, 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 #include "gmxpre.h"
39 #include <math.h>
41 #include "gromacs/gmxlib/nrnb.h"
42 #include "gromacs/math/vec.h"
43 #include "gromacs/mdlib/constr.h"
44 #include "gromacs/mdtypes/md_enums.h"
45 #include "gromacs/utility/smalloc.h"
47 typedef struct gmx_shakedata
49 rvec *rij;
50 real *half_of_reduced_mass;
51 real *distance_squared_tolerance;
52 real *constraint_distance_squared;
53 int nalloc;
54 /* SOR stuff */
55 real delta;
56 real omega;
57 real gamma;
58 } t_gmx_shakedata;
60 gmx_shakedata_t shake_init()
62 gmx_shakedata_t d;
64 snew(d, 1);
66 d->nalloc = 0;
67 d->rij = NULL;
68 d->half_of_reduced_mass = NULL;
69 d->distance_squared_tolerance = NULL;
70 d->constraint_distance_squared = NULL;
72 /* SOR initialization */
73 d->delta = 0.1;
74 d->omega = 1.0;
75 d->gamma = 1000000;
77 return d;
80 /*! \brief Inner kernel for SHAKE constraints
82 * Original implementation from R.C. van Schaik and W.F. van Gunsteren
83 * (ETH Zuerich, June 1992), adapted for GROMACS by David van der
84 * Spoel November 1992.
86 * The algorithm here is based section five of Ryckaert, Ciccotti and
87 * Berendsen, J Comp Phys, 23, 327, 1977.
89 * \param[in] iatom Mini-topology of triples of constraint type (unused in this
90 * function) and indices of the two atoms involved
91 * \param[in] ncon Number of constraints
92 * \param[out] nnit Number of iterations performed
93 * \param[in] maxnit Maximum number of iterations permitted
94 * \param[in] constraint_distance_squared The objective value for each constraint
95 * \param[inout] positions The initial (and final) values of the positions of all atoms
96 * \param[in] initial_displacements The initial displacements of each constraint
97 * \param[in] half_of_reduced_mass Half of the reduced mass for each constraint
98 * \param[in] omega SHAKE over-relaxation factor (set non-1.0 by
99 * using shake-sor=yes in the .mdp, but there is no documentation anywhere)
100 * \param[in] invmass Inverse mass of each atom
101 * \param[in] distance_squared_tolerance Multiplicative tolerance on the difference in the
102 * square of the constrained distance (see code)
103 * \param[out] scaled_lagrange_multiplier Scaled Lagrange multiplier for each constraint (-2 * eta from p. 336
104 * of the paper, divided by the constraint distance)
105 * \param[out] nerror Zero upon success, returns one more than the index of the
106 * problematic constraint if the input was malformed
108 * \todo Make SHAKE use better data structures, in particular for iatom. */
109 void cshake(const int iatom[], int ncon, int *nnit, int maxnit,
110 const real constraint_distance_squared[], real positions[],
111 const real initial_displacements[], const real half_of_reduced_mass[], real omega,
112 const real invmass[], const real distance_squared_tolerance[],
113 real scaled_lagrange_multiplier[], int *nerror)
115 /* default should be increased! MRS 8/4/2009 */
116 const real mytol = 1e-10;
118 int ll, i, j, i3, j3, l3;
119 int ix, iy, iz, jx, jy, jz;
120 real r_dot_r_prime;
121 real constraint_distance_squared_ll;
122 real r_prime_squared;
123 real scaled_lagrange_multiplier_ll;
124 real r_prime_x, r_prime_y, r_prime_z, diff, im, jm;
125 real xh, yh, zh, rijx, rijy, rijz;
126 int nit, error, nconv;
127 real iconvf;
129 // TODO nconv is used solely as a boolean, so we should write the
130 // code like that
131 error = 0;
132 nconv = 1;
133 for (nit = 0; (nit < maxnit) && (nconv != 0) && (error == 0); nit++)
135 nconv = 0;
136 for (ll = 0; (ll < ncon) && (error == 0); ll++)
138 l3 = 3*ll;
139 rijx = initial_displacements[l3+XX];
140 rijy = initial_displacements[l3+YY];
141 rijz = initial_displacements[l3+ZZ];
142 i = iatom[l3+1];
143 j = iatom[l3+2];
144 i3 = 3*i;
145 j3 = 3*j;
146 ix = i3+XX;
147 iy = i3+YY;
148 iz = i3+ZZ;
149 jx = j3+XX;
150 jy = j3+YY;
151 jz = j3+ZZ;
153 /* Compute r prime between atoms i and j, which is the
154 displacement *before* this update stage */
155 r_prime_x = positions[ix]-positions[jx];
156 r_prime_y = positions[iy]-positions[jy];
157 r_prime_z = positions[iz]-positions[jz];
158 r_prime_squared = (r_prime_x * r_prime_x +
159 r_prime_y * r_prime_y +
160 r_prime_z * r_prime_z);
161 constraint_distance_squared_ll = constraint_distance_squared[ll];
162 diff = constraint_distance_squared_ll - r_prime_squared;
164 /* iconvf is less than 1 when the error is smaller than a bound */
165 iconvf = fabs(diff) * distance_squared_tolerance[ll];
167 if (iconvf > 1.0)
169 nconv = static_cast<int>(iconvf);
170 r_dot_r_prime = (rijx * r_prime_x +
171 rijy * r_prime_y +
172 rijz * r_prime_z);
174 if (r_dot_r_prime < constraint_distance_squared_ll * mytol)
176 error = ll+1;
178 else
180 /* The next line solves equation 5.6 (neglecting
181 the term in g^2), for g */
182 scaled_lagrange_multiplier_ll = omega*diff*half_of_reduced_mass[ll]/r_dot_r_prime;
183 scaled_lagrange_multiplier[ll] += scaled_lagrange_multiplier_ll;
184 xh = rijx * scaled_lagrange_multiplier_ll;
185 yh = rijy * scaled_lagrange_multiplier_ll;
186 zh = rijz * scaled_lagrange_multiplier_ll;
187 im = invmass[i];
188 jm = invmass[j];
189 positions[ix] += xh*im;
190 positions[iy] += yh*im;
191 positions[iz] += zh*im;
192 positions[jx] -= xh*jm;
193 positions[jy] -= yh*jm;
194 positions[jz] -= zh*jm;
199 *nnit = nit;
200 *nerror = error;
203 int vec_shakef(FILE *fplog, gmx_shakedata_t shaked,
204 real invmass[], int ncon,
205 t_iparams ip[], t_iatom *iatom,
206 real tol, rvec x[], rvec prime[], real omega,
207 gmx_bool bFEP, real lambda, real scaled_lagrange_multiplier[],
208 real invdt, rvec *v,
209 gmx_bool bCalcVir, tensor vir_r_m_dr, int econq)
211 rvec *rij;
212 real *half_of_reduced_mass, *distance_squared_tolerance, *constraint_distance_squared;
213 int maxnit = 1000;
214 int nit = 0, ll, i, j, d, d2, type;
215 t_iatom *ia;
216 real L1;
217 real mm = 0., tmp;
218 int error = 0;
219 real constraint_distance;
221 if (ncon > shaked->nalloc)
223 shaked->nalloc = over_alloc_dd(ncon);
224 srenew(shaked->rij, shaked->nalloc);
225 srenew(shaked->half_of_reduced_mass, shaked->nalloc);
226 srenew(shaked->distance_squared_tolerance, shaked->nalloc);
227 srenew(shaked->constraint_distance_squared, shaked->nalloc);
229 rij = shaked->rij;
230 half_of_reduced_mass = shaked->half_of_reduced_mass;
231 distance_squared_tolerance = shaked->distance_squared_tolerance;
232 constraint_distance_squared = shaked->constraint_distance_squared;
234 L1 = 1.0-lambda;
235 ia = iatom;
236 for (ll = 0; (ll < ncon); ll++, ia += 3)
238 type = ia[0];
239 i = ia[1];
240 j = ia[2];
242 mm = 2.0*(invmass[i]+invmass[j]);
243 rij[ll][XX] = x[i][XX]-x[j][XX];
244 rij[ll][YY] = x[i][YY]-x[j][YY];
245 rij[ll][ZZ] = x[i][ZZ]-x[j][ZZ];
246 half_of_reduced_mass[ll] = 1.0/mm;
247 if (bFEP)
249 constraint_distance = L1*ip[type].constr.dA + lambda*ip[type].constr.dB;
251 else
253 constraint_distance = ip[type].constr.dA;
255 constraint_distance_squared[ll] = sqr(constraint_distance);
256 distance_squared_tolerance[ll] = 0.5/(constraint_distance_squared[ll]*tol);
259 switch (econq)
261 case econqCoord:
262 cshake(iatom, ncon, &nit, maxnit, constraint_distance_squared, prime[0], rij[0], half_of_reduced_mass, omega, invmass, distance_squared_tolerance, scaled_lagrange_multiplier, &error);
263 break;
264 case econqVeloc:
265 crattle(iatom, ncon, &nit, maxnit, constraint_distance_squared, prime[0], rij[0], half_of_reduced_mass, omega, invmass, distance_squared_tolerance, scaled_lagrange_multiplier, &error, invdt);
266 break;
269 if (nit >= maxnit)
271 if (fplog)
273 fprintf(fplog, "Shake did not converge in %d steps\n", maxnit);
275 fprintf(stderr, "Shake did not converge in %d steps\n", maxnit);
276 nit = 0;
278 else if (error != 0)
280 if (fplog)
282 fprintf(fplog, "Inner product between old and new vector <= 0.0!\n"
283 "constraint #%d atoms %d and %d\n",
284 error-1, iatom[3*(error-1)+1]+1, iatom[3*(error-1)+2]+1);
286 fprintf(stderr, "Inner product between old and new vector <= 0.0!\n"
287 "constraint #%d atoms %d and %d\n",
288 error-1, iatom[3*(error-1)+1]+1, iatom[3*(error-1)+2]+1);
289 nit = 0;
292 /* Constraint virial and correct the Lagrange multipliers for the length */
294 ia = iatom;
296 for (ll = 0; (ll < ncon); ll++, ia += 3)
298 type = ia[0];
299 i = ia[1];
300 j = ia[2];
302 if ((econq == econqCoord) && v != NULL)
304 /* Correct the velocities */
305 mm = scaled_lagrange_multiplier[ll]*invmass[i]*invdt;
306 for (d = 0; d < DIM; d++)
308 v[ia[1]][d] += mm*rij[ll][d];
310 mm = scaled_lagrange_multiplier[ll]*invmass[j]*invdt;
311 for (d = 0; d < DIM; d++)
313 v[ia[2]][d] -= mm*rij[ll][d];
315 /* 16 flops */
318 /* constraint virial */
319 if (bCalcVir)
321 mm = scaled_lagrange_multiplier[ll];
322 for (d = 0; d < DIM; d++)
324 tmp = mm*rij[ll][d];
325 for (d2 = 0; d2 < DIM; d2++)
327 vir_r_m_dr[d][d2] -= tmp*rij[ll][d2];
330 /* 21 flops */
333 /* cshake and crattle produce Lagrange multipliers scaled by
334 the reciprocal of the constraint length, so fix that */
335 if (bFEP)
337 constraint_distance = L1*ip[type].constr.dA + lambda*ip[type].constr.dB;
339 else
341 constraint_distance = ip[type].constr.dA;
343 scaled_lagrange_multiplier[ll] *= constraint_distance;
346 return nit;
349 static void check_cons(FILE *log, int nc, rvec x[], rvec prime[], rvec v[],
350 t_iparams ip[], t_iatom *iatom,
351 real invmass[], int econq)
353 t_iatom *ia;
354 int ai, aj;
355 int i;
356 real d, dp;
357 rvec dx, dv;
359 fprintf(log,
360 " i mi j mj before after should be\n");
361 ia = iatom;
362 for (i = 0; (i < nc); i++, ia += 3)
364 ai = ia[1];
365 aj = ia[2];
366 rvec_sub(x[ai], x[aj], dx);
367 d = norm(dx);
369 switch (econq)
371 case econqCoord:
372 rvec_sub(prime[ai], prime[aj], dx);
373 dp = norm(dx);
374 fprintf(log, "%5d %5.2f %5d %5.2f %10.5f %10.5f %10.5f\n",
375 ai+1, 1.0/invmass[ai],
376 aj+1, 1.0/invmass[aj], d, dp, ip[ia[0]].constr.dA);
377 break;
378 case econqVeloc:
379 rvec_sub(v[ai], v[aj], dv);
380 d = iprod(dx, dv);
381 rvec_sub(prime[ai], prime[aj], dv);
382 dp = iprod(dx, dv);
383 fprintf(log, "%5d %5.2f %5d %5.2f %10.5f %10.5f %10.5f\n",
384 ai+1, 1.0/invmass[ai],
385 aj+1, 1.0/invmass[aj], d, dp, 0.);
386 break;
391 gmx_bool bshakef(FILE *log, gmx_shakedata_t shaked,
392 real invmass[], int nblocks, int sblock[],
393 t_idef *idef, t_inputrec *ir, rvec x_s[], rvec prime[],
394 t_nrnb *nrnb, real *scaled_lagrange_multiplier, real lambda, real *dvdlambda,
395 real invdt, rvec *v, gmx_bool bCalcVir, tensor vir_r_m_dr,
396 gmx_bool bDumpOnError, int econq)
398 t_iatom *iatoms;
399 real dt_2, dvdl;
400 int i, n0, ncon, blen, type, ll;
401 int tnit = 0, trij = 0;
403 #ifdef DEBUG
404 fprintf(log, "nblocks=%d, sblock[0]=%d\n", nblocks, sblock[0]);
405 #endif
407 ncon = idef->il[F_CONSTR].nr/3;
409 for (ll = 0; ll < ncon; ll++)
411 scaled_lagrange_multiplier[ll] = 0;
414 iatoms = &(idef->il[F_CONSTR].iatoms[sblock[0]]);
415 for (i = 0; (i < nblocks); )
417 blen = (sblock[i+1]-sblock[i]);
418 blen /= 3;
419 n0 = vec_shakef(log, shaked, invmass, blen, idef->iparams,
420 iatoms, ir->shake_tol, x_s, prime, shaked->omega,
421 ir->efep != efepNO, lambda, scaled_lagrange_multiplier, invdt, v, bCalcVir, vir_r_m_dr,
422 econq);
424 #ifdef DEBUGSHAKE
425 check_cons(log, blen, x_s, prime, v, idef->iparams, iatoms, invmass, econq);
426 #endif
428 if (n0 == 0)
430 if (bDumpOnError && log)
433 check_cons(log, blen, x_s, prime, v, idef->iparams, iatoms, invmass, econq);
436 return FALSE;
438 tnit += n0*blen;
439 trij += blen;
440 iatoms += 3*blen; /* Increment pointer! */
441 scaled_lagrange_multiplier += blen;
442 i++;
444 /* only for position part? */
445 if (econq == econqCoord)
447 if (ir->efep != efepNO)
449 real bondA, bondB;
450 /* TODO This should probably use invdt, so that sd integrator scaling works properly */
451 dt_2 = 1/sqr(ir->delta_t);
452 dvdl = 0;
453 for (ll = 0; ll < ncon; ll++)
455 type = idef->il[F_CONSTR].iatoms[3*ll];
457 /* Per equations in the manual, dv/dl = -2 \sum_ll lagrangian_ll * r_ll * (d_B - d_A) */
458 /* The vector scaled_lagrange_multiplier[ll] contains the value -2 r_ll eta_ll (eta_ll is the
459 estimate of the Langrangian, definition on page 336 of Ryckaert et al 1977),
460 so the pre-factors are already present. */
461 bondA = idef->iparams[type].constr.dA;
462 bondB = idef->iparams[type].constr.dB;
463 dvdl += scaled_lagrange_multiplier[ll] * dt_2 * (bondB - bondA);
465 *dvdlambda += dvdl;
468 #ifdef DEBUG
469 fprintf(log, "tnit: %5d omega: %10.5f\n", tnit, omega);
470 #endif
471 if (ir->bShakeSOR)
473 if (tnit > shaked->gamma)
475 shaked->delta *= -0.5;
477 shaked->omega += shaked->delta;
478 shaked->gamma = tnit;
480 inc_nrnb(nrnb, eNR_SHAKE, tnit);
481 inc_nrnb(nrnb, eNR_SHAKE_RIJ, trij);
482 if (v)
484 inc_nrnb(nrnb, eNR_CONSTR_V, trij*2);
486 if (bCalcVir)
488 inc_nrnb(nrnb, eNR_CONSTR_VIR, trij);
491 return TRUE;
494 void crattle(int iatom[], int ncon, int *nnit, int maxnit,
495 real constraint_distance_squared[], real vp[], real rij[], real m2[], real omega,
496 real invmass[], real distance_squared_tolerance[], real scaled_lagrange_multiplier[],
497 int *nerror, real invdt)
500 * r.c. van schaik and w.f. van gunsteren
501 * eth zuerich
502 * june 1992
503 * Adapted for use with Gromacs by David van der Spoel november 92 and later.
504 * rattle added by M.R. Shirts, April 2004, from code written by Jay Ponder in TINKER
505 * second part of rattle algorithm
508 int ll, i, j, i3, j3, l3;
509 int ix, iy, iz, jx, jy, jz;
510 real constraint_distance_squared_ll;
511 real vpijd, vx, vy, vz, acor, fac, im, jm;
512 real xh, yh, zh, rijx, rijy, rijz;
513 int nit, error, nconv;
514 real iconvf;
516 // TODO nconv is used solely as a boolean, so we should write the
517 // code like that
518 error = 0;
519 nconv = 1;
520 for (nit = 0; (nit < maxnit) && (nconv != 0) && (error == 0); nit++)
522 nconv = 0;
523 for (ll = 0; (ll < ncon) && (error == 0); ll++)
525 l3 = 3*ll;
526 rijx = rij[l3+XX];
527 rijy = rij[l3+YY];
528 rijz = rij[l3+ZZ];
529 i = iatom[l3+1];
530 j = iatom[l3+2];
531 i3 = 3*i;
532 j3 = 3*j;
533 ix = i3+XX;
534 iy = i3+YY;
535 iz = i3+ZZ;
536 jx = j3+XX;
537 jy = j3+YY;
538 jz = j3+ZZ;
539 vx = vp[ix]-vp[jx];
540 vy = vp[iy]-vp[jy];
541 vz = vp[iz]-vp[jz];
543 vpijd = vx*rijx+vy*rijy+vz*rijz;
544 constraint_distance_squared_ll = constraint_distance_squared[ll];
546 /* iconv is zero when the error is smaller than a bound */
547 iconvf = fabs(vpijd)*(distance_squared_tolerance[ll]/invdt);
549 if (iconvf > 1)
551 nconv = static_cast<int>(iconvf);
552 fac = omega*2.0*m2[ll]/constraint_distance_squared_ll;
553 acor = -fac*vpijd;
554 scaled_lagrange_multiplier[ll] += acor;
555 xh = rijx*acor;
556 yh = rijy*acor;
557 zh = rijz*acor;
559 im = invmass[i];
560 jm = invmass[j];
562 vp[ix] += xh*im;
563 vp[iy] += yh*im;
564 vp[iz] += zh*im;
565 vp[jx] -= xh*jm;
566 vp[jy] -= yh*jm;
567 vp[jz] -= zh*jm;
571 *nnit = nit;
572 *nerror = error;