1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
50 typedef struct gmx_shakedata
63 gmx_shakedata_t shake_init()
75 /* SOR initialization */
83 static void pv(FILE *log,char *s,rvec x)
87 fprintf(log,"%5s:",s);
88 for(m=0; (m<DIM); m++)
89 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;
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;
118 for (nit=0; (nit<maxnit) && (nconv != 0) && (error == 0); nit++) {
120 for(ll=0; (ll<ncon) && (error == 0); ll++) {
139 rpij2 = tx*tx+ty*ty+tz*tz;
143 /* iconvf is less than 1 when the error is smaller than a bound */
144 /* But if tt is too big, then it will result in looping in iconv */
146 iconvf = fabs(diff)*tt[ll];
150 rrpr = rijx*tx+rijy*ty+rijz*tz;
152 if (rrpr < toler*mytol)
155 acor = omega*diff*m2[ll]/rrpr;
176 int vec_shakef(FILE *fplog,gmx_shakedata_t shaked,
177 int natoms,real invmass[],int ncon,
178 t_iparams ip[],t_iatom *iatom,
179 real tol,rvec x[],rvec prime[],real omega,
180 gmx_bool bFEP,real lambda,real lagr[],
182 gmx_bool bCalcVir,tensor vir_r_m_dr,int econq,
188 int nit=0,ll,i,j,type;
193 real g,vscale,rscale,rvscale;
195 if (ncon > shaked->nalloc)
197 shaked->nalloc = over_alloc_dd(ncon);
198 srenew(shaked->rij,shaked->nalloc);
199 srenew(shaked->M2,shaked->nalloc);
200 srenew(shaked->tt,shaked->nalloc);
201 srenew(shaked->dist2,shaked->nalloc);
206 dist2 = shaked->dist2;
211 for(ll=0; (ll<ncon); ll++,ia+=3) {
216 mm=2*(invmass[i]+invmass[j]);
217 rij[ll][XX]=x[i][XX]-x[j][XX];
218 rij[ll][YY]=x[i][YY]-x[j][YY];
219 rij[ll][ZZ]=x[i][ZZ]-x[j][ZZ];
222 toler = sqr(L1*ip[type].constr.dA + lambda*ip[type].constr.dB);
224 toler = sqr(ip[type].constr.dA);
226 tt[ll] = 1.0/(toler*tol2);
231 cshake(iatom,ncon,&nit,maxnit,dist2,prime[0],rij[0],M2,omega,invmass,tt,lagr,&error);
234 crattle(iatom,ncon,&nit,maxnit,dist2,prime[0],rij[0],M2,omega,invmass,tt,lagr,&error,invdt,vetavar);
240 fprintf(fplog,"Shake did not converge in %d steps\n",maxnit);
242 fprintf(stderr,"Shake did not converge in %d steps\n",maxnit);
248 fprintf(fplog,"Inner product between old and new vector <= 0.0!\n"
249 "constraint #%d atoms %u and %u\n",
250 error-1,iatom[3*(error-1)+1]+1,iatom[3*(error-1)+2]+1);
251 fprintf(stderr,"Inner product between old and new vector <= 0.0!\n"
252 "constraint #%d atoms %u and %u\n",
253 error-1,iatom[3*(error-1)+1]+1,iatom[3*(error-1)+2]+1);
257 /* Constraint virial and correct the lagrange multipliers for the length */
261 for(ll=0; (ll<ncon); ll++,ia+=3)
264 if ((econq == econqCoord) && v!=NULL)
266 /* Correct the velocities */
267 mm = lagr[ll]*invmass[ia[1]]*invdt/vetavar->rscale;
270 v[ia[1]][i] += mm*rij[ll][i];
272 mm = lagr[ll]*invmass[ia[2]]*invdt/vetavar->rscale;
275 v[ia[2]][i] -= mm*rij[ll][i];
280 /* constraint virial */
283 if (econq == econqCoord)
285 mm = lagr[ll]/vetavar->rvscale;
287 if (econq == econqVeloc)
289 mm = lagr[ll]/(vetavar->vscale*vetavar->vscale_nhc[0]);
296 vir_r_m_dr[i][j] -= tmp*rij[ll][j];
302 /* Correct the lagrange multipliers for the length */
303 /* (more details would be useful here . . . )*/
308 toler = L1*ip[type].constr.dA + lambda*ip[type].constr.dB;
312 toler = ip[type].constr.dA;
320 static void check_cons(FILE *log,int nc,rvec x[],rvec prime[], rvec v[],
321 t_iparams ip[],t_iatom *iatom,
322 real invmass[], int econq)
331 " i mi j mj before after should be\n");
333 for(i=0; (i<nc); i++,ia+=3) {
336 rvec_sub(x[ai],x[aj],dx);
341 rvec_sub(prime[ai],prime[aj],dx);
343 fprintf(log,"%5d %5.2f %5d %5.2f %10.5f %10.5f %10.5f\n",
344 ai+1,1.0/invmass[ai],
345 aj+1,1.0/invmass[aj],d,dp,ip[ia[0]].constr.dA);
348 rvec_sub(v[ai],v[aj],dv);
350 rvec_sub(prime[ai],prime[aj],dv);
352 fprintf(log,"%5d %5.2f %5d %5.2f %10.5f %10.5f %10.5f\n",
353 ai+1,1.0/invmass[ai],
354 aj+1,1.0/invmass[aj],d,dp,0.);
360 gmx_bool bshakef(FILE *log,gmx_shakedata_t shaked,
361 int natoms,real invmass[],int nblocks,int sblock[],
362 t_idef *idef,t_inputrec *ir,rvec x_s[],rvec prime[],
363 t_nrnb *nrnb,real *lagr,real lambda,real *dvdlambda,
364 real invdt,rvec *v,gmx_bool bCalcVir,tensor vir_r_m_dr,
365 gmx_bool bDumpOnError,int econq,t_vetavars *vetavar)
369 int i,n0,ncons,blen,type;
373 fprintf(log,"nblocks=%d, sblock[0]=%d\n",nblocks,sblock[0]);
376 ncons=idef->il[F_CONSTR].nr/3;
378 for(i=0; i<ncons; i++)
381 iatoms = &(idef->il[F_CONSTR].iatoms[sblock[0]]);
383 for(i=0; (i<nblocks); ) {
384 blen = (sblock[i+1]-sblock[i]);
386 n0 = vec_shakef(log,shaked,natoms,invmass,blen,idef->iparams,
387 iatoms,ir->shake_tol,x_s,prime,shaked->omega,
388 ir->efep!=efepNO,lambda,lam,invdt,v,bCalcVir,vir_r_m_dr,
392 check_cons(log,blen,x_s,prime,v,idef->iparams,iatoms,invmass,econq);
396 if (bDumpOnError && log)
399 check_cons(log,blen,x_s,prime,v,idef->iparams,iatoms,invmass,econq);
406 iatoms += 3*blen; /* Increment pointer! */
410 /* only for position part? */
411 if (econq == econqCoord) {
412 if (ir->efep != efepNO) {
413 dt_2 = 1/sqr(ir->delta_t);
415 for(i=0; i<ncons; i++) {
416 type = idef->il[F_CONSTR].iatoms[3*i];
417 dvdl += lagr[i]*dt_2*
418 (idef->iparams[type].constr.dB-idef->iparams[type].constr.dA);
424 fprintf(log,"tnit: %5d omega: %10.5f\n",tnit,omega);
427 if (tnit > shaked->gamma) {
428 shaked->delta *= -0.5;
430 shaked->omega += shaked->delta;
431 shaked->gamma = tnit;
433 inc_nrnb(nrnb,eNR_SHAKE,tnit);
434 inc_nrnb(nrnb,eNR_SHAKE_RIJ,trij);
436 inc_nrnb(nrnb,eNR_CONSTR_V,trij*2);
438 inc_nrnb(nrnb,eNR_CONSTR_VIR,trij);
443 void crattle(atom_id iatom[],int ncon,int *nnit,int maxnit,
444 real dist2[],real vp[],real rij[],real m2[],real omega,
445 real invmass[],real tt[],real lagr[],int *nerror,real invdt,t_vetavars *vetavar)
448 * r.c. van schaik and w.f. van gunsteren
451 * Adapted for use with Gromacs by David van der Spoel november 92 and later.
452 * rattle added by M.R. Shirts, April 2004, from code written by Jay Ponder in TINKER
453 * second part of rattle algorithm
456 const real mytol=1e-10;
458 int ll,i,j,i3,j3,l3,ii;
459 int ix,iy,iz,jx,jy,jz;
460 real toler,rijd,vpijd, vx,vy,vz,diff,acor,xdotd,fac,im,jm,imdt,jmdt;
461 real xh,yh,zh,rijx,rijy,rijz;
465 real veta,vscale_nhc,iconvf;
467 veta = vetavar->veta;
468 vscale_nhc = vetavar->vscale_nhc[0]; /* for now, just use the first state */
472 for (nit=0; (nit<maxnit) && (nconv != 0) && (error == 0); nit++) {
474 for(ll=0; (ll<ncon) && (error == 0); ll++) {
493 vpijd = vx*rijx+vy*rijy+vz*rijz;
495 /* this is r(t+dt) \dotproduct \dot{r}(t+dt) */
496 xdotd = vpijd*vscale_nhc + veta*toler;
498 /* iconv is zero when the error is smaller than a bound */
499 iconvf = fabs(xdotd)*(tt[ll]/invdt);
503 fac = omega*2.0*m2[ll]/toler;
511 im = invmass[i]/vscale_nhc;
512 jm = invmass[j]/vscale_nhc;