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
54 #include "gmx_random.h"
68 #include "gmx_wallcycle.h"
70 /*For debugging, start at v(-dt/2) for velolcity verlet -- uncomment next line */
71 /*#define STARTFROMDT2*/
91 /* The random state */
97 gmx_sd_sigma_t *sdsig;
102 typedef struct gmx_update
107 /* Variables for the deform algorithm */
108 gmx_large_int_t deformref_step;
109 matrix deformref_box;
113 static void do_update_md(int start,int nrend,double dt,
114 t_grp_tcstat *tcstat,t_grp_acc *gstat,double nh_vxi[],
115 rvec accel[],ivec nFreeze[],real invmass[],
116 unsigned short ptype[],unsigned short cFREEZE[],
117 unsigned short cACC[],unsigned short cTC[],
118 rvec x[],rvec xprime[],rvec v[],
120 gmx_bool bNH,gmx_bool bPR)
125 real vn,vv,va,vb,vnrel;
131 /* Update with coupling to extended ensembles, used for
132 * Nose-Hoover and Parrinello-Rahman coupling
133 * Nose-Hoover uses the reversible leap-frog integrator from
134 * Holian et al. Phys Rev E 52(3) : 2338, 1995
136 for(n=start; n<nrend; n++)
151 lg = tcstat[gt].lambda;
155 rvec_sub(v[n],gstat[ga].u,vrel);
159 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
161 vnrel = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
162 - iprod(M[d],vrel)))/(1 + 0.5*vxi*dt);
163 /* do not scale the mean velocities u */
164 vn = gstat[ga].u[d] + accel[ga][d]*dt + vnrel;
166 xprime[n][d] = x[n][d]+vn*dt;
171 xprime[n][d] = x[n][d];
178 /* Classic version of update, used with berendsen coupling */
179 for(n=start; n<nrend; n++)
181 w_dt = invmass[n]*dt;
194 lg = tcstat[gt].lambda;
199 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
201 vv = lg*vn + f[n][d]*w_dt;
203 /* do not scale the mean velocities u */
205 va = vv + accel[ga][d]*dt;
206 vb = va + (1.0-lg)*u;
208 xprime[n][d] = x[n][d]+vb*dt;
213 xprime[n][d] = x[n][d];
220 static void do_update_vv_vel(int start,int nrend,double dt,
221 t_grp_tcstat *tcstat,t_grp_acc *gstat,
222 rvec accel[],ivec nFreeze[],real invmass[],
223 unsigned short ptype[],
224 unsigned short cFREEZE[],unsigned short cACC[],
226 gmx_bool bExtended, real veta, real alpha)
231 real u,vn,vv,va,vb,vnrel;
237 g = 0.25*dt*veta*alpha;
239 mv2 = series_sinhx(g);
246 for(n=start; n<nrend; n++)
248 w_dt = invmass[n]*dt;
260 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
262 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
270 } /* do_update_vv_vel */
272 static void do_update_vv_pos(int start,int nrend,double dt,
273 t_grp_tcstat *tcstat,t_grp_acc *gstat,
274 rvec accel[],ivec nFreeze[],real invmass[],
275 unsigned short ptype[],
276 unsigned short cFREEZE[],
277 rvec x[],rvec xprime[],rvec v[],
278 rvec f[],gmx_bool bExtended, real veta, real alpha)
288 mr2 = series_sinhx(g);
296 for(n=start; n<nrend; n++) {
304 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
306 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
310 xprime[n][d] = x[n][d];
314 }/* do_update_vv_pos */
316 static void do_update_visc(int start,int nrend,double dt,
317 t_grp_tcstat *tcstat,real invmass[],double nh_vxi[],
318 unsigned short ptype[],unsigned short cTC[],
319 rvec x[],rvec xprime[],rvec v[],
320 rvec f[],matrix M,matrix box,real
322 gmx_bool bNH,gmx_bool bPR)
332 fac = 2*M_PI/(box[ZZ][ZZ]);
335 /* Update with coupling to extended ensembles, used for
336 * Nose-Hoover and Parrinello-Rahman coupling
338 for(n=start; n<nrend; n++) {
344 lg = tcstat[gt].lambda;
345 cosz = cos(fac*x[n][ZZ]);
347 copy_rvec(v[n],vrel);
359 if((ptype[n] != eptVSite) && (ptype[n] != eptShell))
361 vn = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
362 - iprod(M[d],vrel)))/(1 + 0.5*vxi*dt);
365 vn += vc + dt*cosz*cos_accel;
368 xprime[n][d] = x[n][d]+vn*dt;
372 xprime[n][d] = x[n][d];
379 /* Classic version of update, used with berendsen coupling */
380 for(n=start; n<nrend; n++)
382 w_dt = invmass[n]*dt;
387 lg = tcstat[gt].lambda;
388 cosz = cos(fac*x[n][ZZ]);
394 if((ptype[n] != eptVSite) && (ptype[n] != eptShell))
399 /* Do not scale the cosine velocity profile */
400 vv = vc + lg*(vn - vc + f[n][d]*w_dt);
401 /* Add the cosine accelaration profile */
402 vv += dt*cosz*cos_accel;
406 vv = lg*(vn + f[n][d]*w_dt);
409 xprime[n][d] = x[n][d]+vv*dt;
414 xprime[n][d] = x[n][d];
421 static gmx_stochd_t *init_stochd(FILE *fplog,t_inputrec *ir)
430 /* Initiate random number generator for langevin type dynamics,
431 * for BD, SD or velocity rescaling temperature coupling.
433 sd->gaussrand = gmx_rng_init(ir->ld_seed);
435 ngtc = ir->opts.ngtc;
439 snew(sd->bd_rf,ngtc);
441 else if (EI_SD(ir->eI))
444 snew(sd->sdsig,ngtc);
447 for(n=0; n<ngtc; n++)
449 if (ir->opts.tau_t[n] > 0)
451 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
452 sdc[n].eph = exp(sdc[n].gdt/2);
453 sdc[n].emh = exp(-sdc[n].gdt/2);
454 sdc[n].em = exp(-sdc[n].gdt);
458 /* No friction and noise on this group */
464 if (sdc[n].gdt >= 0.05)
466 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
467 - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
468 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
469 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
474 /* Seventh order expansions for small y */
475 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
476 sdc[n].c = y*y*y*(2/3.0+y*(-1/2.0+y*(7/30.0+y*(-1/12.0+y*31/1260.0))));
477 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
480 fprintf(debug,"SD const tc-grp %d: b %g c %g d %g\n",
481 n,sdc[n].b,sdc[n].c,sdc[n].d);
488 void get_stochd_state(gmx_update_t upd,t_state *state)
490 gmx_rng_get_state(upd->sd->gaussrand,state->ld_rng,state->ld_rngi);
493 void set_stochd_state(gmx_update_t upd,t_state *state)
495 gmx_rng_set_state(upd->sd->gaussrand,state->ld_rng,state->ld_rngi[0]);
498 gmx_update_t init_update(FILE *fplog,t_inputrec *ir)
504 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE)
506 upd->sd = init_stochd(fplog,ir);
515 static void do_update_sd1(gmx_stochd_t *sd,
516 int start,int homenr,double dt,
517 rvec accel[],ivec nFreeze[],
518 real invmass[],unsigned short ptype[],
519 unsigned short cFREEZE[],unsigned short cACC[],
520 unsigned short cTC[],
521 rvec x[],rvec xprime[],rvec v[],rvec f[],
523 int ngtc,real tau_t[],real ref_t[])
535 if (homenr > sd->sd_V_nalloc)
537 sd->sd_V_nalloc = over_alloc_dd(homenr);
538 srenew(sd->sd_V,sd->sd_V_nalloc);
540 gaussrand = sd->gaussrand;
542 for(n=0; n<ngtc; n++)
545 /* The mass is encounted for later, since this differs per atom */
546 sig[n].V = sqrt(kT*(1 - sdc[n].em*sdc[n].em));
549 for(n=start; n<start+homenr; n++)
551 ism = sqrt(invmass[n]);
567 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
569 sd_V = ism*sig[gt].V*gmx_rng_gaussian_table(gaussrand);
571 v[n][d] = v[n][d]*sdc[gt].em
572 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
575 xprime[n][d] = x[n][d] + v[n][d]*dt;
580 xprime[n][d] = x[n][d];
586 static void do_update_sd2(gmx_stochd_t *sd,gmx_bool bInitStep,
587 int start,int homenr,
588 rvec accel[],ivec nFreeze[],
589 real invmass[],unsigned short ptype[],
590 unsigned short cFREEZE[],unsigned short cACC[],
591 unsigned short cTC[],
592 rvec x[],rvec xprime[],rvec v[],rvec f[],
594 int ngtc,real tau_t[],real ref_t[],
599 /* The random part of the velocity update, generated in the first
600 * half of the update, needs to be remembered for the second half.
612 if (homenr > sd->sd_V_nalloc)
614 sd->sd_V_nalloc = over_alloc_dd(homenr);
615 srenew(sd->sd_V,sd->sd_V_nalloc);
618 gaussrand = sd->gaussrand;
622 for (n=0; n<ngtc; n++)
625 /* The mass is encounted for later, since this differs per atom */
626 sig[n].V = sqrt(kT*(1-sdc[n].em));
627 sig[n].X = sqrt(kT*sqr(tau_t[n])*sdc[n].c);
628 sig[n].Yv = sqrt(kT*sdc[n].b/sdc[n].c);
629 sig[n].Yx = sqrt(kT*sqr(tau_t[n])*sdc[n].b/(1-sdc[n].em));
633 for (n=start; n<start+homenr; n++)
635 ism = sqrt(invmass[n]);
655 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
661 sd_X[n][d] = ism*sig[gt].X*gmx_rng_gaussian_table(gaussrand);
663 Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
664 + ism*sig[gt].Yv*gmx_rng_gaussian_table(gaussrand);
665 sd_V[n-start][d] = ism*sig[gt].V*gmx_rng_gaussian_table(gaussrand);
667 v[n][d] = vn*sdc[gt].em
668 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
669 + sd_V[n-start][d] - sdc[gt].em*Vmh;
671 xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
676 /* Correct the velocities for the constraints.
677 * This operation introduces some inaccuracy,
678 * since the velocity is determined from differences in coordinates.
681 (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
683 Xmh = sd_V[n-start][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
684 + ism*sig[gt].Yx*gmx_rng_gaussian_table(gaussrand);
685 sd_X[n][d] = ism*sig[gt].X*gmx_rng_gaussian_table(gaussrand);
687 xprime[n][d] += sd_X[n][d] - Xmh;
696 xprime[n][d] = x[n][d];
703 static void do_update_bd(int start,int nrend,double dt,
705 real invmass[],unsigned short ptype[],
706 unsigned short cFREEZE[],unsigned short cTC[],
707 rvec x[],rvec xprime[],rvec v[],
708 rvec f[],real friction_coefficient,
709 int ngtc,real tau_t[],real ref_t[],
710 real *rf,gmx_rng_t gaussrand)
712 /* note -- these appear to be full step velocities . . . */
718 if (friction_coefficient != 0)
720 invfr = 1.0/friction_coefficient;
721 for(n=0; n<ngtc; n++)
723 rf[n] = sqrt(2.0*BOLTZ*ref_t[n]/(friction_coefficient*dt));
728 for(n=0; n<ngtc; n++)
730 rf[n] = sqrt(2.0*BOLTZ*ref_t[n]);
733 for(n=start; (n<nrend); n++)
743 for(d=0; (d<DIM); d++)
745 if((ptype[n]!=eptVSite) && (ptype[n]!=eptShell) && !nFreeze[gf][d])
747 if (friction_coefficient != 0) {
748 vn = invfr*f[n][d] + rf[gt]*gmx_rng_gaussian_table(gaussrand);
752 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
753 vn = 0.5*invmass[n]*f[n][d]*dt
754 + sqrt(0.5*invmass[n])*rf[gt]*gmx_rng_gaussian_table(gaussrand);
758 xprime[n][d] = x[n][d]+vn*dt;
763 xprime[n][d] = x[n][d];
769 static void dump_it_all(FILE *fp,const char *title,
770 int natoms,rvec x[],rvec xp[],rvec v[],rvec f[])
775 fprintf(fp,"%s\n",title);
776 pr_rvecs(fp,0,"x",x,natoms);
777 pr_rvecs(fp,0,"xp",xp,natoms);
778 pr_rvecs(fp,0,"v",v,natoms);
779 pr_rvecs(fp,0,"f",f,natoms);
784 static void calc_ke_part_normal(rvec v[], t_grpopts *opts,t_mdatoms *md,
785 gmx_ekindata_t *ekind,t_nrnb *nrnb,gmx_bool bEkinAveVel,
786 gmx_bool bSaveEkinOld)
788 int start=md->start,homenr=md->homenr;
789 int g,d,n,m,ga=0,gt=0;
792 t_grp_tcstat *tcstat=ekind->tcstat;
793 t_grp_acc *grpstat=ekind->grpstat;
796 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
797 an option, but not supported now. Additionally, if we are doing iterations.
798 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
799 bSavEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't copy over the ekinh_old.
800 If FALSE, we overrwrite it.
803 /* group velocities are calculated in update_ekindata and
804 * accumulated in acumulate_groups.
805 * Now the partial global and groups ekin.
807 for(g=0; (g<opts->ngtc); g++)
811 copy_mat(tcstat[g].ekinh,tcstat[g].ekinh_old);
814 clear_mat(tcstat[g].ekinf);
816 clear_mat(tcstat[g].ekinh);
819 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
822 ekind->dekindl_old = ekind->dekindl;
825 for(n=start; (n<start+homenr); n++)
835 hm = 0.5*md->massT[n];
837 for(d=0; (d<DIM); d++)
839 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
841 for(d=0; (d<DIM); d++)
843 for (m=0;(m<DIM); m++)
845 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
848 tcstat[gt].ekinf[m][d]+=hm*v_corrt[m]*v_corrt[d];
852 tcstat[gt].ekinh[m][d]+=hm*v_corrt[m]*v_corrt[d];
856 if (md->nMassPerturbed && md->bPerturbed[n])
858 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt,v_corrt);
861 ekind->dekindl = dekindl;
862 inc_nrnb(nrnb,eNR_EKIN,homenr);
865 static void calc_ke_part_visc(matrix box,rvec x[],rvec v[],
866 t_grpopts *opts,t_mdatoms *md,
867 gmx_ekindata_t *ekind,
868 t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
870 int start=md->start,homenr=md->homenr;
874 t_grp_tcstat *tcstat=ekind->tcstat;
875 t_cos_acc *cosacc=&(ekind->cosacc);
880 for(g=0; g<opts->ngtc; g++)
882 copy_mat(ekind->tcstat[g].ekinh,ekind->tcstat[g].ekinh_old);
883 clear_mat(ekind->tcstat[g].ekinh);
885 ekind->dekindl_old = ekind->dekindl;
887 fac = 2*M_PI/box[ZZ][ZZ];
890 for(n=start; n<start+homenr; n++)
896 hm = 0.5*md->massT[n];
898 /* Note that the times of x and v differ by half a step */
899 /* MRS -- would have to be changed for VV */
900 cosz = cos(fac*x[n][ZZ]);
901 /* Calculate the amplitude of the new velocity profile */
902 mvcos += 2*cosz*md->massT[n]*v[n][XX];
904 copy_rvec(v[n],v_corrt);
905 /* Subtract the profile for the kinetic energy */
906 v_corrt[XX] -= cosz*cosacc->vcos;
907 for (d=0; (d<DIM); d++)
909 for (m=0; (m<DIM); m++)
911 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
914 tcstat[gt].ekinf[m][d]+=hm*v_corrt[m]*v_corrt[d];
918 tcstat[gt].ekinh[m][d]+=hm*v_corrt[m]*v_corrt[d];
922 if(md->nPerturbed && md->bPerturbed[n])
924 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt,v_corrt);
927 ekind->dekindl = dekindl;
928 cosacc->mvcos = mvcos;
930 inc_nrnb(nrnb,eNR_EKIN,homenr);
933 void calc_ke_part(t_state *state,t_grpopts *opts,t_mdatoms *md,
934 gmx_ekindata_t *ekind,t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
936 if (ekind->cosacc.cos_accel == 0)
938 calc_ke_part_normal(state->v,opts,md,ekind,nrnb,bEkinAveVel,bSaveEkinOld);
942 calc_ke_part_visc(state->box,state->x,state->v,opts,md,ekind,nrnb,bEkinAveVel,bSaveEkinOld);
946 void init_ekinstate(ekinstate_t *ekinstate,const t_inputrec *ir)
948 ekinstate->ekin_n = ir->opts.ngtc;
949 snew(ekinstate->ekinh,ekinstate->ekin_n);
950 snew(ekinstate->ekinf,ekinstate->ekin_n);
951 snew(ekinstate->ekinh_old,ekinstate->ekin_n);
952 snew(ekinstate->ekinscalef_nhc,ekinstate->ekin_n);
953 snew(ekinstate->ekinscaleh_nhc,ekinstate->ekin_n);
954 snew(ekinstate->vscale_nhc,ekinstate->ekin_n);
955 ekinstate->dekindl = 0;
956 ekinstate->mvcos = 0;
959 void update_ekinstate(ekinstate_t *ekinstate,gmx_ekindata_t *ekind)
963 for(i=0;i<ekinstate->ekin_n;i++)
965 copy_mat(ekind->tcstat[i].ekinh,ekinstate->ekinh[i]);
966 copy_mat(ekind->tcstat[i].ekinf,ekinstate->ekinf[i]);
967 copy_mat(ekind->tcstat[i].ekinh_old,ekinstate->ekinh_old[i]);
968 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
969 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
970 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
973 copy_mat(ekind->ekin,ekinstate->ekin_total);
974 ekinstate->dekindl = ekind->dekindl;
975 ekinstate->mvcos = ekind->cosacc.mvcos;
979 void restore_ekinstate_from_state(t_commrec *cr,
980 gmx_ekindata_t *ekind,ekinstate_t *ekinstate)
986 for(i=0;i<ekinstate->ekin_n;i++)
988 copy_mat(ekinstate->ekinh[i],ekind->tcstat[i].ekinh);
989 copy_mat(ekinstate->ekinf[i],ekind->tcstat[i].ekinf);
990 copy_mat(ekinstate->ekinh_old[i],ekind->tcstat[i].ekinh_old);
991 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
992 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
993 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
996 copy_mat(ekinstate->ekin_total,ekind->ekin);
998 ekind->dekindl = ekinstate->dekindl;
999 ekind->cosacc.mvcos = ekinstate->mvcos;
1000 n = ekinstate->ekin_n;
1005 gmx_bcast(sizeof(n),&n,cr);
1008 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1009 ekind->tcstat[i].ekinh[0],cr);
1010 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1011 ekind->tcstat[i].ekinf[0],cr);
1012 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1013 ekind->tcstat[i].ekinh_old[0],cr);
1015 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1016 &(ekind->tcstat[i].ekinscalef_nhc),cr);
1017 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1018 &(ekind->tcstat[i].ekinscaleh_nhc),cr);
1019 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1020 &(ekind->tcstat[i].vscale_nhc),cr);
1022 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1025 gmx_bcast(sizeof(ekind->dekindl),&ekind->dekindl,cr);
1026 gmx_bcast(sizeof(ekind->cosacc.mvcos),&ekind->cosacc.mvcos,cr);
1030 void set_deform_reference_box(gmx_update_t upd,gmx_large_int_t step,matrix box)
1032 upd->deformref_step = step;
1033 copy_mat(box,upd->deformref_box);
1036 static void deform(gmx_update_t upd,
1037 int start,int homenr,rvec x[],matrix box,matrix *scale_tot,
1038 const t_inputrec *ir,gmx_large_int_t step)
1040 matrix bnew,invbox,mu;
1044 elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1046 for(i=0; i<DIM; i++)
1048 for(j=0; j<DIM; j++)
1050 if (ir->deform[i][j] != 0)
1053 upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1057 /* We correct the off-diagonal elements,
1058 * which can grow indefinitely during shearing,
1059 * so the shifts do not get messed up.
1061 for(i=1; i<DIM; i++)
1063 for(j=i-1; j>=0; j--)
1065 while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1067 rvec_dec(bnew[i],bnew[j]);
1069 while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1071 rvec_inc(bnew[i],bnew[j]);
1075 m_inv_ur0(box,invbox);
1077 mmul_ur0(box,invbox,mu);
1079 for(i=start; i<start+homenr; i++)
1081 x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1082 x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1083 x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1087 /* The transposes of the scaling matrices are stored,
1088 * so we need to do matrix multiplication in the inverse order.
1090 mmul_ur0(*scale_tot,mu,*scale_tot);
1094 static void combine_forces(int nstlist,
1095 gmx_constr_t constr,
1096 t_inputrec *ir,t_mdatoms *md,t_idef *idef,
1097 t_commrec *cr,gmx_large_int_t step,t_state *state,
1098 int start,int nrend,
1099 rvec f[],rvec f_lr[],
1104 /* f contains the short-range forces + the long range forces
1105 * which are stored separately in f_lr.
1108 if (constr != NULL && !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1110 /* We need to constrain the LR forces separately,
1111 * because due to the different pre-factor for the SR and LR
1112 * forces in the update algorithm, we can not determine
1113 * the constraint force for the coordinate constraining.
1114 * Constrain only the additional LR part of the force.
1116 /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1117 constrain(NULL,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
1118 state->x,f_lr,f_lr,state->box,state->lambda,NULL,
1119 NULL,NULL,nrnb,econqForce,ir->epc==epcMTTK,state->veta,state->veta);
1122 /* Add nstlist-1 times the LR force to the sum of both forces
1123 * and store the result in forces_lr.
1126 for(i=start; i<nrend; i++)
1128 for(d=0; d<DIM; d++)
1130 f_lr[i][d] = f[i][d] + nm1*f_lr[i][d];
1135 void update_tcouple(FILE *fplog,
1136 gmx_large_int_t step,
1137 t_inputrec *inputrec,
1139 gmx_ekindata_t *ekind,
1140 gmx_wallcycle_t wcycle,
1146 gmx_bool bTCouple=FALSE;
1148 int i,start,end,homenr;
1150 /* if using vv, we do this elsewhere in the code */
1151 if (inputrec->etc != etcNO &&
1152 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec)))
1154 /* We should only couple after a step where energies were determined */
1155 bTCouple = (inputrec->nsttcouple == 1 ||
1156 do_per_step(step+inputrec->nsttcouple-1,
1157 inputrec->nsttcouple));
1162 dttc = inputrec->nsttcouple*inputrec->delta_t;
1164 switch (inputrec->etc)
1169 berendsen_tcoupl(inputrec,ekind,dttc);
1172 nosehoover_tcoupl(&(inputrec->opts),ekind,dttc,
1173 state->nosehoover_xi,state->nosehoover_vxi,MassQ);
1176 vrescale_tcoupl(inputrec,ekind,dttc,
1177 state->therm_integral,upd->sd->gaussrand);
1180 /* rescale in place here */
1181 if (EI_VV(inputrec->eI))
1183 rescale_velocities(ekind,md,md->start,md->start+md->homenr,state->v);
1188 /* Set the T scaling lambda to 1 to have no scaling */
1189 for(i=0; (i<inputrec->opts.ngtc); i++)
1191 ekind->tcstat[i].lambda = 1.0;
1196 void update_pcouple(FILE *fplog,
1197 gmx_large_int_t step,
1198 t_inputrec *inputrec,
1202 gmx_wallcycle_t wcycle,
1206 gmx_bool bPCouple=FALSE;
1210 /* if using vv, we do this elsewhere in the code */
1211 if (inputrec->epc != epcNO &&
1212 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec)))
1214 /* We should only couple after a step where energies were determined */
1215 bPCouple = (inputrec->nstpcouple == 1 ||
1216 do_per_step(step+inputrec->nstpcouple-1,
1217 inputrec->nstpcouple));
1220 clear_mat(pcoupl_mu);
1221 for(i=0; i<DIM; i++)
1223 pcoupl_mu[i][i] = 1.0;
1230 dtpc = inputrec->nstpcouple*inputrec->delta_t;
1232 switch (inputrec->epc)
1234 /* We can always pcoupl, even if we did not sum the energies
1235 * the previous step, since state->pres_prev is only updated
1236 * when the energies have been summed.
1240 case (epcBERENDSEN):
1243 berendsen_pcoupl(fplog,step,inputrec,dtpc,state->pres_prev,state->box,
1247 case (epcPARRINELLORAHMAN):
1248 parrinellorahman_pcoupl(fplog,step,inputrec,dtpc,state->pres_prev,
1249 state->box,state->box_rel,state->boxv,
1250 M,pcoupl_mu,bInitStep);
1258 static rvec *get_xprime(const t_state *state,gmx_update_t upd)
1260 if (state->nalloc > upd->xp_nalloc)
1262 upd->xp_nalloc = state->nalloc;
1263 srenew(upd->xp,upd->xp_nalloc);
1269 void update_constraints(FILE *fplog,
1270 gmx_large_int_t step,
1271 real *dvdlambda, /* FEP stuff */
1272 t_inputrec *inputrec, /* input record and box stuff */
1273 gmx_ekindata_t *ekind,
1277 rvec force[], /* forces on home particles */
1280 tensor vir, /* tensors for virial and ekin, needed for computing */
1283 gmx_wallcycle_t wcycle,
1285 gmx_constr_t constr,
1287 gmx_bool bFirstHalf,
1291 gmx_bool bExtended,bTrotter,bLastStep,bLog=FALSE,bEner=FALSE,bDoConstr=FALSE;
1294 int start,homenr,nrend,i,n,m,g,d;
1296 rvec *vbuf,*xprime=NULL;
1298 if (constr) {bDoConstr=TRUE;}
1299 if (bFirstHalf && !EI_VV(inputrec->eI)) {bDoConstr=FALSE;}
1301 /* for now, SD update is here -- though it really seems like it
1302 should be reformulated as a velocity verlet method, since it has two parts */
1305 homenr = md->homenr;
1306 nrend = start+homenr;
1308 dt = inputrec->delta_t;
1313 * APPLY CONSTRAINTS:
1316 * When doing PR pressure coupling we have to constrain the
1317 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1318 * it is enough to do this once though, since the relative velocities
1319 * after this will be normal to the bond vector
1324 /* clear out constraints before applying */
1325 clear_mat(vir_part);
1327 xprime = get_xprime(state,upd);
1329 bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1330 bLog = (do_per_step(step,inputrec->nstlog) || bLastStep || (step < 0));
1331 bEner = (do_per_step(step,inputrec->nstenergy) || bLastStep);
1332 /* Constrain the coordinates xprime */
1333 wallcycle_start(wcycle,ewcCONSTR);
1334 if (EI_VV(inputrec->eI) && bFirstHalf)
1336 constrain(NULL,bLog,bEner,constr,idef,
1337 inputrec,ekind,cr,step,1,md,
1338 state->x,state->v,state->v,
1339 state->box,state->lambda,dvdlambda,
1340 NULL,bCalcVir ? &vir_con : NULL,nrnb,econqVeloc,
1341 inputrec->epc==epcMTTK,state->veta,vetanew);
1345 constrain(NULL,bLog,bEner,constr,idef,
1346 inputrec,ekind,cr,step,1,md,
1347 state->x,xprime,NULL,
1348 state->box,state->lambda,dvdlambda,
1349 state->v,bCalcVir ? &vir_con : NULL ,nrnb,econqCoord,
1350 inputrec->epc==epcMTTK,state->veta,state->veta);
1352 wallcycle_stop(wcycle,ewcCONSTR);
1356 dump_it_all(fplog,"After Shake",
1357 state->natoms,state->x,xprime,state->v,force);
1361 if (inputrec->eI == eiSD2)
1363 /* A correction factor eph is needed for the SD constraint force */
1364 /* Here we can, unfortunately, not have proper corrections
1365 * for different friction constants, so we use the first one.
1367 for(i=0; i<DIM; i++)
1369 for(m=0; m<DIM; m++)
1371 vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1377 m_add(vir_part,vir_con,vir_part);
1381 pr_rvecs(debug,0,"constraint virial",vir_part,DIM);
1387 if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1389 xprime = get_xprime(state,upd);
1391 /* The second part of the SD integration */
1392 do_update_sd2(upd->sd,FALSE,start,homenr,
1393 inputrec->opts.acc,inputrec->opts.nFreeze,
1394 md->invmass,md->ptype,
1395 md->cFREEZE,md->cACC,md->cTC,
1396 state->x,xprime,state->v,force,state->sd_X,
1397 inputrec->opts.ngtc,inputrec->opts.tau_t,
1398 inputrec->opts.ref_t,FALSE);
1399 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1403 /* Constrain the coordinates xprime */
1404 wallcycle_start(wcycle,ewcCONSTR);
1405 constrain(NULL,bLog,bEner,constr,idef,
1406 inputrec,NULL,cr,step,1,md,
1407 state->x,xprime,NULL,
1408 state->box,state->lambda,dvdlambda,
1409 NULL,NULL,nrnb,econqCoord,FALSE,0,0);
1410 wallcycle_stop(wcycle,ewcCONSTR);
1414 /* We must always unshift after updating coordinates; if we did not shake
1415 x was shifted in do_force */
1417 if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1419 if (graph && (graph->nnodes > 0))
1421 unshift_x(graph,state->box,state->x,upd->xp);
1422 if (TRICLINIC(state->box))
1424 inc_nrnb(nrnb,eNR_SHIFTX,2*graph->nnodes);
1428 inc_nrnb(nrnb,eNR_SHIFTX,graph->nnodes);
1433 copy_rvecn(upd->xp,state->x,start,nrend);
1436 dump_it_all(fplog,"After unshift",
1437 state->natoms,state->x,upd->xp,state->v,force);
1439 /* ############# END the update of velocities and positions ######### */
1442 void update_box(FILE *fplog,
1443 gmx_large_int_t step,
1444 t_inputrec *inputrec, /* input record and box stuff */
1448 rvec force[], /* forces on home particles */
1452 gmx_wallcycle_t wcycle,
1455 gmx_bool bFirstHalf)
1457 gmx_bool bExtended,bTrotter,bLastStep,bLog=FALSE,bEner=FALSE;
1460 int start,homenr,nrend,i,n,m,g;
1464 homenr = md->homenr;
1465 nrend = start+homenr;
1468 (inputrec->etc == etcNOSEHOOVER) ||
1469 (inputrec->epc == epcPARRINELLORAHMAN) ||
1470 (inputrec->epc == epcMTTK);
1472 dt = inputrec->delta_t;
1476 /* now update boxes */
1477 switch (inputrec->epc) {
1480 case (epcBERENDSEN):
1481 berendsen_pscale(inputrec,pcoupl_mu,state->box,state->box_rel,
1482 start,homenr,state->x,md->cFREEZE,nrnb);
1484 case (epcPARRINELLORAHMAN):
1485 /* The box velocities were updated in do_pr_pcoupl in the update
1486 * iteration, but we dont change the box vectors until we get here
1487 * since we need to be able to shift/unshift above.
1493 state->box[i][m] += dt*state->boxv[i][m];
1496 preserve_box_shape(inputrec,state->box_rel,state->box);
1498 /* Scale the coordinates */
1499 for(n=start; (n<start+homenr); n++)
1501 tmvmul_ur0(pcoupl_mu,state->x[n],state->x[n]);
1505 switch (inputrec->epct)
1507 case (epctISOTROPIC):
1508 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1509 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1510 Side length scales as exp(veta*dt) */
1512 msmul(state->box,exp(state->veta*dt),state->box);
1514 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1515 o If we assume isotropic scaling, and box length scaling
1516 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1517 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1518 determinant of B is L^DIM det(M), and the determinant
1519 of dB/dt is (dL/dT)^DIM det (M). veta will be
1520 (det(dB/dT)/det(B))^(1/3). Then since M =
1521 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1523 msmul(state->box,state->veta,state->boxv);
1533 if ((!(IR_NPT_TROTTER(inputrec))) && scale_tot)
1535 /* The transposes of the scaling matrices are stored,
1536 * therefore we need to reverse the order in the multiplication.
1538 mmul_ur0(*scale_tot,pcoupl_mu,*scale_tot);
1541 if (DEFORM(*inputrec))
1543 deform(upd,start,homenr,state->x,state->box,scale_tot,inputrec,step);
1546 dump_it_all(fplog,"After update",
1547 state->natoms,state->x,upd->xp,state->v,force);
1550 void update_coords(FILE *fplog,
1551 gmx_large_int_t step,
1552 t_inputrec *inputrec, /* input record and box stuff */
1555 rvec *f, /* forces on home particles */
1559 gmx_ekindata_t *ekind,
1561 gmx_wallcycle_t wcycle,
1565 t_commrec *cr, /* these shouldn't be here -- need to think about it */
1567 gmx_constr_t constr,
1570 gmx_bool bExtended,bNH,bPR,bTrotter,bLastStep,bLog=FALSE,bEner=FALSE;
1572 real *imass,*imassin;
1575 int start,homenr,nrend,i,j,d,n,m,g;
1576 int blen0,blen1,iatom,jatom,nshake,nsettle,nconstr,nexpand;
1579 rvec *vcom,*xcom,*vall,*xall,*xin,*vin,*forcein,*fall,*xpall,*xprimein,*xprime;
1582 /* Running the velocity half does nothing except for velocity verlet */
1583 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1584 !EI_VV(inputrec->eI))
1586 gmx_incons("update_coords called for velocity without VV integrator");
1590 homenr = md->homenr;
1591 nrend = start+homenr;
1593 xprime = get_xprime(state,upd);
1595 dt = inputrec->delta_t;
1598 /* We need to update the NMR restraint history when time averaging is used */
1599 if (state->flags & (1<<estDISRE_RM3TAV))
1601 update_disres_history(fcd,&state->hist);
1603 if (state->flags & (1<<estORIRE_DTAV))
1605 update_orires_history(fcd,&state->hist);
1608 bNH = inputrec->etc == etcNOSEHOOVER;
1609 bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1611 bExtended = bNH || bPR;
1613 if (bDoLR && inputrec->nstlist > 1 && !EI_VV(inputrec->eI)) /* get this working with VV? */
1615 /* Store the total force + nstlist-1 times the LR force
1616 * in forces_lr, so it can be used in a normal update algorithm
1617 * to produce twin time stepping.
1619 /* is this correct in the new construction? MRS */
1620 combine_forces(inputrec->nstlist,constr,inputrec,md,idef,cr,step,state,
1621 start,nrend,f,f_lr,nrnb);
1629 /* ############# START The update of velocities and positions ######### */
1631 dump_it_all(fplog,"Before update",
1632 state->natoms,state->x,xprime,state->v,force);
1634 switch (inputrec->eI) {
1636 if (ekind->cosacc.cos_accel == 0) {
1637 /* use normal version of update */
1638 do_update_md(start,nrend,dt,
1639 ekind->tcstat,ekind->grpstat,state->nosehoover_vxi,
1640 inputrec->opts.acc,inputrec->opts.nFreeze,md->invmass,md->ptype,
1641 md->cFREEZE,md->cACC,md->cTC,
1642 state->x,xprime,state->v,force,M,
1647 do_update_visc(start,nrend,dt,
1648 ekind->tcstat,md->invmass,state->nosehoover_vxi,
1649 md->ptype,md->cTC,state->x,xprime,state->v,force,M,
1650 state->box,ekind->cosacc.cos_accel,ekind->cosacc.vcos,bNH,bPR);
1654 do_update_sd1(upd->sd,start,homenr,dt,
1655 inputrec->opts.acc,inputrec->opts.nFreeze,
1656 md->invmass,md->ptype,
1657 md->cFREEZE,md->cACC,md->cTC,
1658 state->x,xprime,state->v,force,state->sd_X,
1659 inputrec->opts.ngtc,inputrec->opts.tau_t,inputrec->opts.ref_t);
1662 /* The SD update is done in 2 parts, because an extra constraint step
1665 do_update_sd2(upd->sd,bInitStep,start,homenr,
1666 inputrec->opts.acc,inputrec->opts.nFreeze,
1667 md->invmass,md->ptype,
1668 md->cFREEZE,md->cACC,md->cTC,
1669 state->x,xprime,state->v,force,state->sd_X,
1670 inputrec->opts.ngtc,inputrec->opts.tau_t,inputrec->opts.ref_t,
1674 do_update_bd(start,nrend,dt,
1675 inputrec->opts.nFreeze,md->invmass,md->ptype,
1676 md->cFREEZE,md->cTC,
1677 state->x,xprime,state->v,force,
1679 inputrec->opts.ngtc,inputrec->opts.tau_t,inputrec->opts.ref_t,
1680 upd->sd->bd_rf,upd->sd->gaussrand);
1684 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
1685 switch (UpdatePart) {
1688 do_update_vv_vel(start,nrend,dt,
1689 ekind->tcstat,ekind->grpstat,
1690 inputrec->opts.acc,inputrec->opts.nFreeze,
1691 md->invmass,md->ptype,
1692 md->cFREEZE,md->cACC,
1694 bExtended,state->veta,alpha);
1697 do_update_vv_pos(start,nrend,dt,
1698 ekind->tcstat,ekind->grpstat,
1699 inputrec->opts.acc,inputrec->opts.nFreeze,
1700 md->invmass,md->ptype,md->cFREEZE,
1701 state->x,xprime,state->v,force,
1702 bExtended,state->veta,alpha);
1707 gmx_fatal(FARGS,"Don't know how to update coordinates");
1713 void correct_ekin(FILE *log,int start,int end,rvec v[],rvec vcm,real mass[],
1714 real tmass,tensor ekin)
1717 * This is a debugging routine. It should not be called for production code
1719 * The kinetic energy should calculated according to:
1720 * Ekin = 1/2 m (v-vcm)^2
1721 * However the correction is not always applied, since vcm may not be
1722 * known in time and we compute
1723 * Ekin' = 1/2 m v^2 instead
1724 * This can be corrected afterwards by computing
1725 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
1727 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
1734 /* Local particles */
1737 /* Processor dependent part. */
1739 for(i=start; (i<end); i++)
1743 for(j=0; (j<DIM); j++)
1749 svmul(1/tmass,vcm,vcm);
1750 svmul(0.5,vcm,hvcm);
1752 for(j=0; (j<DIM); j++)
1754 for(k=0; (k<DIM); k++)
1756 dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
1759 pr_rvecs(log,0,"dekin",dekin,DIM);
1760 pr_rvecs(log,0," ekin", ekin,DIM);
1761 fprintf(log,"dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
1762 trace(dekin),trace(ekin),vcm[XX],vcm[YY],vcm[ZZ]);
1763 fprintf(log,"mv = (%8.4f %8.4f %8.4f)\n",
1764 mv[XX],mv[YY],mv[ZZ]);