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 void store_rvec(rvec *from, rvec *to, int n) {
116 copy_rvec(from[i],to[i]);
120 static void do_update_md(int start,int nrend,double dt,
121 t_grp_tcstat *tcstat,t_grp_acc *gstat,double nh_vxi[],
122 rvec accel[],ivec nFreeze[],real invmass[],
123 unsigned short ptype[],unsigned short cFREEZE[],
124 unsigned short cACC[],unsigned short cTC[],
125 rvec x[],rvec xprime[],rvec v[],
127 gmx_bool bNH,gmx_bool bPR)
132 real vn,vv,va,vb,vnrel;
138 /* Update with coupling to extended ensembles, used for
139 * Nose-Hoover and Parrinello-Rahman coupling
140 * Nose-Hoover uses the reversible leap-frog integrator from
141 * Holian et al. Phys Rev E 52(3) : 2338, 1995
143 for(n=start; n<nrend; n++)
158 lg = tcstat[gt].lambda;
162 rvec_sub(v[n],gstat[ga].u,vrel);
166 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
168 vnrel = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
169 - iprod(M[d],vrel)))/(1 + 0.5*vxi*dt);
170 /* do not scale the mean velocities u */
171 vn = gstat[ga].u[d] + accel[ga][d]*dt + vnrel;
173 xprime[n][d] = x[n][d]+vn*dt;
178 xprime[n][d] = x[n][d];
185 /* Classic version of update, used with berendsen coupling */
186 for(n=start; n<nrend; n++)
188 w_dt = invmass[n]*dt;
201 lg = tcstat[gt].lambda;
206 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
208 vv = lg*vn + f[n][d]*w_dt;
210 /* do not scale the mean velocities u */
212 va = vv + accel[ga][d]*dt;
213 vb = va + (1.0-lg)*u;
215 xprime[n][d] = x[n][d]+vb*dt;
220 xprime[n][d] = x[n][d];
227 static void do_update_vv_vel(int start,int nrend,double dt,
228 t_grp_tcstat *tcstat,t_grp_acc *gstat,
229 rvec accel[],ivec nFreeze[],real invmass[],
230 unsigned short ptype[],
231 unsigned short cFREEZE[],unsigned short cACC[],
233 gmx_bool bExtended, real veta, real alpha)
238 real u,vn,vv,va,vb,vnrel;
244 g = 0.25*dt*veta*alpha;
246 mv2 = series_sinhx(g);
253 for(n=start; n<nrend; n++)
255 w_dt = invmass[n]*dt;
267 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
269 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
277 } /* do_update_vv_vel */
279 static void do_update_vv_pos(int start,int nrend,double dt,
280 t_grp_tcstat *tcstat,t_grp_acc *gstat,
281 rvec accel[],ivec nFreeze[],real invmass[],
282 unsigned short ptype[],
283 unsigned short cFREEZE[],
284 rvec x[],rvec xprime[],rvec v[],
285 rvec f[],gmx_bool bExtended, real veta, real alpha)
295 mr2 = series_sinhx(g);
303 for(n=start; n<nrend; n++) {
311 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
313 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
317 xprime[n][d] = x[n][d];
321 }/* do_update_vv_pos */
323 static void do_update_visc(int start,int nrend,double dt,
324 t_grp_tcstat *tcstat,real invmass[],double nh_vxi[],
325 unsigned short ptype[],unsigned short cTC[],
326 rvec x[],rvec xprime[],rvec v[],
327 rvec f[],matrix M,matrix box,real
329 gmx_bool bNH,gmx_bool bPR)
339 fac = 2*M_PI/(box[ZZ][ZZ]);
342 /* Update with coupling to extended ensembles, used for
343 * Nose-Hoover and Parrinello-Rahman coupling
345 for(n=start; n<nrend; n++) {
351 lg = tcstat[gt].lambda;
352 cosz = cos(fac*x[n][ZZ]);
354 copy_rvec(v[n],vrel);
366 if((ptype[n] != eptVSite) && (ptype[n] != eptShell))
368 vn = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
369 - iprod(M[d],vrel)))/(1 + 0.5*vxi*dt);
372 vn += vc + dt*cosz*cos_accel;
375 xprime[n][d] = x[n][d]+vn*dt;
379 xprime[n][d] = x[n][d];
386 /* Classic version of update, used with berendsen coupling */
387 for(n=start; n<nrend; n++)
389 w_dt = invmass[n]*dt;
394 lg = tcstat[gt].lambda;
395 cosz = cos(fac*x[n][ZZ]);
401 if((ptype[n] != eptVSite) && (ptype[n] != eptShell))
406 /* Do not scale the cosine velocity profile */
407 vv = vc + lg*(vn - vc + f[n][d]*w_dt);
408 /* Add the cosine accelaration profile */
409 vv += dt*cosz*cos_accel;
413 vv = lg*(vn + f[n][d]*w_dt);
416 xprime[n][d] = x[n][d]+vv*dt;
421 xprime[n][d] = x[n][d];
428 static gmx_stochd_t *init_stochd(FILE *fplog,t_inputrec *ir)
437 /* Initiate random number generator for langevin type dynamics,
438 * for BD, SD or velocity rescaling temperature coupling.
440 sd->gaussrand = gmx_rng_init(ir->ld_seed);
442 ngtc = ir->opts.ngtc;
444 if (ir->eI == eiBD) {
445 snew(sd->bd_rf,ngtc);
446 } else if (EI_SD(ir->eI)) {
448 snew(sd->sdsig,ngtc);
451 for(n=0; n<ngtc; n++) {
452 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
453 sdc[n].eph = exp(sdc[n].gdt/2);
454 sdc[n].emh = exp(-sdc[n].gdt/2);
455 sdc[n].em = exp(-sdc[n].gdt);
456 if (sdc[n].gdt >= 0.05) {
457 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
458 - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
459 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
460 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
463 /* Seventh order expansions for small y */
464 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
465 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))));
466 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
469 fprintf(debug,"SD const tc-grp %d: b %g c %g d %g\n",
470 n,sdc[n].b,sdc[n].c,sdc[n].d);
477 void get_stochd_state(gmx_update_t upd,t_state *state)
479 gmx_rng_get_state(upd->sd->gaussrand,state->ld_rng,state->ld_rngi);
482 void set_stochd_state(gmx_update_t upd,t_state *state)
484 gmx_rng_set_state(upd->sd->gaussrand,state->ld_rng,state->ld_rngi[0]);
487 gmx_update_t init_update(FILE *fplog,t_inputrec *ir)
493 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE)
495 upd->sd = init_stochd(fplog,ir);
504 static void do_update_sd1(gmx_stochd_t *sd,
505 int start,int homenr,double dt,
506 rvec accel[],ivec nFreeze[],
507 real invmass[],unsigned short ptype[],
508 unsigned short cFREEZE[],unsigned short cACC[],
509 unsigned short cTC[],
510 rvec x[],rvec xprime[],rvec v[],rvec f[],
512 int ngtc,real tau_t[],real ref_t[])
524 if (homenr > sd->sd_V_nalloc)
526 sd->sd_V_nalloc = over_alloc_dd(homenr);
527 srenew(sd->sd_V,sd->sd_V_nalloc);
529 gaussrand = sd->gaussrand;
531 for(n=0; n<ngtc; n++)
534 /* The mass is encounted for later, since this differs per atom */
535 sig[n].V = sqrt(2*kT*(1 - sdc[n].em));
538 for(n=start; n<start+homenr; n++)
540 ism = sqrt(invmass[n]);
556 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
558 sd_V = ism*sig[gt].V*gmx_rng_gaussian_table(gaussrand);
560 v[n][d] = v[n][d]*sdc[gt].em
561 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
564 xprime[n][d] = x[n][d] + v[n][d]*dt;
569 xprime[n][d] = x[n][d];
575 static void do_update_sd2(gmx_stochd_t *sd,gmx_bool bInitStep,
576 int start,int homenr,
577 rvec accel[],ivec nFreeze[],
578 real invmass[],unsigned short ptype[],
579 unsigned short cFREEZE[],unsigned short cACC[],
580 unsigned short cTC[],
581 rvec x[],rvec xprime[],rvec v[],rvec f[],
583 int ngtc,real tau_t[],real ref_t[],
588 /* The random part of the velocity update, generated in the first
589 * half of the update, needs to be remembered for the second half.
601 if (homenr > sd->sd_V_nalloc)
603 sd->sd_V_nalloc = over_alloc_dd(homenr);
604 srenew(sd->sd_V,sd->sd_V_nalloc);
607 gaussrand = sd->gaussrand;
611 for (n=0; n<ngtc; n++)
614 /* The mass is encounted for later, since this differs per atom */
615 sig[n].V = sqrt(kT*(1-sdc[n].em));
616 sig[n].X = sqrt(kT*sqr(tau_t[n])*sdc[n].c);
617 sig[n].Yv = sqrt(kT*sdc[n].b/sdc[n].c);
618 sig[n].Yx = sqrt(kT*sqr(tau_t[n])*sdc[n].b/(1-sdc[n].em));
622 for (n=start; n<start+homenr; n++)
624 ism = sqrt(invmass[n]);
644 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
650 sd_X[n][d] = ism*sig[gt].X*gmx_rng_gaussian_table(gaussrand);
652 Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
653 + ism*sig[gt].Yv*gmx_rng_gaussian_table(gaussrand);
654 sd_V[n-start][d] = ism*sig[gt].V*gmx_rng_gaussian_table(gaussrand);
656 v[n][d] = vn*sdc[gt].em
657 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
658 + sd_V[n-start][d] - sdc[gt].em*Vmh;
660 xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
665 /* Correct the velocities for the constraints.
666 * This operation introduces some inaccuracy,
667 * since the velocity is determined from differences in coordinates.
670 (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
672 Xmh = sd_V[n-start][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
673 + ism*sig[gt].Yx*gmx_rng_gaussian_table(gaussrand);
674 sd_X[n][d] = ism*sig[gt].X*gmx_rng_gaussian_table(gaussrand);
676 xprime[n][d] += sd_X[n][d] - Xmh;
685 xprime[n][d] = x[n][d];
692 static void do_update_bd(int start,int nrend,double dt,
694 real invmass[],unsigned short ptype[],
695 unsigned short cFREEZE[],unsigned short cTC[],
696 rvec x[],rvec xprime[],rvec v[],
697 rvec f[],real friction_coefficient,
698 int ngtc,real tau_t[],real ref_t[],
699 real *rf,gmx_rng_t gaussrand)
701 /* note -- these appear to be full step velocities . . . */
707 if (friction_coefficient != 0)
709 invfr = 1.0/friction_coefficient;
710 for(n=0; n<ngtc; n++)
712 rf[n] = sqrt(2.0*BOLTZ*ref_t[n]/(friction_coefficient*dt));
717 for(n=0; n<ngtc; n++)
719 rf[n] = sqrt(2.0*BOLTZ*ref_t[n]);
722 for(n=start; (n<nrend); n++)
732 for(d=0; (d<DIM); d++)
734 if((ptype[n]!=eptVSite) && (ptype[n]!=eptShell) && !nFreeze[gf][d])
736 if (friction_coefficient != 0) {
737 vn = invfr*f[n][d] + rf[gt]*gmx_rng_gaussian_table(gaussrand);
741 /* NOTE: invmass = 1/(mass*friction_constant*dt) */
742 vn = invmass[n]*f[n][d]*dt
743 + sqrt(invmass[n])*rf[gt]*gmx_rng_gaussian_table(gaussrand);
746 xprime[n][d] = x[n][d]+vn*dt;
752 xprime[n][d] = x[n][d];
758 static void dump_it_all(FILE *fp,const char *title,
759 int natoms,rvec x[],rvec xp[],rvec v[],rvec f[])
764 fprintf(fp,"%s\n",title);
765 pr_rvecs(fp,0,"x",x,natoms);
766 pr_rvecs(fp,0,"xp",xp,natoms);
767 pr_rvecs(fp,0,"v",v,natoms);
768 pr_rvecs(fp,0,"f",f,natoms);
773 static void calc_ke_part_normal(rvec v[], t_grpopts *opts,t_mdatoms *md,
774 gmx_ekindata_t *ekind,t_nrnb *nrnb,gmx_bool bEkinAveVel,
775 gmx_bool bSaveEkinOld)
777 int start=md->start,homenr=md->homenr;
778 int g,d,n,m,ga=0,gt=0;
781 t_grp_tcstat *tcstat=ekind->tcstat;
782 t_grp_acc *grpstat=ekind->grpstat;
785 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
786 an option, but not supported now. Additionally, if we are doing iterations.
787 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
788 bSavEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't copy over the ekinh_old.
789 If FALSE, we overrwrite it.
792 /* group velocities are calculated in update_ekindata and
793 * accumulated in acumulate_groups.
794 * Now the partial global and groups ekin.
796 for(g=0; (g<opts->ngtc); g++)
800 copy_mat(tcstat[g].ekinh,tcstat[g].ekinh_old);
803 clear_mat(tcstat[g].ekinf);
805 clear_mat(tcstat[g].ekinh);
808 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
811 ekind->dekindl_old = ekind->dekindl;
814 for(n=start; (n<start+homenr); n++)
824 hm = 0.5*md->massT[n];
826 for(d=0; (d<DIM); d++)
828 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
830 for(d=0; (d<DIM); d++)
832 for (m=0;(m<DIM); m++)
834 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
837 tcstat[gt].ekinf[m][d]+=hm*v_corrt[m]*v_corrt[d];
841 tcstat[gt].ekinh[m][d]+=hm*v_corrt[m]*v_corrt[d];
845 if (md->nMassPerturbed && md->bPerturbed[n])
847 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt,v_corrt);
850 ekind->dekindl = dekindl;
851 inc_nrnb(nrnb,eNR_EKIN,homenr);
854 static void calc_ke_part_visc(matrix box,rvec x[],rvec v[],
855 t_grpopts *opts,t_mdatoms *md,
856 gmx_ekindata_t *ekind,
857 t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
859 int start=md->start,homenr=md->homenr;
863 t_grp_tcstat *tcstat=ekind->tcstat;
864 t_cos_acc *cosacc=&(ekind->cosacc);
869 for(g=0; g<opts->ngtc; g++)
871 copy_mat(ekind->tcstat[g].ekinh,ekind->tcstat[g].ekinh_old);
872 clear_mat(ekind->tcstat[g].ekinh);
874 ekind->dekindl_old = ekind->dekindl;
876 fac = 2*M_PI/box[ZZ][ZZ];
879 for(n=start; n<start+homenr; n++)
885 hm = 0.5*md->massT[n];
887 /* Note that the times of x and v differ by half a step */
888 /* MRS -- would have to be changed for VV */
889 cosz = cos(fac*x[n][ZZ]);
890 /* Calculate the amplitude of the new velocity profile */
891 mvcos += 2*cosz*md->massT[n]*v[n][XX];
893 copy_rvec(v[n],v_corrt);
894 /* Subtract the profile for the kinetic energy */
895 v_corrt[XX] -= cosz*cosacc->vcos;
896 for (d=0; (d<DIM); d++)
898 for (m=0; (m<DIM); m++)
900 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
903 tcstat[gt].ekinf[m][d]+=hm*v_corrt[m]*v_corrt[d];
907 tcstat[gt].ekinh[m][d]+=hm*v_corrt[m]*v_corrt[d];
911 if(md->nPerturbed && md->bPerturbed[n])
913 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt,v_corrt);
916 ekind->dekindl = dekindl;
917 cosacc->mvcos = mvcos;
919 inc_nrnb(nrnb,eNR_EKIN,homenr);
922 void calc_ke_part(t_state *state,t_grpopts *opts,t_mdatoms *md,
923 gmx_ekindata_t *ekind,t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
925 if (ekind->cosacc.cos_accel == 0)
927 calc_ke_part_normal(state->v,opts,md,ekind,nrnb,bEkinAveVel,bSaveEkinOld);
931 calc_ke_part_visc(state->box,state->x,state->v,opts,md,ekind,nrnb,bEkinAveVel,bSaveEkinOld);
935 void init_ekinstate(ekinstate_t *ekinstate,const t_inputrec *ir)
937 ekinstate->ekin_n = ir->opts.ngtc;
938 snew(ekinstate->ekinh,ekinstate->ekin_n);
939 snew(ekinstate->ekinf,ekinstate->ekin_n);
940 snew(ekinstate->ekinh_old,ekinstate->ekin_n);
941 snew(ekinstate->ekinscalef_nhc,ekinstate->ekin_n);
942 snew(ekinstate->ekinscaleh_nhc,ekinstate->ekin_n);
943 snew(ekinstate->vscale_nhc,ekinstate->ekin_n);
944 ekinstate->dekindl = 0;
945 ekinstate->mvcos = 0;
948 void update_ekinstate(ekinstate_t *ekinstate,gmx_ekindata_t *ekind)
952 for(i=0;i<ekinstate->ekin_n;i++)
954 copy_mat(ekind->tcstat[i].ekinh,ekinstate->ekinh[i]);
955 copy_mat(ekind->tcstat[i].ekinf,ekinstate->ekinf[i]);
956 copy_mat(ekind->tcstat[i].ekinh_old,ekinstate->ekinh_old[i]);
957 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
958 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
959 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
962 copy_mat(ekind->ekin,ekinstate->ekin_total);
963 ekinstate->dekindl = ekind->dekindl;
964 ekinstate->mvcos = ekind->cosacc.mvcos;
968 void restore_ekinstate_from_state(t_commrec *cr,
969 gmx_ekindata_t *ekind,ekinstate_t *ekinstate)
975 for(i=0;i<ekinstate->ekin_n;i++)
977 copy_mat(ekinstate->ekinh[i],ekind->tcstat[i].ekinh);
978 copy_mat(ekinstate->ekinf[i],ekind->tcstat[i].ekinf);
979 copy_mat(ekinstate->ekinh_old[i],ekind->tcstat[i].ekinh_old);
980 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
981 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
982 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
985 copy_mat(ekinstate->ekin_total,ekind->ekin);
987 ekind->dekindl = ekinstate->dekindl;
988 ekind->cosacc.mvcos = ekinstate->mvcos;
989 n = ekinstate->ekin_n;
994 gmx_bcast(sizeof(n),&n,cr);
997 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
998 ekind->tcstat[i].ekinh[0],cr);
999 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1000 ekind->tcstat[i].ekinf[0],cr);
1001 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1002 ekind->tcstat[i].ekinh_old[0],cr);
1004 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1005 &(ekind->tcstat[i].ekinscalef_nhc),cr);
1006 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1007 &(ekind->tcstat[i].ekinscaleh_nhc),cr);
1008 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1009 &(ekind->tcstat[i].vscale_nhc),cr);
1011 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1014 gmx_bcast(sizeof(ekind->dekindl),&ekind->dekindl,cr);
1015 gmx_bcast(sizeof(ekind->cosacc.mvcos),&ekind->cosacc.mvcos,cr);
1019 void set_deform_reference_box(gmx_update_t upd,gmx_large_int_t step,matrix box)
1021 upd->deformref_step = step;
1022 copy_mat(box,upd->deformref_box);
1025 static void deform(gmx_update_t upd,
1026 int start,int homenr,rvec x[],matrix box,matrix *scale_tot,
1027 const t_inputrec *ir,gmx_large_int_t step)
1029 matrix bnew,invbox,mu;
1033 elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1035 for(i=0; i<DIM; i++)
1037 for(j=0; j<DIM; j++)
1039 if (ir->deform[i][j] != 0)
1042 upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1046 /* We correct the off-diagonal elements,
1047 * which can grow indefinitely during shearing,
1048 * so the shifts do not get messed up.
1050 for(i=1; i<DIM; i++)
1052 for(j=i-1; j>=0; j--)
1054 while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1056 rvec_dec(bnew[i],bnew[j]);
1058 while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1060 rvec_inc(bnew[i],bnew[j]);
1064 m_inv_ur0(box,invbox);
1066 mmul_ur0(box,invbox,mu);
1068 for(i=start; i<start+homenr; i++)
1070 x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1071 x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1072 x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1076 /* The transposes of the scaling matrices are stored,
1077 * so we need to do matrix multiplication in the inverse order.
1079 mmul_ur0(*scale_tot,mu,*scale_tot);
1083 static void combine_forces(int nstlist,
1084 gmx_constr_t constr,
1085 t_inputrec *ir,t_mdatoms *md,t_idef *idef,
1086 t_commrec *cr,gmx_large_int_t step,t_state *state,
1087 int start,int nrend,
1088 rvec f[],rvec f_lr[],
1093 /* f contains the short-range forces + the long range forces
1094 * which are stored separately in f_lr.
1097 if (constr != NULL && !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1099 /* We need to constrain the LR forces separately,
1100 * because due to the different pre-factor for the SR and LR
1101 * forces in the update algorithm, we can not determine
1102 * the constraint force for the coordinate constraining.
1103 * Constrain only the additional LR part of the force.
1105 /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1106 constrain(NULL,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
1107 state->x,f_lr,f_lr,state->box,state->lambda,NULL,
1108 NULL,NULL,nrnb,econqForce,ir->epc==epcMTTK,state->veta,state->veta);
1111 /* Add nstlist-1 times the LR force to the sum of both forces
1112 * and store the result in forces_lr.
1115 for(i=start; i<nrend; i++)
1117 for(d=0; d<DIM; d++)
1119 f_lr[i][d] = f[i][d] + nm1*f_lr[i][d];
1124 void update_tcouple(FILE *fplog,
1125 gmx_large_int_t step,
1126 t_inputrec *inputrec,
1128 gmx_ekindata_t *ekind,
1129 gmx_wallcycle_t wcycle,
1135 gmx_bool bTCouple=FALSE;
1137 int i,start,end,homenr;
1139 /* if using vv, we do this elsewhere in the code */
1140 if (inputrec->etc != etcNO &&
1141 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec)))
1143 /* We should only couple after a step where energies were determined */
1144 bTCouple = (inputrec->nsttcouple == 1 ||
1145 do_per_step(step+inputrec->nsttcouple-1,
1146 inputrec->nsttcouple));
1151 dttc = inputrec->nsttcouple*inputrec->delta_t;
1153 switch (inputrec->etc)
1158 berendsen_tcoupl(inputrec,ekind,dttc);
1161 nosehoover_tcoupl(&(inputrec->opts),ekind,dttc,
1162 state->nosehoover_xi,state->nosehoover_vxi,MassQ);
1165 vrescale_tcoupl(inputrec,ekind,dttc,
1166 state->therm_integral,upd->sd->gaussrand);
1169 /* rescale in place here */
1170 if (EI_VV(inputrec->eI))
1172 rescale_velocities(ekind,md,md->start,md->start+md->homenr,state->v);
1177 /* Set the T scaling lambda to 1 to have no scaling */
1178 for(i=0; (i<inputrec->opts.ngtc); i++)
1180 ekind->tcstat[i].lambda = 1.0;
1185 void update_pcouple(FILE *fplog,
1186 gmx_large_int_t step,
1187 t_inputrec *inputrec,
1191 gmx_wallcycle_t wcycle,
1195 gmx_bool bPCouple=FALSE;
1199 /* if using vv, we do this elsewhere in the code */
1200 if (inputrec->epc != epcNO &&
1201 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec)))
1203 /* We should only couple after a step where energies were determined */
1204 bPCouple = (inputrec->nstpcouple == 1 ||
1205 do_per_step(step+inputrec->nstpcouple-1,
1206 inputrec->nstpcouple));
1209 clear_mat(pcoupl_mu);
1210 for(i=0; i<DIM; i++)
1212 pcoupl_mu[i][i] = 1.0;
1219 dtpc = inputrec->nstpcouple*inputrec->delta_t;
1221 switch (inputrec->epc)
1223 /* We can always pcoupl, even if we did not sum the energies
1224 * the previous step, since state->pres_prev is only updated
1225 * when the energies have been summed.
1229 case (epcBERENDSEN):
1232 berendsen_pcoupl(fplog,step,inputrec,dtpc,state->pres_prev,state->box,
1236 case (epcPARRINELLORAHMAN):
1237 parrinellorahman_pcoupl(fplog,step,inputrec,dtpc,state->pres_prev,
1238 state->box,state->box_rel,state->boxv,
1239 M,pcoupl_mu,bInitStep);
1247 static rvec *get_xprime(const t_state *state,gmx_update_t upd)
1249 if (state->nalloc > upd->xp_nalloc)
1251 upd->xp_nalloc = state->nalloc;
1252 srenew(upd->xp,upd->xp_nalloc);
1258 void update_constraints(FILE *fplog,
1259 gmx_large_int_t step,
1260 real *dvdlambda, /* FEP stuff */
1261 t_inputrec *inputrec, /* input record and box stuff */
1262 gmx_ekindata_t *ekind,
1266 rvec force[], /* forces on home particles */
1269 tensor vir, /* tensors for virial and ekin, needed for computing */
1272 gmx_wallcycle_t wcycle,
1274 gmx_constr_t constr,
1276 gmx_bool bFirstHalf,
1280 gmx_bool bExtended,bTrotter,bLastStep,bLog=FALSE,bEner=FALSE,bDoConstr=FALSE;
1283 int start,homenr,nrend,i,n,m,g,d;
1285 rvec *vbuf,*xprime=NULL;
1287 if (constr) {bDoConstr=TRUE;}
1288 if (bFirstHalf && !EI_VV(inputrec->eI)) {bDoConstr=FALSE;}
1290 /* for now, SD update is here -- though it really seems like it
1291 should be reformulated as a velocity verlet method, since it has two parts */
1294 homenr = md->homenr;
1295 nrend = start+homenr;
1297 dt = inputrec->delta_t;
1302 * APPLY CONSTRAINTS:
1305 * When doing PR pressure coupling we have to constrain the
1306 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1307 * it is enough to do this once though, since the relative velocities
1308 * after this will be normal to the bond vector
1313 /* clear out constraints before applying */
1314 clear_mat(vir_part);
1316 xprime = get_xprime(state,upd);
1318 bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1319 bLog = (do_per_step(step,inputrec->nstlog) || bLastStep || (step < 0));
1320 bEner = (do_per_step(step,inputrec->nstenergy) || bLastStep);
1321 /* Constrain the coordinates xprime */
1322 wallcycle_start(wcycle,ewcCONSTR);
1323 if (EI_VV(inputrec->eI) && bFirstHalf)
1325 constrain(NULL,bLog,bEner,constr,idef,
1326 inputrec,ekind,cr,step,1,md,
1327 state->x,state->v,state->v,
1328 state->box,state->lambda,dvdlambda,
1329 NULL,bCalcVir ? &vir_con : NULL,nrnb,econqVeloc,
1330 inputrec->epc==epcMTTK,state->veta,vetanew);
1334 constrain(NULL,bLog,bEner,constr,idef,
1335 inputrec,ekind,cr,step,1,md,
1336 state->x,xprime,NULL,
1337 state->box,state->lambda,dvdlambda,
1338 state->v,bCalcVir ? &vir_con : NULL ,nrnb,econqCoord,
1339 inputrec->epc==epcMTTK,state->veta,state->veta);
1341 wallcycle_stop(wcycle,ewcCONSTR);
1345 dump_it_all(fplog,"After Shake",
1346 state->natoms,state->x,xprime,state->v,force);
1350 if (inputrec->eI == eiSD2)
1352 /* A correction factor eph is needed for the SD constraint force */
1353 /* Here we can, unfortunately, not have proper corrections
1354 * for different friction constants, so we use the first one.
1356 for(i=0; i<DIM; i++)
1358 for(m=0; m<DIM; m++)
1360 vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1366 m_add(vir_part,vir_con,vir_part);
1370 pr_rvecs(debug,0,"constraint virial",vir_part,DIM);
1376 if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1378 xprime = get_xprime(state,upd);
1380 /* The second part of the SD integration */
1381 do_update_sd2(upd->sd,FALSE,start,homenr,
1382 inputrec->opts.acc,inputrec->opts.nFreeze,
1383 md->invmass,md->ptype,
1384 md->cFREEZE,md->cACC,md->cTC,
1385 state->x,xprime,state->v,force,state->sd_X,
1386 inputrec->opts.ngtc,inputrec->opts.tau_t,
1387 inputrec->opts.ref_t,FALSE);
1388 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1392 /* Constrain the coordinates xprime */
1393 wallcycle_start(wcycle,ewcCONSTR);
1394 constrain(NULL,bLog,bEner,constr,idef,
1395 inputrec,NULL,cr,step,1,md,
1396 state->x,xprime,NULL,
1397 state->box,state->lambda,dvdlambda,
1398 NULL,NULL,nrnb,econqCoord,FALSE,0,0);
1399 wallcycle_stop(wcycle,ewcCONSTR);
1403 /* We must always unshift after updating coordinates; if we did not shake
1404 x was shifted in do_force */
1406 if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1408 if (graph && (graph->nnodes > 0))
1410 unshift_x(graph,state->box,state->x,upd->xp);
1411 if (TRICLINIC(state->box))
1413 inc_nrnb(nrnb,eNR_SHIFTX,2*graph->nnodes);
1417 inc_nrnb(nrnb,eNR_SHIFTX,graph->nnodes);
1419 copy_rvecn(upd->xp,state->x,start,graph->start);
1420 copy_rvecn(upd->xp,state->x,graph->start+graph->nnodes,nrend);
1424 copy_rvecn(upd->xp,state->x,start,nrend);
1427 dump_it_all(fplog,"After unshift",
1428 state->natoms,state->x,upd->xp,state->v,force);
1430 /* ############# END the update of velocities and positions ######### */
1433 void update_box(FILE *fplog,
1434 gmx_large_int_t step,
1435 t_inputrec *inputrec, /* input record and box stuff */
1439 rvec force[], /* forces on home particles */
1443 gmx_wallcycle_t wcycle,
1446 gmx_bool bFirstHalf)
1448 gmx_bool bExtended,bTrotter,bLastStep,bLog=FALSE,bEner=FALSE;
1451 int start,homenr,nrend,i,n,m,g;
1455 homenr = md->homenr;
1456 nrend = start+homenr;
1459 (inputrec->etc == etcNOSEHOOVER) ||
1460 (inputrec->epc == epcPARRINELLORAHMAN) ||
1461 (inputrec->epc == epcMTTK);
1463 dt = inputrec->delta_t;
1467 /* now update boxes */
1468 switch (inputrec->epc) {
1471 case (epcBERENDSEN):
1472 berendsen_pscale(inputrec,pcoupl_mu,state->box,state->box_rel,
1473 start,homenr,state->x,md->cFREEZE,nrnb);
1475 case (epcPARRINELLORAHMAN):
1476 /* The box velocities were updated in do_pr_pcoupl in the update
1477 * iteration, but we dont change the box vectors until we get here
1478 * since we need to be able to shift/unshift above.
1484 state->box[i][m] += dt*state->boxv[i][m];
1487 preserve_box_shape(inputrec,state->box_rel,state->box);
1489 /* Scale the coordinates */
1490 for(n=start; (n<start+homenr); n++)
1492 tmvmul_ur0(pcoupl_mu,state->x[n],state->x[n]);
1496 switch (inputrec->epct)
1498 case (epctISOTROPIC):
1499 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1500 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1501 Side length scales as exp(veta*dt) */
1503 msmul(state->box,exp(state->veta*dt),state->box);
1505 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1506 o If we assume isotropic scaling, and box length scaling
1507 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1508 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1509 determinant of B is L^DIM det(M), and the determinant
1510 of dB/dt is (dL/dT)^DIM det (M). veta will be
1511 (det(dB/dT)/det(B))^(1/3). Then since M =
1512 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1514 msmul(state->box,state->veta,state->boxv);
1524 if ((!(IR_NPT_TROTTER(inputrec))) && scale_tot)
1526 /* The transposes of the scaling matrices are stored,
1527 * therefore we need to reverse the order in the multiplication.
1529 mmul_ur0(*scale_tot,pcoupl_mu,*scale_tot);
1532 if (DEFORM(*inputrec))
1534 deform(upd,start,homenr,state->x,state->box,scale_tot,inputrec,step);
1537 dump_it_all(fplog,"After update",
1538 state->natoms,state->x,upd->xp,state->v,force);
1541 void update_coords(FILE *fplog,
1542 gmx_large_int_t step,
1543 t_inputrec *inputrec, /* input record and box stuff */
1546 rvec *f, /* forces on home particles */
1550 gmx_ekindata_t *ekind,
1552 gmx_wallcycle_t wcycle,
1556 t_commrec *cr, /* these shouldn't be here -- need to think about it */
1558 gmx_constr_t constr,
1561 gmx_bool bExtended,bNH,bPR,bTrotter,bLastStep,bLog=FALSE,bEner=FALSE;
1563 real *imass,*imassin;
1566 int start,homenr,nrend,i,j,d,n,m,g;
1567 int blen0,blen1,iatom,jatom,nshake,nsettle,nconstr,nexpand;
1570 rvec *vcom,*xcom,*vall,*xall,*xin,*vin,*forcein,*fall,*xpall,*xprimein,*xprime;
1573 /* Running the velocity half does nothing except for velocity verlet */
1574 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1575 !EI_VV(inputrec->eI))
1577 gmx_incons("update_coords called for velocity without VV integrator");
1581 homenr = md->homenr;
1582 nrend = start+homenr;
1584 xprime = get_xprime(state,upd);
1586 dt = inputrec->delta_t;
1589 /* We need to update the NMR restraint history when time averaging is used */
1590 if (state->flags & (1<<estDISRE_RM3TAV))
1592 update_disres_history(fcd,&state->hist);
1594 if (state->flags & (1<<estORIRE_DTAV))
1596 update_orires_history(fcd,&state->hist);
1599 bNH = inputrec->etc == etcNOSEHOOVER;
1600 bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1602 bExtended = bNH || bPR;
1604 if (bDoLR && inputrec->nstlist > 1 && !EI_VV(inputrec->eI)) /* get this working with VV? */
1606 /* Store the total force + nstlist-1 times the LR force
1607 * in forces_lr, so it can be used in a normal update algorithm
1608 * to produce twin time stepping.
1610 /* is this correct in the new construction? MRS */
1611 combine_forces(inputrec->nstlist,constr,inputrec,md,idef,cr,step,state,
1612 start,nrend,f,f_lr,nrnb);
1620 /* ############# START The update of velocities and positions ######### */
1622 dump_it_all(fplog,"Before update",
1623 state->natoms,state->x,xprime,state->v,force);
1625 switch (inputrec->eI) {
1627 if (ekind->cosacc.cos_accel == 0) {
1628 /* use normal version of update */
1629 do_update_md(start,nrend,dt,
1630 ekind->tcstat,ekind->grpstat,state->nosehoover_vxi,
1631 inputrec->opts.acc,inputrec->opts.nFreeze,md->invmass,md->ptype,
1632 md->cFREEZE,md->cACC,md->cTC,
1633 state->x,xprime,state->v,force,M,
1638 do_update_visc(start,nrend,dt,
1639 ekind->tcstat,md->invmass,state->nosehoover_vxi,
1640 md->ptype,md->cTC,state->x,xprime,state->v,force,M,
1641 state->box,ekind->cosacc.cos_accel,ekind->cosacc.vcos,bNH,bPR);
1645 do_update_sd1(upd->sd,start,homenr,dt,
1646 inputrec->opts.acc,inputrec->opts.nFreeze,
1647 md->invmass,md->ptype,
1648 md->cFREEZE,md->cACC,md->cTC,
1649 state->x,xprime,state->v,force,state->sd_X,
1650 inputrec->opts.ngtc,inputrec->opts.tau_t,inputrec->opts.ref_t);
1653 /* The SD update is done in 2 parts, because an extra constraint step
1656 do_update_sd2(upd->sd,bInitStep,start,homenr,
1657 inputrec->opts.acc,inputrec->opts.nFreeze,
1658 md->invmass,md->ptype,
1659 md->cFREEZE,md->cACC,md->cTC,
1660 state->x,xprime,state->v,force,state->sd_X,
1661 inputrec->opts.ngtc,inputrec->opts.tau_t,inputrec->opts.ref_t,
1665 do_update_bd(start,nrend,dt,
1666 inputrec->opts.nFreeze,md->invmass,md->ptype,
1667 md->cFREEZE,md->cTC,
1668 state->x,xprime,state->v,force,
1670 inputrec->opts.ngtc,inputrec->opts.tau_t,inputrec->opts.ref_t,
1671 upd->sd->bd_rf,upd->sd->gaussrand);
1675 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
1676 switch (UpdatePart) {
1679 do_update_vv_vel(start,nrend,dt,
1680 ekind->tcstat,ekind->grpstat,
1681 inputrec->opts.acc,inputrec->opts.nFreeze,
1682 md->invmass,md->ptype,
1683 md->cFREEZE,md->cACC,
1685 bExtended,state->veta,alpha);
1688 do_update_vv_pos(start,nrend,dt,
1689 ekind->tcstat,ekind->grpstat,
1690 inputrec->opts.acc,inputrec->opts.nFreeze,
1691 md->invmass,md->ptype,md->cFREEZE,
1692 state->x,xprime,state->v,force,
1693 bExtended,state->veta,alpha);
1698 gmx_fatal(FARGS,"Don't know how to update coordinates");
1704 void correct_ekin(FILE *log,int start,int end,rvec v[],rvec vcm,real mass[],
1705 real tmass,tensor ekin)
1708 * This is a debugging routine. It should not be called for production code
1710 * The kinetic energy should calculated according to:
1711 * Ekin = 1/2 m (v-vcm)^2
1712 * However the correction is not always applied, since vcm may not be
1713 * known in time and we compute
1714 * Ekin' = 1/2 m v^2 instead
1715 * This can be corrected afterwards by computing
1716 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
1718 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
1725 /* Local particles */
1728 /* Processor dependent part. */
1730 for(i=start; (i<end); i++)
1734 for(j=0; (j<DIM); j++)
1740 svmul(1/tmass,vcm,vcm);
1741 svmul(0.5,vcm,hvcm);
1743 for(j=0; (j<DIM); j++)
1745 for(k=0; (k<DIM); k++)
1747 dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
1750 pr_rvecs(log,0,"dekin",dekin,DIM);
1751 pr_rvecs(log,0," ekin", ekin,DIM);
1752 fprintf(log,"dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
1753 trace(dekin),trace(ekin),vcm[XX],vcm[YY],vcm[ZZ]);
1754 fprintf(log,"mv = (%8.4f %8.4f %8.4f)\n",
1755 mv[XX],mv[YY],mv[ZZ]);