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;
100 /* andersen temperature control stuff */
101 gmx_bool *randomize_group;
105 typedef struct gmx_update
108 /* xprime for constraint algorithms */
112 /* variable size arrays for andersen */
115 gmx_bool randatom_list_init;
117 /* Variables for the deform algorithm */
118 gmx_large_int_t deformref_step;
119 matrix deformref_box;
123 static void do_update_md(int start,int nrend,double dt,
124 t_grp_tcstat *tcstat,t_grp_acc *gstat,double nh_vxi[],
125 rvec accel[],ivec nFreeze[],real invmass[],
126 unsigned short ptype[],unsigned short cFREEZE[],
127 unsigned short cACC[],unsigned short cTC[],
128 rvec x[],rvec xprime[],rvec v[],
130 gmx_bool bNH,gmx_bool bPR)
135 real vn,vv,va,vb,vnrel;
141 /* Update with coupling to extended ensembles, used for
142 * Nose-Hoover and Parrinello-Rahman coupling
143 * Nose-Hoover uses the reversible leap-frog integrator from
144 * Holian et al. Phys Rev E 52(3) : 2338, 1995
146 for(n=start; n<nrend; n++)
161 lg = tcstat[gt].lambda;
165 rvec_sub(v[n],gstat[ga].u,vrel);
169 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
171 vnrel = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
172 - iprod(M[d],vrel)))/(1 + 0.5*vxi*dt);
173 /* do not scale the mean velocities u */
174 vn = gstat[ga].u[d] + accel[ga][d]*dt + vnrel;
176 xprime[n][d] = x[n][d]+vn*dt;
181 xprime[n][d] = x[n][d];
188 /* Classic version of update, used with berendsen coupling */
189 for(n=start; n<nrend; n++)
191 w_dt = invmass[n]*dt;
204 lg = tcstat[gt].lambda;
209 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
211 vv = lg*vn + f[n][d]*w_dt;
213 /* do not scale the mean velocities u */
215 va = vv + accel[ga][d]*dt;
216 vb = va + (1.0-lg)*u;
218 xprime[n][d] = x[n][d]+vb*dt;
223 xprime[n][d] = x[n][d];
230 static void do_update_vv_vel(int start,int nrend,double dt,
231 t_grp_tcstat *tcstat,t_grp_acc *gstat,
232 rvec accel[],ivec nFreeze[],real invmass[],
233 unsigned short ptype[],unsigned short cFREEZE[],
234 unsigned short cACC[],rvec v[],rvec f[],
235 gmx_bool bExtended, real veta, real alpha)
240 real u,vn,vv,va,vb,vnrel;
246 g = 0.25*dt*veta*alpha;
248 mv2 = series_sinhx(g);
255 for(n=start; n<nrend; n++)
257 w_dt = invmass[n]*dt;
269 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
271 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
279 } /* do_update_vv_vel */
281 static void do_update_vv_pos(int start,int nrend,double dt,
282 t_grp_tcstat *tcstat,t_grp_acc *gstat,
283 rvec accel[],ivec nFreeze[],real invmass[],
284 unsigned short ptype[],unsigned short cFREEZE[],
285 rvec x[],rvec xprime[],rvec v[],
286 rvec f[],gmx_bool bExtended, real veta, real alpha)
293 /* Would it make more sense if Parrinello-Rahman was put here? */
298 mr2 = series_sinhx(g);
304 for(n=start; n<nrend; n++) {
313 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
315 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
319 xprime[n][d] = x[n][d];
323 }/* do_update_vv_pos */
325 static void do_update_visc(int start,int nrend,double dt,
326 t_grp_tcstat *tcstat,real invmass[],double nh_vxi[],
327 unsigned short ptype[],unsigned short cTC[],
328 rvec x[],rvec xprime[],rvec v[],
329 rvec f[],matrix M,matrix box,real
331 gmx_bool bNH,gmx_bool bPR)
341 fac = 2*M_PI/(box[ZZ][ZZ]);
344 /* Update with coupling to extended ensembles, used for
345 * Nose-Hoover and Parrinello-Rahman coupling
347 for(n=start; n<nrend; n++) {
353 lg = tcstat[gt].lambda;
354 cosz = cos(fac*x[n][ZZ]);
356 copy_rvec(v[n],vrel);
368 if((ptype[n] != eptVSite) && (ptype[n] != eptShell))
370 vn = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
371 - iprod(M[d],vrel)))/(1 + 0.5*vxi*dt);
374 vn += vc + dt*cosz*cos_accel;
377 xprime[n][d] = x[n][d]+vn*dt;
381 xprime[n][d] = x[n][d];
388 /* Classic version of update, used with berendsen coupling */
389 for(n=start; n<nrend; n++)
391 w_dt = invmass[n]*dt;
396 lg = tcstat[gt].lambda;
397 cosz = cos(fac*x[n][ZZ]);
403 if((ptype[n] != eptVSite) && (ptype[n] != eptShell))
408 /* Do not scale the cosine velocity profile */
409 vv = vc + lg*(vn - vc + f[n][d]*w_dt);
410 /* Add the cosine accelaration profile */
411 vv += dt*cosz*cos_accel;
415 vv = lg*(vn + f[n][d]*w_dt);
418 xprime[n][d] = x[n][d]+vv*dt;
423 xprime[n][d] = x[n][d];
430 static gmx_stochd_t *init_stochd(FILE *fplog,t_inputrec *ir)
439 /* Initiate random number generator for langevin type dynamics,
440 * for BD, SD or velocity rescaling temperature coupling.
442 sd->gaussrand = gmx_rng_init(ir->ld_seed);
444 ngtc = ir->opts.ngtc;
448 snew(sd->bd_rf,ngtc);
450 else if (EI_SD(ir->eI))
453 snew(sd->sdsig,ngtc);
456 for(n=0; n<ngtc; n++)
458 if (ir->opts.tau_t[n] > 0)
460 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
461 sdc[n].eph = exp(sdc[n].gdt/2);
462 sdc[n].emh = exp(-sdc[n].gdt/2);
463 sdc[n].em = exp(-sdc[n].gdt);
467 /* No friction and noise on this group */
473 if (sdc[n].gdt >= 0.05)
475 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
476 - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
477 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
478 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
483 /* Seventh order expansions for small y */
484 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
485 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))));
486 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
489 fprintf(debug,"SD const tc-grp %d: b %g c %g d %g\n",
490 n,sdc[n].b,sdc[n].c,sdc[n].d);
493 else if (ETC_ANDERSEN(ir->etc))
502 snew(sd->randomize_group,ngtc);
503 snew(sd->boltzfac,ngtc);
505 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
506 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
508 for (n=0;n<ngtc;n++) {
509 reft = max(0.0,opts->ref_t[n]);
510 if ((opts->tau_t[n] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
512 sd->randomize_group[n] = TRUE;
513 sd->boltzfac[n] = BOLTZ*opts->ref_t[n];
515 sd->randomize_group[n] = FALSE;
522 void get_stochd_state(gmx_update_t upd,t_state *state)
524 gmx_rng_get_state(upd->sd->gaussrand,state->ld_rng,state->ld_rngi);
527 void set_stochd_state(gmx_update_t upd,t_state *state)
529 gmx_rng_set_state(upd->sd->gaussrand,state->ld_rng,state->ld_rngi[0]);
532 gmx_update_t init_update(FILE *fplog,t_inputrec *ir)
538 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
540 upd->sd = init_stochd(fplog,ir);
545 upd->randatom = NULL;
546 upd->randatom_list = NULL;
547 upd->randatom_list_init = FALSE; /* we have not yet cleared the data structure at this point */
552 static void do_update_sd1(gmx_stochd_t *sd,
553 int start,int homenr,double dt,
554 rvec accel[],ivec nFreeze[],
555 real invmass[],unsigned short ptype[],
556 unsigned short cFREEZE[],unsigned short cACC[],
557 unsigned short cTC[],
558 rvec x[],rvec xprime[],rvec v[],rvec f[],
560 int ngtc,real tau_t[],real ref_t[])
572 if (homenr > sd->sd_V_nalloc)
574 sd->sd_V_nalloc = over_alloc_dd(homenr);
575 srenew(sd->sd_V,sd->sd_V_nalloc);
577 gaussrand = sd->gaussrand;
579 for(n=0; n<ngtc; n++)
582 /* The mass is encounted for later, since this differs per atom */
583 sig[n].V = sqrt(kT*(1 - sdc[n].em*sdc[n].em));
586 for(n=start; n<start+homenr; n++)
588 ism = sqrt(invmass[n]);
604 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
606 sd_V = ism*sig[gt].V*gmx_rng_gaussian_table(gaussrand);
608 v[n][d] = v[n][d]*sdc[gt].em
609 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
612 xprime[n][d] = x[n][d] + v[n][d]*dt;
617 xprime[n][d] = x[n][d];
623 static void do_update_sd2(gmx_stochd_t *sd,gmx_bool bInitStep,
624 int start,int homenr,
625 rvec accel[],ivec nFreeze[],
626 real invmass[],unsigned short ptype[],
627 unsigned short cFREEZE[],unsigned short cACC[],
628 unsigned short cTC[],
629 rvec x[],rvec xprime[],rvec v[],rvec f[],
631 int ngtc,real tau_t[],real ref_t[],
636 /* The random part of the velocity update, generated in the first
637 * half of the update, needs to be remembered for the second half.
649 if (homenr > sd->sd_V_nalloc)
651 sd->sd_V_nalloc = over_alloc_dd(homenr);
652 srenew(sd->sd_V,sd->sd_V_nalloc);
655 gaussrand = sd->gaussrand;
659 for (n=0; n<ngtc; n++)
662 /* The mass is encounted for later, since this differs per atom */
663 sig[n].V = sqrt(kT*(1-sdc[n].em));
664 sig[n].X = sqrt(kT*sqr(tau_t[n])*sdc[n].c);
665 sig[n].Yv = sqrt(kT*sdc[n].b/sdc[n].c);
666 sig[n].Yx = sqrt(kT*sqr(tau_t[n])*sdc[n].b/(1-sdc[n].em));
670 for (n=start; n<start+homenr; n++)
672 ism = sqrt(invmass[n]);
692 if((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
698 sd_X[n][d] = ism*sig[gt].X*gmx_rng_gaussian_table(gaussrand);
700 Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
701 + ism*sig[gt].Yv*gmx_rng_gaussian_table(gaussrand);
702 sd_V[n-start][d] = ism*sig[gt].V*gmx_rng_gaussian_table(gaussrand);
704 v[n][d] = vn*sdc[gt].em
705 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
706 + sd_V[n-start][d] - sdc[gt].em*Vmh;
708 xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
713 /* Correct the velocities for the constraints.
714 * This operation introduces some inaccuracy,
715 * since the velocity is determined from differences in coordinates.
718 (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
720 Xmh = sd_V[n-start][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
721 + ism*sig[gt].Yx*gmx_rng_gaussian_table(gaussrand);
722 sd_X[n][d] = ism*sig[gt].X*gmx_rng_gaussian_table(gaussrand);
724 xprime[n][d] += sd_X[n][d] - Xmh;
733 xprime[n][d] = x[n][d];
740 static void do_update_bd(int start,int nrend,double dt,
742 real invmass[],unsigned short ptype[],
743 unsigned short cFREEZE[],unsigned short cTC[],
744 rvec x[],rvec xprime[],rvec v[],
745 rvec f[],real friction_coefficient,
746 int ngtc,real tau_t[],real ref_t[],
747 real *rf,gmx_rng_t gaussrand)
749 /* note -- these appear to be full step velocities . . . */
755 if (friction_coefficient != 0)
757 invfr = 1.0/friction_coefficient;
758 for(n=0; n<ngtc; n++)
760 rf[n] = sqrt(2.0*BOLTZ*ref_t[n]/(friction_coefficient*dt));
765 for(n=0; n<ngtc; n++)
767 rf[n] = sqrt(2.0*BOLTZ*ref_t[n]);
770 for(n=start; (n<nrend); n++)
780 for(d=0; (d<DIM); d++)
782 if((ptype[n]!=eptVSite) && (ptype[n]!=eptShell) && !nFreeze[gf][d])
784 if (friction_coefficient != 0) {
785 vn = invfr*f[n][d] + rf[gt]*gmx_rng_gaussian_table(gaussrand);
789 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
790 vn = 0.5*invmass[n]*f[n][d]*dt
791 + sqrt(0.5*invmass[n])*rf[gt]*gmx_rng_gaussian_table(gaussrand);
795 xprime[n][d] = x[n][d]+vn*dt;
800 xprime[n][d] = x[n][d];
806 static void dump_it_all(FILE *fp,const char *title,
807 int natoms,rvec x[],rvec xp[],rvec v[],rvec f[])
812 fprintf(fp,"%s\n",title);
813 pr_rvecs(fp,0,"x",x,natoms);
814 pr_rvecs(fp,0,"xp",xp,natoms);
815 pr_rvecs(fp,0,"v",v,natoms);
816 pr_rvecs(fp,0,"f",f,natoms);
821 static void calc_ke_part_normal(rvec v[], t_grpopts *opts,t_mdatoms *md,
822 gmx_ekindata_t *ekind,t_nrnb *nrnb,gmx_bool bEkinAveVel,
823 gmx_bool bSaveEkinOld)
825 int start=md->start,homenr=md->homenr;
826 int g,d,n,m,ga=0,gt=0;
829 t_grp_tcstat *tcstat=ekind->tcstat;
830 t_grp_acc *grpstat=ekind->grpstat;
833 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
834 an option, but not supported now. Additionally, if we are doing iterations.
835 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
836 bSavEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't copy over the ekinh_old.
837 If FALSE, we overrwrite it.
840 /* group velocities are calculated in update_ekindata and
841 * accumulated in acumulate_groups.
842 * Now the partial global and groups ekin.
844 for(g=0; (g<opts->ngtc); g++)
848 copy_mat(tcstat[g].ekinh,tcstat[g].ekinh_old);
851 clear_mat(tcstat[g].ekinf);
853 clear_mat(tcstat[g].ekinh);
856 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
859 ekind->dekindl_old = ekind->dekindl;
862 for(n=start; (n<start+homenr); n++)
872 hm = 0.5*md->massT[n];
874 for(d=0; (d<DIM); d++)
876 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
878 for(d=0; (d<DIM); d++)
880 for (m=0;(m<DIM); m++)
882 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
885 tcstat[gt].ekinf[m][d]+=hm*v_corrt[m]*v_corrt[d];
889 tcstat[gt].ekinh[m][d]+=hm*v_corrt[m]*v_corrt[d];
893 if (md->nMassPerturbed && md->bPerturbed[n])
895 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt,v_corrt);
898 ekind->dekindl = dekindl;
899 inc_nrnb(nrnb,eNR_EKIN,homenr);
902 static void calc_ke_part_visc(matrix box,rvec x[],rvec v[],
903 t_grpopts *opts,t_mdatoms *md,
904 gmx_ekindata_t *ekind,
905 t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
907 int start=md->start,homenr=md->homenr;
911 t_grp_tcstat *tcstat=ekind->tcstat;
912 t_cos_acc *cosacc=&(ekind->cosacc);
917 for(g=0; g<opts->ngtc; g++)
919 copy_mat(ekind->tcstat[g].ekinh,ekind->tcstat[g].ekinh_old);
920 clear_mat(ekind->tcstat[g].ekinh);
922 ekind->dekindl_old = ekind->dekindl;
924 fac = 2*M_PI/box[ZZ][ZZ];
927 for(n=start; n<start+homenr; n++)
933 hm = 0.5*md->massT[n];
935 /* Note that the times of x and v differ by half a step */
936 /* MRS -- would have to be changed for VV */
937 cosz = cos(fac*x[n][ZZ]);
938 /* Calculate the amplitude of the new velocity profile */
939 mvcos += 2*cosz*md->massT[n]*v[n][XX];
941 copy_rvec(v[n],v_corrt);
942 /* Subtract the profile for the kinetic energy */
943 v_corrt[XX] -= cosz*cosacc->vcos;
944 for (d=0; (d<DIM); d++)
946 for (m=0; (m<DIM); m++)
948 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
951 tcstat[gt].ekinf[m][d]+=hm*v_corrt[m]*v_corrt[d];
955 tcstat[gt].ekinh[m][d]+=hm*v_corrt[m]*v_corrt[d];
959 if(md->nPerturbed && md->bPerturbed[n])
961 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt,v_corrt);
964 ekind->dekindl = dekindl;
965 cosacc->mvcos = mvcos;
967 inc_nrnb(nrnb,eNR_EKIN,homenr);
970 void calc_ke_part(t_state *state,t_grpopts *opts,t_mdatoms *md,
971 gmx_ekindata_t *ekind,t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
973 if (ekind->cosacc.cos_accel == 0)
975 calc_ke_part_normal(state->v,opts,md,ekind,nrnb,bEkinAveVel,bSaveEkinOld);
979 calc_ke_part_visc(state->box,state->x,state->v,opts,md,ekind,nrnb,bEkinAveVel,bSaveEkinOld);
983 extern void init_ekinstate(ekinstate_t *ekinstate,const t_inputrec *ir)
985 ekinstate->ekin_n = ir->opts.ngtc;
986 snew(ekinstate->ekinh,ekinstate->ekin_n);
987 snew(ekinstate->ekinf,ekinstate->ekin_n);
988 snew(ekinstate->ekinh_old,ekinstate->ekin_n);
989 snew(ekinstate->ekinscalef_nhc,ekinstate->ekin_n);
990 snew(ekinstate->ekinscaleh_nhc,ekinstate->ekin_n);
991 snew(ekinstate->vscale_nhc,ekinstate->ekin_n);
992 ekinstate->dekindl = 0;
993 ekinstate->mvcos = 0;
996 void update_ekinstate(ekinstate_t *ekinstate,gmx_ekindata_t *ekind)
1000 for(i=0;i<ekinstate->ekin_n;i++)
1002 copy_mat(ekind->tcstat[i].ekinh,ekinstate->ekinh[i]);
1003 copy_mat(ekind->tcstat[i].ekinf,ekinstate->ekinf[i]);
1004 copy_mat(ekind->tcstat[i].ekinh_old,ekinstate->ekinh_old[i]);
1005 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1006 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1007 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1010 copy_mat(ekind->ekin,ekinstate->ekin_total);
1011 ekinstate->dekindl = ekind->dekindl;
1012 ekinstate->mvcos = ekind->cosacc.mvcos;
1016 void restore_ekinstate_from_state(t_commrec *cr,
1017 gmx_ekindata_t *ekind,ekinstate_t *ekinstate)
1023 for(i=0;i<ekinstate->ekin_n;i++)
1025 copy_mat(ekinstate->ekinh[i],ekind->tcstat[i].ekinh);
1026 copy_mat(ekinstate->ekinf[i],ekind->tcstat[i].ekinf);
1027 copy_mat(ekinstate->ekinh_old[i],ekind->tcstat[i].ekinh_old);
1028 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1029 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1030 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1033 copy_mat(ekinstate->ekin_total,ekind->ekin);
1035 ekind->dekindl = ekinstate->dekindl;
1036 ekind->cosacc.mvcos = ekinstate->mvcos;
1037 n = ekinstate->ekin_n;
1042 gmx_bcast(sizeof(n),&n,cr);
1045 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1046 ekind->tcstat[i].ekinh[0],cr);
1047 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1048 ekind->tcstat[i].ekinf[0],cr);
1049 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1050 ekind->tcstat[i].ekinh_old[0],cr);
1052 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1053 &(ekind->tcstat[i].ekinscalef_nhc),cr);
1054 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1055 &(ekind->tcstat[i].ekinscaleh_nhc),cr);
1056 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1057 &(ekind->tcstat[i].vscale_nhc),cr);
1059 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1062 gmx_bcast(sizeof(ekind->dekindl),&ekind->dekindl,cr);
1063 gmx_bcast(sizeof(ekind->cosacc.mvcos),&ekind->cosacc.mvcos,cr);
1067 void set_deform_reference_box(gmx_update_t upd,gmx_large_int_t step,matrix box)
1069 upd->deformref_step = step;
1070 copy_mat(box,upd->deformref_box);
1073 static void deform(gmx_update_t upd,
1074 int start,int homenr,rvec x[],matrix box,matrix *scale_tot,
1075 const t_inputrec *ir,gmx_large_int_t step)
1077 matrix bnew,invbox,mu;
1081 elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1083 for(i=0; i<DIM; i++)
1085 for(j=0; j<DIM; j++)
1087 if (ir->deform[i][j] != 0)
1090 upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1094 /* We correct the off-diagonal elements,
1095 * which can grow indefinitely during shearing,
1096 * so the shifts do not get messed up.
1098 for(i=1; i<DIM; i++)
1100 for(j=i-1; j>=0; j--)
1102 while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1104 rvec_dec(bnew[i],bnew[j]);
1106 while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1108 rvec_inc(bnew[i],bnew[j]);
1112 m_inv_ur0(box,invbox);
1114 mmul_ur0(box,invbox,mu);
1116 for(i=start; i<start+homenr; i++)
1118 x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1119 x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1120 x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1124 /* The transposes of the scaling matrices are stored,
1125 * so we need to do matrix multiplication in the inverse order.
1127 mmul_ur0(*scale_tot,mu,*scale_tot);
1131 static void combine_forces(int nstlist,
1132 gmx_constr_t constr,
1133 t_inputrec *ir,t_mdatoms *md,t_idef *idef,
1134 t_commrec *cr,gmx_large_int_t step,t_state *state,
1135 int start,int nrend,
1136 rvec f[],rvec f_lr[],
1141 /* f contains the short-range forces + the long range forces
1142 * which are stored separately in f_lr.
1145 if (constr != NULL && !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1147 /* We need to constrain the LR forces separately,
1148 * because due to the different pre-factor for the SR and LR
1149 * forces in the update algorithm, we can not determine
1150 * the constraint force for the coordinate constraining.
1151 * Constrain only the additional LR part of the force.
1153 /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1154 constrain(NULL,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
1155 state->x,f_lr,f_lr,state->box,state->lambda[efptBONDED],NULL,
1156 NULL,NULL,nrnb,econqForce,ir->epc==epcMTTK,state->veta,state->veta);
1159 /* Add nstlist-1 times the LR force to the sum of both forces
1160 * and store the result in forces_lr.
1163 for(i=start; i<nrend; i++)
1165 for(d=0; d<DIM; d++)
1167 f_lr[i][d] = f[i][d] + nm1*f_lr[i][d];
1172 void update_tcouple(FILE *fplog,
1173 gmx_large_int_t step,
1174 t_inputrec *inputrec,
1176 gmx_ekindata_t *ekind,
1177 gmx_wallcycle_t wcycle,
1183 gmx_bool bTCouple=FALSE;
1185 int i,start,end,homenr,offset;
1187 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1188 if (inputrec->etc != etcNO &&
1189 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec)))
1191 /* We should only couple after a step where energies were determined (for leapfrog versions)
1192 or the step energies are determined, for velocity verlet versions */
1194 if (EI_VV(inputrec->eI)) {
1199 bTCouple = (inputrec->nsttcouple == 1 ||
1200 do_per_step(step+inputrec->nsttcouple-offset,
1201 inputrec->nsttcouple));
1206 dttc = inputrec->nsttcouple*inputrec->delta_t;
1208 switch (inputrec->etc)
1213 berendsen_tcoupl(inputrec,ekind,dttc);
1216 nosehoover_tcoupl(&(inputrec->opts),ekind,dttc,
1217 state->nosehoover_xi,state->nosehoover_vxi,MassQ);
1220 vrescale_tcoupl(inputrec,ekind,dttc,
1221 state->therm_integral,upd->sd->gaussrand);
1224 /* rescale in place here */
1225 if (EI_VV(inputrec->eI))
1227 rescale_velocities(ekind,md,md->start,md->start+md->homenr,state->v);
1232 /* Set the T scaling lambda to 1 to have no scaling */
1233 for(i=0; (i<inputrec->opts.ngtc); i++)
1235 ekind->tcstat[i].lambda = 1.0;
1240 void update_pcouple(FILE *fplog,
1241 gmx_large_int_t step,
1242 t_inputrec *inputrec,
1246 gmx_wallcycle_t wcycle,
1250 gmx_bool bPCouple=FALSE;
1254 /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1255 if (inputrec->epc != epcNO && (!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))))
1257 /* We should only couple after a step where energies were determined */
1258 bPCouple = (inputrec->nstpcouple == 1 ||
1259 do_per_step(step+inputrec->nstpcouple-1,
1260 inputrec->nstpcouple));
1263 clear_mat(pcoupl_mu);
1264 for(i=0; i<DIM; i++)
1266 pcoupl_mu[i][i] = 1.0;
1273 dtpc = inputrec->nstpcouple*inputrec->delta_t;
1275 switch (inputrec->epc)
1277 /* We can always pcoupl, even if we did not sum the energies
1278 * the previous step, since state->pres_prev is only updated
1279 * when the energies have been summed.
1283 case (epcBERENDSEN):
1286 berendsen_pcoupl(fplog,step,inputrec,dtpc,state->pres_prev,state->box,
1290 case (epcPARRINELLORAHMAN):
1291 parrinellorahman_pcoupl(fplog,step,inputrec,dtpc,state->pres_prev,
1292 state->box,state->box_rel,state->boxv,
1293 M,pcoupl_mu,bInitStep);
1301 static rvec *get_xprime(const t_state *state,gmx_update_t upd)
1303 if (state->nalloc > upd->xp_nalloc)
1305 upd->xp_nalloc = state->nalloc;
1306 srenew(upd->xp,upd->xp_nalloc);
1312 void update_constraints(FILE *fplog,
1313 gmx_large_int_t step,
1314 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1315 t_inputrec *inputrec, /* input record and box stuff */
1316 gmx_ekindata_t *ekind,
1320 rvec force[], /* forces on home particles */
1323 tensor vir, /* tensors for virial and ekin, needed for computing */
1326 gmx_wallcycle_t wcycle,
1328 gmx_constr_t constr,
1330 gmx_bool bFirstHalf,
1334 gmx_bool bExtended,bLastStep,bLog=FALSE,bEner=FALSE,bDoConstr=FALSE;
1337 int start,homenr,nrend,i,n,m,g,d;
1339 rvec *vbuf,*xprime=NULL;
1341 if (constr) {bDoConstr=TRUE;}
1342 if (bFirstHalf && !EI_VV(inputrec->eI)) {bDoConstr=FALSE;}
1344 /* for now, SD update is here -- though it really seems like it
1345 should be reformulated as a velocity verlet method, since it has two parts */
1348 homenr = md->homenr;
1349 nrend = start+homenr;
1351 dt = inputrec->delta_t;
1356 * APPLY CONSTRAINTS:
1359 * When doing PR pressure coupling we have to constrain the
1360 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1361 * it is enough to do this once though, since the relative velocities
1362 * after this will be normal to the bond vector
1367 /* clear out constraints before applying */
1368 clear_mat(vir_part);
1370 xprime = get_xprime(state,upd);
1372 bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1373 bLog = (do_per_step(step,inputrec->nstlog) || bLastStep || (step < 0));
1374 bEner = (do_per_step(step,inputrec->nstenergy) || bLastStep);
1375 /* Constrain the coordinates xprime */
1376 wallcycle_start(wcycle,ewcCONSTR);
1377 if (EI_VV(inputrec->eI) && bFirstHalf)
1379 constrain(NULL,bLog,bEner,constr,idef,
1380 inputrec,ekind,cr,step,1,md,
1381 state->x,state->v,state->v,
1382 state->box,state->lambda[efptBONDED],dvdlambda,
1383 NULL,bCalcVir ? &vir_con : NULL,nrnb,econqVeloc,
1384 inputrec->epc==epcMTTK,state->veta,vetanew);
1388 constrain(NULL,bLog,bEner,constr,idef,
1389 inputrec,ekind,cr,step,1,md,
1390 state->x,xprime,NULL,
1391 state->box,state->lambda[efptBONDED],dvdlambda,
1392 state->v,bCalcVir ? &vir_con : NULL ,nrnb,econqCoord,
1393 inputrec->epc==epcMTTK,state->veta,state->veta);
1395 wallcycle_stop(wcycle,ewcCONSTR);
1399 dump_it_all(fplog,"After Shake",
1400 state->natoms,state->x,xprime,state->v,force);
1404 if (inputrec->eI == eiSD2)
1406 /* A correction factor eph is needed for the SD constraint force */
1407 /* Here we can, unfortunately, not have proper corrections
1408 * for different friction constants, so we use the first one.
1410 for(i=0; i<DIM; i++)
1412 for(m=0; m<DIM; m++)
1414 vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1420 m_add(vir_part,vir_con,vir_part);
1424 pr_rvecs(debug,0,"constraint virial",vir_part,DIM);
1430 if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1432 xprime = get_xprime(state,upd);
1434 /* The second part of the SD integration */
1435 do_update_sd2(upd->sd,FALSE,start,homenr,
1436 inputrec->opts.acc,inputrec->opts.nFreeze,
1437 md->invmass,md->ptype,
1438 md->cFREEZE,md->cACC,md->cTC,
1439 state->x,xprime,state->v,force,state->sd_X,
1440 inputrec->opts.ngtc,inputrec->opts.tau_t,
1441 inputrec->opts.ref_t,FALSE);
1442 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1446 /* Constrain the coordinates xprime */
1447 wallcycle_start(wcycle,ewcCONSTR);
1448 constrain(NULL,bLog,bEner,constr,idef,
1449 inputrec,NULL,cr,step,1,md,
1450 state->x,xprime,NULL,
1451 state->box,state->lambda[efptBONDED],dvdlambda,
1452 NULL,NULL,nrnb,econqCoord,FALSE,0,0);
1453 wallcycle_stop(wcycle,ewcCONSTR);
1457 /* We must always unshift after updating coordinates; if we did not shake
1458 x was shifted in do_force */
1460 if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1462 if (graph && (graph->nnodes > 0))
1464 unshift_x(graph,state->box,state->x,upd->xp);
1465 if (TRICLINIC(state->box))
1467 inc_nrnb(nrnb,eNR_SHIFTX,2*graph->nnodes);
1471 inc_nrnb(nrnb,eNR_SHIFTX,graph->nnodes);
1476 copy_rvecn(upd->xp,state->x,start,nrend);
1479 dump_it_all(fplog,"After unshift",
1480 state->natoms,state->x,upd->xp,state->v,force);
1482 /* ############# END the update of velocities and positions ######### */
1485 void update_box(FILE *fplog,
1486 gmx_large_int_t step,
1487 t_inputrec *inputrec, /* input record and box stuff */
1491 rvec force[], /* forces on home particles */
1495 gmx_wallcycle_t wcycle,
1498 gmx_bool bFirstHalf)
1500 gmx_bool bExtended,bLastStep,bLog=FALSE,bEner=FALSE;
1503 int start,homenr,nrend,i,n,m,g;
1507 homenr = md->homenr;
1508 nrend = start+homenr;
1511 (inputrec->etc == etcNOSEHOOVER) ||
1512 (inputrec->epc == epcPARRINELLORAHMAN) ||
1513 (inputrec->epc == epcMTTK);
1515 dt = inputrec->delta_t;
1519 /* now update boxes */
1520 switch (inputrec->epc) {
1523 case (epcBERENDSEN):
1524 berendsen_pscale(inputrec,pcoupl_mu,state->box,state->box_rel,
1525 start,homenr,state->x,md->cFREEZE,nrnb);
1527 case (epcPARRINELLORAHMAN):
1528 /* The box velocities were updated in do_pr_pcoupl in the update
1529 * iteration, but we dont change the box vectors until we get here
1530 * since we need to be able to shift/unshift above.
1536 state->box[i][m] += dt*state->boxv[i][m];
1539 preserve_box_shape(inputrec,state->box_rel,state->box);
1541 /* Scale the coordinates */
1542 for(n=start; (n<start+homenr); n++)
1544 tmvmul_ur0(pcoupl_mu,state->x[n],state->x[n]);
1548 switch (inputrec->epct)
1550 case (epctISOTROPIC):
1551 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1552 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1553 Side length scales as exp(veta*dt) */
1555 msmul(state->box,exp(state->veta*dt),state->box);
1557 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1558 o If we assume isotropic scaling, and box length scaling
1559 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1560 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1561 determinant of B is L^DIM det(M), and the determinant
1562 of dB/dt is (dL/dT)^DIM det (M). veta will be
1563 (det(dB/dT)/det(B))^(1/3). Then since M =
1564 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1566 msmul(state->box,state->veta,state->boxv);
1576 if ((!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))) && scale_tot)
1578 /* The transposes of the scaling matrices are stored,
1579 * therefore we need to reverse the order in the multiplication.
1581 mmul_ur0(*scale_tot,pcoupl_mu,*scale_tot);
1584 if (DEFORM(*inputrec))
1586 deform(upd,start,homenr,state->x,state->box,scale_tot,inputrec,step);
1589 dump_it_all(fplog,"After update",
1590 state->natoms,state->x,upd->xp,state->v,force);
1593 void update_coords(FILE *fplog,
1594 gmx_large_int_t step,
1595 t_inputrec *inputrec, /* input record and box stuff */
1598 rvec *f, /* forces on home particles */
1602 gmx_ekindata_t *ekind,
1604 gmx_wallcycle_t wcycle,
1608 t_commrec *cr, /* these shouldn't be here -- need to think about it */
1610 gmx_constr_t constr,
1613 gmx_bool bNH,bPR,bLastStep,bLog=FALSE,bEner=FALSE;
1615 real *imass,*imassin;
1618 int start,homenr,nrend,i,j,d,n,m,g;
1619 int blen0,blen1,iatom,jatom,nshake,nsettle,nconstr,nexpand;
1622 rvec *vcom,*xcom,*vall,*xall,*xin,*vin,*forcein,*fall,*xpall,*xprimein,*xprime;
1625 /* Running the velocity half does nothing except for velocity verlet */
1626 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1627 !EI_VV(inputrec->eI))
1629 gmx_incons("update_coords called for velocity without VV integrator");
1633 homenr = md->homenr;
1634 nrend = start+homenr;
1636 xprime = get_xprime(state,upd);
1638 dt = inputrec->delta_t;
1641 /* We need to update the NMR restraint history when time averaging is used */
1642 if (state->flags & (1<<estDISRE_RM3TAV))
1644 update_disres_history(fcd,&state->hist);
1646 if (state->flags & (1<<estORIRE_DTAV))
1648 update_orires_history(fcd,&state->hist);
1652 bNH = inputrec->etc == etcNOSEHOOVER;
1653 bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1655 if (bDoLR && inputrec->nstlist > 1 && !EI_VV(inputrec->eI)) /* get this working with VV? */
1657 /* Store the total force + nstlist-1 times the LR force
1658 * in forces_lr, so it can be used in a normal update algorithm
1659 * to produce twin time stepping.
1661 /* is this correct in the new construction? MRS */
1662 combine_forces(inputrec->nstlist,constr,inputrec,md,idef,cr,step,state,
1663 start,nrend,f,f_lr,nrnb);
1671 /* ############# START The update of velocities and positions ######### */
1673 dump_it_all(fplog,"Before update",
1674 state->natoms,state->x,xprime,state->v,force);
1676 switch (inputrec->eI) {
1678 if (ekind->cosacc.cos_accel == 0) {
1679 /* use normal version of update */
1680 do_update_md(start,nrend,dt,
1681 ekind->tcstat,ekind->grpstat,state->nosehoover_vxi,
1682 inputrec->opts.acc,inputrec->opts.nFreeze,md->invmass,md->ptype,
1683 md->cFREEZE,md->cACC,md->cTC,
1684 state->x,xprime,state->v,force,M,
1689 do_update_visc(start,nrend,dt,
1690 ekind->tcstat,md->invmass,state->nosehoover_vxi,
1691 md->ptype,md->cTC,state->x,xprime,state->v,force,M,
1692 state->box,ekind->cosacc.cos_accel,ekind->cosacc.vcos,bNH,bPR);
1696 do_update_sd1(upd->sd,start,homenr,dt,
1697 inputrec->opts.acc,inputrec->opts.nFreeze,
1698 md->invmass,md->ptype,
1699 md->cFREEZE,md->cACC,md->cTC,
1700 state->x,xprime,state->v,force,state->sd_X,
1701 inputrec->opts.ngtc,inputrec->opts.tau_t,inputrec->opts.ref_t);
1704 /* The SD update is done in 2 parts, because an extra constraint step
1707 do_update_sd2(upd->sd,bInitStep,start,homenr,
1708 inputrec->opts.acc,inputrec->opts.nFreeze,
1709 md->invmass,md->ptype,
1710 md->cFREEZE,md->cACC,md->cTC,
1711 state->x,xprime,state->v,force,state->sd_X,
1712 inputrec->opts.ngtc,inputrec->opts.tau_t,inputrec->opts.ref_t,
1716 do_update_bd(start,nrend,dt,
1717 inputrec->opts.nFreeze,md->invmass,md->ptype,
1718 md->cFREEZE,md->cTC,
1719 state->x,xprime,state->v,force,
1721 inputrec->opts.ngtc,inputrec->opts.tau_t,inputrec->opts.ref_t,
1722 upd->sd->bd_rf,upd->sd->gaussrand);
1726 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
1727 switch (UpdatePart) {
1730 do_update_vv_vel(start,nrend,dt,
1731 ekind->tcstat,ekind->grpstat,
1732 inputrec->opts.acc,inputrec->opts.nFreeze,
1733 md->invmass,md->ptype,md->cFREEZE,md->cACC,
1734 state->v,force,(bNH || bPR),state->veta,alpha);
1737 do_update_vv_pos(start,nrend,dt,
1738 ekind->tcstat,ekind->grpstat,
1739 inputrec->opts.acc,inputrec->opts.nFreeze,
1740 md->invmass,md->ptype,md->cFREEZE,
1741 state->x,xprime,state->v,force,
1742 (bNH || bPR) ,state->veta,alpha);
1747 gmx_fatal(FARGS,"Don't know how to update coordinates");
1754 void correct_ekin(FILE *log,int start,int end,rvec v[],rvec vcm,real mass[],
1755 real tmass,tensor ekin)
1758 * This is a debugging routine. It should not be called for production code
1760 * The kinetic energy should calculated according to:
1761 * Ekin = 1/2 m (v-vcm)^2
1762 * However the correction is not always applied, since vcm may not be
1763 * known in time and we compute
1764 * Ekin' = 1/2 m v^2 instead
1765 * This can be corrected afterwards by computing
1766 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
1768 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
1775 /* Local particles */
1778 /* Processor dependent part. */
1780 for(i=start; (i<end); i++)
1784 for(j=0; (j<DIM); j++)
1790 svmul(1/tmass,vcm,vcm);
1791 svmul(0.5,vcm,hvcm);
1793 for(j=0; (j<DIM); j++)
1795 for(k=0; (k<DIM); k++)
1797 dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
1800 pr_rvecs(log,0,"dekin",dekin,DIM);
1801 pr_rvecs(log,0," ekin", ekin,DIM);
1802 fprintf(log,"dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
1803 trace(dekin),trace(ekin),vcm[XX],vcm[YY],vcm[ZZ]);
1804 fprintf(log,"mv = (%8.4f %8.4f %8.4f)\n",
1805 mv[XX],mv[YY],mv[ZZ]);
1808 extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_large_int_t step, t_mdatoms *md, t_state *state, gmx_update_t upd, t_idef *idef, gmx_constr_t constr) {
1811 real rate = (ir->delta_t)/ir->opts.tau_t[0];
1812 /* proceed with andersen if 1) it's fixed probability per
1813 particle andersen or 2) it's massive andersen and it's tau_t/dt */
1814 if ((ir->etc==etcANDERSEN) || do_per_step(step,(int)(1.0/rate)))
1816 srenew(upd->randatom,state->nalloc);
1817 srenew(upd->randatom_list,state->nalloc);
1818 if (upd->randatom_list_init == FALSE) {
1819 for (i=0;i<state->nalloc;i++) {
1820 upd->randatom[i] = FALSE;
1821 upd->randatom_list[i] = 0;
1823 upd->randatom_list_init = TRUE;
1825 andersen_tcoupl(ir,md,state,upd->sd->gaussrand,rate,
1826 (ir->etc==etcANDERSEN)?idef:NULL,
1827 constr?get_nblocks(constr):0,
1828 constr?get_sblock(constr):NULL,
1829 upd->randatom,upd->randatom_list,
1830 upd->sd->randomize_group,upd->sd->boltzfac);