2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "types/commrec.h"
47 #include "gromacs/utility/smalloc.h"
55 #include "gromacs/random/random.h"
65 #include "gmx_omp_nthreads.h"
67 #include "gromacs/fileio/confio.h"
68 #include "gromacs/fileio/futil.h"
69 #include "gromacs/timing/wallcycle.h"
70 #include "gromacs/utility/gmxomp.h"
71 #include "gromacs/pulling/pull.h"
73 /*For debugging, start at v(-dt/2) for velolcity verlet -- uncomment next line */
74 /*#define STARTFROMDT2*/
98 gmx_sd_sigma_t *sdsig;
101 /* andersen temperature control stuff */
102 gmx_bool *randomize_group;
106 typedef struct gmx_update
109 /* xprime for constraint algorithms */
113 /* Variables for the deform algorithm */
114 gmx_int64_t deformref_step;
115 matrix deformref_box;
119 static void do_update_md(int start, int nrend, double dt,
120 t_grp_tcstat *tcstat,
122 gmx_bool bNEMD, t_grp_acc *gstat, rvec accel[],
125 unsigned short ptype[], unsigned short cFREEZE[],
126 unsigned short cACC[], unsigned short cTC[],
127 rvec x[], rvec xprime[], rvec v[],
129 gmx_bool bNH, gmx_bool bPR)
132 int gf = 0, ga = 0, gt = 0;
134 real vn, vv, va, vb, vnrel;
140 /* Update with coupling to extended ensembles, used for
141 * Nose-Hoover and Parrinello-Rahman coupling
142 * Nose-Hoover uses the reversible leap-frog integrator from
143 * Holian et al. Phys Rev E 52(3) : 2338, 1995
145 for (n = start; n < nrend; n++)
160 lg = tcstat[gt].lambda;
165 rvec_sub(v[n], gstat[ga].u, vrel);
167 for (d = 0; d < DIM; d++)
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];
186 else if (cFREEZE != NULL ||
187 nFreeze[0][XX] || nFreeze[0][YY] || nFreeze[0][ZZ] ||
190 /* Update with Berendsen/v-rescale coupling and freeze or NEMD */
191 for (n = start; n < nrend; n++)
193 w_dt = invmass[n]*dt;
206 lg = tcstat[gt].lambda;
208 for (d = 0; d < DIM; d++)
211 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
213 vv = lg*vn + f[n][d]*w_dt;
215 /* do not scale the mean velocities u */
217 va = vv + accel[ga][d]*dt;
218 vb = va + (1.0-lg)*u;
220 xprime[n][d] = x[n][d]+vb*dt;
225 xprime[n][d] = x[n][d];
232 /* Plain update with Berendsen/v-rescale coupling */
233 for (n = start; n < nrend; n++)
235 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
237 w_dt = invmass[n]*dt;
242 lg = tcstat[gt].lambda;
244 for (d = 0; d < DIM; d++)
246 vn = lg*v[n][d] + f[n][d]*w_dt;
248 xprime[n][d] = x[n][d] + vn*dt;
253 for (d = 0; d < DIM; d++)
256 xprime[n][d] = x[n][d];
263 static void do_update_vv_vel(int start, int nrend, double dt,
264 rvec accel[], ivec nFreeze[], real invmass[],
265 unsigned short ptype[], unsigned short cFREEZE[],
266 unsigned short cACC[], rvec v[], rvec f[],
267 gmx_bool bExtended, real veta, real alpha)
272 real u, vn, vv, va, vb, vnrel;
278 g = 0.25*dt*veta*alpha;
280 mv2 = series_sinhx(g);
287 for (n = start; n < nrend; n++)
289 w_dt = invmass[n]*dt;
299 for (d = 0; d < DIM; d++)
301 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
303 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
311 } /* do_update_vv_vel */
313 static void do_update_vv_pos(int start, int nrend, double dt,
315 unsigned short ptype[], unsigned short cFREEZE[],
316 rvec x[], rvec xprime[], rvec v[],
317 gmx_bool bExtended, real veta)
324 /* Would it make more sense if Parrinello-Rahman was put here? */
329 mr2 = series_sinhx(g);
337 for (n = start; n < nrend; n++)
345 for (d = 0; d < DIM; d++)
347 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
349 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
353 xprime[n][d] = x[n][d];
357 } /* do_update_vv_pos */
359 static void do_update_visc(int start, int nrend, double dt,
360 t_grp_tcstat *tcstat,
363 unsigned short ptype[], unsigned short cTC[],
364 rvec x[], rvec xprime[], rvec v[],
365 rvec f[], matrix M, matrix box, real
366 cos_accel, real vcos,
367 gmx_bool bNH, gmx_bool bPR)
372 real lg, vxi = 0, vv;
377 fac = 2*M_PI/(box[ZZ][ZZ]);
381 /* Update with coupling to extended ensembles, used for
382 * Nose-Hoover and Parrinello-Rahman coupling
384 for (n = start; n < nrend; n++)
391 lg = tcstat[gt].lambda;
392 cosz = cos(fac*x[n][ZZ]);
394 copy_rvec(v[n], vrel);
402 for (d = 0; d < DIM; d++)
406 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
408 vn = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
409 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
412 vn += vc + dt*cosz*cos_accel;
415 xprime[n][d] = x[n][d]+vn*dt;
419 xprime[n][d] = x[n][d];
426 /* Classic version of update, used with berendsen coupling */
427 for (n = start; n < nrend; n++)
429 w_dt = invmass[n]*dt;
434 lg = tcstat[gt].lambda;
435 cosz = cos(fac*x[n][ZZ]);
437 for (d = 0; d < DIM; d++)
441 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
446 /* Do not scale the cosine velocity profile */
447 vv = vc + lg*(vn - vc + f[n][d]*w_dt);
448 /* Add the cosine accelaration profile */
449 vv += dt*cosz*cos_accel;
453 vv = lg*(vn + f[n][d]*w_dt);
456 xprime[n][d] = x[n][d]+vv*dt;
461 xprime[n][d] = x[n][d];
468 static gmx_stochd_t *init_stochd(t_inputrec *ir)
477 ngtc = ir->opts.ngtc;
481 snew(sd->bd_rf, ngtc);
483 else if (EI_SD(ir->eI))
486 snew(sd->sdsig, ngtc);
489 for (n = 0; n < ngtc; n++)
491 if (ir->opts.tau_t[n] > 0)
493 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
494 sdc[n].eph = exp(sdc[n].gdt/2);
495 sdc[n].emh = exp(-sdc[n].gdt/2);
496 sdc[n].em = exp(-sdc[n].gdt);
500 /* No friction and noise on this group */
506 if (sdc[n].gdt >= 0.05)
508 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
509 - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
510 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
511 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
516 /* Seventh order expansions for small y */
517 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
518 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))));
519 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
523 fprintf(debug, "SD const tc-grp %d: b %g c %g d %g\n",
524 n, sdc[n].b, sdc[n].c, sdc[n].d);
528 else if (ETC_ANDERSEN(ir->etc))
537 snew(sd->randomize_group, ngtc);
538 snew(sd->boltzfac, ngtc);
540 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
541 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
543 for (n = 0; n < ngtc; n++)
545 reft = max(0.0, opts->ref_t[n]);
546 if ((opts->tau_t[n] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
548 sd->randomize_group[n] = TRUE;
549 sd->boltzfac[n] = BOLTZ*opts->ref_t[n];
553 sd->randomize_group[n] = FALSE;
560 gmx_update_t init_update(t_inputrec *ir)
566 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
568 upd->sd = init_stochd(ir);
577 static void do_update_sd1(gmx_stochd_t *sd,
578 int start, int nrend, double dt,
579 rvec accel[], ivec nFreeze[],
580 real invmass[], unsigned short ptype[],
581 unsigned short cFREEZE[], unsigned short cACC[],
582 unsigned short cTC[],
583 rvec x[], rvec xprime[], rvec v[], rvec f[],
584 int ngtc, real tau_t[], real ref_t[],
585 gmx_int64_t step, int seed, int* gatindex)
590 int gf = 0, ga = 0, gt = 0;
597 for (n = 0; n < ngtc; n++)
600 /* The mass is encounted for later, since this differs per atom */
601 sig[n].V = sqrt(kT*(1 - sdc[n].em*sdc[n].em));
604 for (n = start; n < nrend; n++)
607 int ng = gatindex ? gatindex[n] : n;
608 ism = sqrt(invmass[n]);
622 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
623 for (d = 0; d < DIM; d++)
625 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
627 sd_V = ism*sig[gt].V*rnd[d];
629 v[n][d] = v[n][d]*sdc[gt].em
630 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
633 xprime[n][d] = x[n][d] + v[n][d]*dt;
638 xprime[n][d] = x[n][d];
644 static void check_sd2_work_data_allocation(gmx_stochd_t *sd, int nrend)
646 if (nrend > sd->sd_V_nalloc)
648 sd->sd_V_nalloc = over_alloc_dd(nrend);
649 srenew(sd->sd_V, sd->sd_V_nalloc);
653 static void do_update_sd2_Tconsts(gmx_stochd_t *sd,
658 /* This is separated from the update below, because it is single threaded */
667 for (gt = 0; gt < ngtc; gt++)
669 kT = BOLTZ*ref_t[gt];
670 /* The mass is encounted for later, since this differs per atom */
671 sig[gt].V = sqrt(kT*(1-sdc[gt].em));
672 sig[gt].X = sqrt(kT*sqr(tau_t[gt])*sdc[gt].c);
673 sig[gt].Yv = sqrt(kT*sdc[gt].b/sdc[gt].c);
674 sig[gt].Yx = sqrt(kT*sqr(tau_t[gt])*sdc[gt].b/(1-sdc[gt].em));
678 static void do_update_sd2(gmx_stochd_t *sd,
680 int start, int nrend,
681 rvec accel[], ivec nFreeze[],
682 real invmass[], unsigned short ptype[],
683 unsigned short cFREEZE[], unsigned short cACC[],
684 unsigned short cTC[],
685 rvec x[], rvec xprime[], rvec v[], rvec f[],
688 gmx_bool bFirstHalf, gmx_int64_t step, int seed,
693 /* The random part of the velocity update, generated in the first
694 * half of the update, needs to be remembered for the second half.
698 int gf = 0, ga = 0, gt = 0;
699 real vn = 0, Vmh, Xmh;
707 for (n = start; n < nrend; n++)
709 real rnd[6], rndi[3];
710 ng = gatindex ? gatindex[n] : n;
711 ism = sqrt(invmass[n]);
725 gmx_rng_cycle_6gaussian_table(step*2+(bFirstHalf ? 1 : 2), ng, seed, RND_SEED_UPDATE, rnd);
728 gmx_rng_cycle_3gaussian_table(step*2, ng, seed, RND_SEED_UPDATE, rndi);
730 for (d = 0; d < DIM; d++)
736 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
742 sd_X[n][d] = ism*sig[gt].X*rndi[d];
744 Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
745 + ism*sig[gt].Yv*rnd[d*2];
746 sd_V[n][d] = ism*sig[gt].V*rnd[d*2+1];
748 v[n][d] = vn*sdc[gt].em
749 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
750 + sd_V[n][d] - sdc[gt].em*Vmh;
752 xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
756 /* Correct the velocities for the constraints.
757 * This operation introduces some inaccuracy,
758 * since the velocity is determined from differences in coordinates.
761 (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
763 Xmh = sd_V[n][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
764 + ism*sig[gt].Yx*rnd[d*2];
765 sd_X[n][d] = ism*sig[gt].X*rnd[d*2+1];
767 xprime[n][d] += sd_X[n][d] - Xmh;
776 xprime[n][d] = x[n][d];
783 static void do_update_bd_Tconsts(double dt, real friction_coefficient,
784 int ngtc, const real ref_t[],
787 /* This is separated from the update below, because it is single threaded */
790 if (friction_coefficient != 0)
792 for (gt = 0; gt < ngtc; gt++)
794 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]/(friction_coefficient*dt));
799 for (gt = 0; gt < ngtc; gt++)
801 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]);
806 static void do_update_bd(int start, int nrend, double dt,
808 real invmass[], unsigned short ptype[],
809 unsigned short cFREEZE[], unsigned short cTC[],
810 rvec x[], rvec xprime[], rvec v[],
811 rvec f[], real friction_coefficient,
812 real *rf, gmx_int64_t step, int seed,
815 /* note -- these appear to be full step velocities . . . */
821 if (friction_coefficient != 0)
823 invfr = 1.0/friction_coefficient;
826 for (n = start; (n < nrend); n++)
829 int ng = gatindex ? gatindex[n] : n;
839 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
840 for (d = 0; (d < DIM); d++)
842 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
844 if (friction_coefficient != 0)
846 vn = invfr*f[n][d] + rf[gt]*rnd[d];
850 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
851 vn = 0.5*invmass[n]*f[n][d]*dt
852 + sqrt(0.5*invmass[n])*rf[gt]*rnd[d];
856 xprime[n][d] = x[n][d]+vn*dt;
861 xprime[n][d] = x[n][d];
867 static void dump_it_all(FILE gmx_unused *fp, const char gmx_unused *title,
868 int gmx_unused natoms, rvec gmx_unused x[], rvec gmx_unused xp[],
869 rvec gmx_unused v[], rvec gmx_unused f[])
874 fprintf(fp, "%s\n", title);
875 pr_rvecs(fp, 0, "x", x, natoms);
876 pr_rvecs(fp, 0, "xp", xp, natoms);
877 pr_rvecs(fp, 0, "v", v, natoms);
878 pr_rvecs(fp, 0, "f", f, natoms);
883 static void calc_ke_part_normal(rvec v[], t_grpopts *opts, t_mdatoms *md,
884 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel,
885 gmx_bool bSaveEkinOld)
888 t_grp_tcstat *tcstat = ekind->tcstat;
889 t_grp_acc *grpstat = ekind->grpstat;
892 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
893 an option, but not supported now. Additionally, if we are doing iterations.
894 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
895 bSavEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't copy over the ekinh_old.
896 If FALSE, we overrwrite it.
899 /* group velocities are calculated in update_ekindata and
900 * accumulated in acumulate_groups.
901 * Now the partial global and groups ekin.
903 for (g = 0; (g < opts->ngtc); g++)
908 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
912 clear_mat(tcstat[g].ekinf);
916 clear_mat(tcstat[g].ekinh);
920 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
923 ekind->dekindl_old = ekind->dekindl;
925 nthread = gmx_omp_nthreads_get(emntUpdate);
927 #pragma omp parallel for num_threads(nthread) schedule(static)
928 for (thread = 0; thread < nthread; thread++)
930 int start_t, end_t, n;
938 start_t = ((thread+0)*md->homenr)/nthread;
939 end_t = ((thread+1)*md->homenr)/nthread;
941 ekin_sum = ekind->ekin_work[thread];
942 dekindl_sum = ekind->dekindl_work[thread];
944 for (gt = 0; gt < opts->ngtc; gt++)
946 clear_mat(ekin_sum[gt]);
952 for (n = start_t; n < end_t; n++)
962 hm = 0.5*md->massT[n];
964 for (d = 0; (d < DIM); d++)
966 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
968 for (d = 0; (d < DIM); d++)
970 for (m = 0; (m < DIM); m++)
972 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
973 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
976 if (md->nMassPerturbed && md->bPerturbed[n])
979 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
985 for (thread = 0; thread < nthread; thread++)
987 for (g = 0; g < opts->ngtc; g++)
991 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
996 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1001 ekind->dekindl += *ekind->dekindl_work[thread];
1004 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1007 static void calc_ke_part_visc(matrix box, rvec x[], rvec v[],
1008 t_grpopts *opts, t_mdatoms *md,
1009 gmx_ekindata_t *ekind,
1010 t_nrnb *nrnb, gmx_bool bEkinAveVel)
1012 int start = 0, homenr = md->homenr;
1013 int g, d, n, m, gt = 0;
1016 t_grp_tcstat *tcstat = ekind->tcstat;
1017 t_cos_acc *cosacc = &(ekind->cosacc);
1022 for (g = 0; g < opts->ngtc; g++)
1024 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1025 clear_mat(ekind->tcstat[g].ekinh);
1027 ekind->dekindl_old = ekind->dekindl;
1029 fac = 2*M_PI/box[ZZ][ZZ];
1032 for (n = start; n < start+homenr; n++)
1038 hm = 0.5*md->massT[n];
1040 /* Note that the times of x and v differ by half a step */
1041 /* MRS -- would have to be changed for VV */
1042 cosz = cos(fac*x[n][ZZ]);
1043 /* Calculate the amplitude of the new velocity profile */
1044 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1046 copy_rvec(v[n], v_corrt);
1047 /* Subtract the profile for the kinetic energy */
1048 v_corrt[XX] -= cosz*cosacc->vcos;
1049 for (d = 0; (d < DIM); d++)
1051 for (m = 0; (m < DIM); m++)
1053 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1056 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1060 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1064 if (md->nPerturbed && md->bPerturbed[n])
1066 dekindl += 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1069 ekind->dekindl = dekindl;
1070 cosacc->mvcos = mvcos;
1072 inc_nrnb(nrnb, eNR_EKIN, homenr);
1075 void calc_ke_part(t_state *state, t_grpopts *opts, t_mdatoms *md,
1076 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
1078 if (ekind->cosacc.cos_accel == 0)
1080 calc_ke_part_normal(state->v, opts, md, ekind, nrnb, bEkinAveVel, bSaveEkinOld);
1084 calc_ke_part_visc(state->box, state->x, state->v, opts, md, ekind, nrnb, bEkinAveVel);
1088 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1090 ekinstate->ekin_n = ir->opts.ngtc;
1091 snew(ekinstate->ekinh, ekinstate->ekin_n);
1092 snew(ekinstate->ekinf, ekinstate->ekin_n);
1093 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1094 snew(ekinstate->ekinscalef_nhc, ekinstate->ekin_n);
1095 snew(ekinstate->ekinscaleh_nhc, ekinstate->ekin_n);
1096 snew(ekinstate->vscale_nhc, ekinstate->ekin_n);
1097 ekinstate->dekindl = 0;
1098 ekinstate->mvcos = 0;
1101 void update_ekinstate(ekinstate_t *ekinstate, gmx_ekindata_t *ekind)
1105 for (i = 0; i < ekinstate->ekin_n; i++)
1107 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1108 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1109 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1110 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1111 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1112 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1115 copy_mat(ekind->ekin, ekinstate->ekin_total);
1116 ekinstate->dekindl = ekind->dekindl;
1117 ekinstate->mvcos = ekind->cosacc.mvcos;
1121 void restore_ekinstate_from_state(t_commrec *cr,
1122 gmx_ekindata_t *ekind, ekinstate_t *ekinstate)
1128 for (i = 0; i < ekinstate->ekin_n; i++)
1130 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1131 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1132 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1133 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1134 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1135 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1138 copy_mat(ekinstate->ekin_total, ekind->ekin);
1140 ekind->dekindl = ekinstate->dekindl;
1141 ekind->cosacc.mvcos = ekinstate->mvcos;
1142 n = ekinstate->ekin_n;
1147 gmx_bcast(sizeof(n), &n, cr);
1148 for (i = 0; i < n; i++)
1150 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1151 ekind->tcstat[i].ekinh[0], cr);
1152 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1153 ekind->tcstat[i].ekinf[0], cr);
1154 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1155 ekind->tcstat[i].ekinh_old[0], cr);
1157 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1158 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1159 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1160 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1161 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1162 &(ekind->tcstat[i].vscale_nhc), cr);
1164 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1165 ekind->ekin[0], cr);
1167 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1168 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1172 void set_deform_reference_box(gmx_update_t upd, gmx_int64_t step, matrix box)
1174 upd->deformref_step = step;
1175 copy_mat(box, upd->deformref_box);
1178 static void deform(gmx_update_t upd,
1179 int start, int homenr, rvec x[], matrix box, matrix *scale_tot,
1180 const t_inputrec *ir, gmx_int64_t step)
1182 matrix bnew, invbox, mu;
1186 elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1187 copy_mat(box, bnew);
1188 for (i = 0; i < DIM; i++)
1190 for (j = 0; j < DIM; j++)
1192 if (ir->deform[i][j] != 0)
1195 upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1199 /* We correct the off-diagonal elements,
1200 * which can grow indefinitely during shearing,
1201 * so the shifts do not get messed up.
1203 for (i = 1; i < DIM; i++)
1205 for (j = i-1; j >= 0; j--)
1207 while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1209 rvec_dec(bnew[i], bnew[j]);
1211 while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1213 rvec_inc(bnew[i], bnew[j]);
1217 m_inv_ur0(box, invbox);
1218 copy_mat(bnew, box);
1219 mmul_ur0(box, invbox, mu);
1221 for (i = start; i < start+homenr; i++)
1223 x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1224 x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1225 x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1229 /* The transposes of the scaling matrices are stored,
1230 * so we need to do matrix multiplication in the inverse order.
1232 mmul_ur0(*scale_tot, mu, *scale_tot);
1236 static void combine_forces(int nstcalclr,
1237 gmx_constr_t constr,
1238 t_inputrec *ir, t_mdatoms *md, t_idef *idef,
1241 t_state *state, gmx_bool bMolPBC,
1242 int start, int nrend,
1243 rvec f[], rvec f_lr[],
1248 /* f contains the short-range forces + the long range forces
1249 * which are stored separately in f_lr.
1252 if (constr != NULL && !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1254 /* We need to constrain the LR forces separately,
1255 * because due to the different pre-factor for the SR and LR
1256 * forces in the update algorithm, we can not determine
1257 * the constraint force for the coordinate constraining.
1258 * Constrain only the additional LR part of the force.
1260 /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1261 constrain(NULL, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
1262 state->x, f_lr, f_lr, bMolPBC, state->box, state->lambda[efptBONDED], NULL,
1263 NULL, NULL, nrnb, econqForce, ir->epc == epcMTTK, state->veta, state->veta);
1266 /* Add nstcalclr-1 times the LR force to the sum of both forces
1267 * and store the result in forces_lr.
1269 nm1 = nstcalclr - 1;
1270 for (i = start; i < nrend; i++)
1272 for (d = 0; d < DIM; d++)
1274 f_lr[i][d] = f[i][d] + nm1*f_lr[i][d];
1279 void update_tcouple(gmx_int64_t step,
1280 t_inputrec *inputrec,
1282 gmx_ekindata_t *ekind,
1287 gmx_bool bTCouple = FALSE;
1289 int i, start, end, homenr, offset;
1291 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1292 if (inputrec->etc != etcNO &&
1293 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec)))
1295 /* We should only couple after a step where energies were determined (for leapfrog versions)
1296 or the step energies are determined, for velocity verlet versions */
1298 if (EI_VV(inputrec->eI))
1306 bTCouple = (inputrec->nsttcouple == 1 ||
1307 do_per_step(step+inputrec->nsttcouple-offset,
1308 inputrec->nsttcouple));
1313 dttc = inputrec->nsttcouple*inputrec->delta_t;
1315 switch (inputrec->etc)
1320 berendsen_tcoupl(inputrec, ekind, dttc);
1323 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1324 state->nosehoover_xi, state->nosehoover_vxi, MassQ);
1327 vrescale_tcoupl(inputrec, step, ekind, dttc,
1328 state->therm_integral);
1331 /* rescale in place here */
1332 if (EI_VV(inputrec->eI))
1334 rescale_velocities(ekind, md, 0, md->homenr, state->v);
1339 /* Set the T scaling lambda to 1 to have no scaling */
1340 for (i = 0; (i < inputrec->opts.ngtc); i++)
1342 ekind->tcstat[i].lambda = 1.0;
1347 void update_pcouple(FILE *fplog,
1349 t_inputrec *inputrec,
1355 gmx_bool bPCouple = FALSE;
1359 /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1360 if (inputrec->epc != epcNO && (!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))))
1362 /* We should only couple after a step where energies were determined */
1363 bPCouple = (inputrec->nstpcouple == 1 ||
1364 do_per_step(step+inputrec->nstpcouple-1,
1365 inputrec->nstpcouple));
1368 clear_mat(pcoupl_mu);
1369 for (i = 0; i < DIM; i++)
1371 pcoupl_mu[i][i] = 1.0;
1378 dtpc = inputrec->nstpcouple*inputrec->delta_t;
1380 switch (inputrec->epc)
1382 /* We can always pcoupl, even if we did not sum the energies
1383 * the previous step, since state->pres_prev is only updated
1384 * when the energies have been summed.
1388 case (epcBERENDSEN):
1391 berendsen_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev, state->box,
1395 case (epcPARRINELLORAHMAN):
1396 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1397 state->box, state->box_rel, state->boxv,
1398 M, pcoupl_mu, bInitStep);
1406 static rvec *get_xprime(const t_state *state, gmx_update_t upd)
1408 if (state->nalloc > upd->xp_nalloc)
1410 upd->xp_nalloc = state->nalloc;
1411 srenew(upd->xp, upd->xp_nalloc);
1417 void update_constraints(FILE *fplog,
1419 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1420 t_inputrec *inputrec, /* input record and box stuff */
1421 gmx_ekindata_t *ekind,
1426 rvec force[], /* forces on home particles */
1431 gmx_wallcycle_t wcycle,
1433 gmx_constr_t constr,
1434 gmx_bool bFirstHalf,
1438 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1441 int start, homenr, nrend, i, n, m, g, d;
1443 rvec *vbuf, *xprime = NULL;
1450 if (bFirstHalf && !EI_VV(inputrec->eI))
1455 /* for now, SD update is here -- though it really seems like it
1456 should be reformulated as a velocity verlet method, since it has two parts */
1459 homenr = md->homenr;
1460 nrend = start+homenr;
1462 dt = inputrec->delta_t;
1467 * APPLY CONSTRAINTS:
1470 * When doing PR pressure coupling we have to constrain the
1471 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1472 * it is enough to do this once though, since the relative velocities
1473 * after this will be normal to the bond vector
1478 /* clear out constraints before applying */
1479 clear_mat(vir_part);
1481 xprime = get_xprime(state, upd);
1483 bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1484 bLog = (do_per_step(step, inputrec->nstlog) || bLastStep || (step < 0));
1485 bEner = (do_per_step(step, inputrec->nstenergy) || bLastStep);
1486 /* Constrain the coordinates xprime */
1487 wallcycle_start(wcycle, ewcCONSTR);
1488 if (EI_VV(inputrec->eI) && bFirstHalf)
1490 constrain(NULL, bLog, bEner, constr, idef,
1491 inputrec, ekind, cr, step, 1, md,
1492 state->x, state->v, state->v,
1493 bMolPBC, state->box,
1494 state->lambda[efptBONDED], dvdlambda,
1495 NULL, bCalcVir ? &vir_con : NULL, nrnb, econqVeloc,
1496 inputrec->epc == epcMTTK, state->veta, vetanew);
1500 constrain(NULL, bLog, bEner, constr, idef,
1501 inputrec, ekind, cr, step, 1, md,
1502 state->x, xprime, NULL,
1503 bMolPBC, state->box,
1504 state->lambda[efptBONDED], dvdlambda,
1505 state->v, bCalcVir ? &vir_con : NULL, nrnb, econqCoord,
1506 inputrec->epc == epcMTTK, state->veta, state->veta);
1508 wallcycle_stop(wcycle, ewcCONSTR);
1512 dump_it_all(fplog, "After Shake",
1513 state->natoms, state->x, xprime, state->v, force);
1517 if (inputrec->eI == eiSD2)
1519 /* A correction factor eph is needed for the SD constraint force */
1520 /* Here we can, unfortunately, not have proper corrections
1521 * for different friction constants, so we use the first one.
1523 for (i = 0; i < DIM; i++)
1525 for (m = 0; m < DIM; m++)
1527 vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1533 m_add(vir_part, vir_con, vir_part);
1537 pr_rvecs(debug, 0, "constraint virial", vir_part, DIM);
1543 if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1545 xprime = get_xprime(state, upd);
1547 nth = gmx_omp_nthreads_get(emntUpdate);
1549 #pragma omp parallel for num_threads(nth) schedule(static)
1550 for (th = 0; th < nth; th++)
1552 int start_th, end_th;
1554 start_th = start + ((nrend-start)* th )/nth;
1555 end_th = start + ((nrend-start)*(th+1))/nth;
1557 /* The second part of the SD integration */
1558 do_update_sd2(upd->sd,
1559 FALSE, start_th, end_th,
1560 inputrec->opts.acc, inputrec->opts.nFreeze,
1561 md->invmass, md->ptype,
1562 md->cFREEZE, md->cACC, md->cTC,
1563 state->x, xprime, state->v, force, state->sd_X,
1564 inputrec->opts.tau_t,
1565 FALSE, step, inputrec->ld_seed,
1566 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1568 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1572 /* Constrain the coordinates xprime */
1573 wallcycle_start(wcycle, ewcCONSTR);
1574 constrain(NULL, bLog, bEner, constr, idef,
1575 inputrec, NULL, cr, step, 1, md,
1576 state->x, xprime, NULL,
1577 bMolPBC, state->box,
1578 state->lambda[efptBONDED], dvdlambda,
1579 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
1580 wallcycle_stop(wcycle, ewcCONSTR);
1584 /* We must always unshift after updating coordinates; if we did not shake
1585 x was shifted in do_force */
1587 if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1589 if (graph && (graph->nnodes > 0))
1591 unshift_x(graph, state->box, state->x, upd->xp);
1592 if (TRICLINIC(state->box))
1594 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1598 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1603 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntUpdate)) schedule(static)
1604 for (i = start; i < nrend; i++)
1606 copy_rvec(upd->xp[i], state->x[i]);
1610 dump_it_all(fplog, "After unshift",
1611 state->natoms, state->x, upd->xp, state->v, force);
1613 /* ############# END the update of velocities and positions ######### */
1616 void update_box(FILE *fplog,
1618 t_inputrec *inputrec, /* input record and box stuff */
1621 rvec force[], /* forces on home particles */
1627 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE;
1630 int start, homenr, nrend, i, n, m, g;
1634 homenr = md->homenr;
1635 nrend = start+homenr;
1638 (inputrec->etc == etcNOSEHOOVER) ||
1639 (inputrec->epc == epcPARRINELLORAHMAN) ||
1640 (inputrec->epc == epcMTTK);
1642 dt = inputrec->delta_t;
1646 /* now update boxes */
1647 switch (inputrec->epc)
1651 case (epcBERENDSEN):
1652 berendsen_pscale(inputrec, pcoupl_mu, state->box, state->box_rel,
1653 start, homenr, state->x, md->cFREEZE, nrnb);
1655 case (epcPARRINELLORAHMAN):
1656 /* The box velocities were updated in do_pr_pcoupl in the update
1657 * iteration, but we dont change the box vectors until we get here
1658 * since we need to be able to shift/unshift above.
1660 for (i = 0; i < DIM; i++)
1662 for (m = 0; m <= i; m++)
1664 state->box[i][m] += dt*state->boxv[i][m];
1667 preserve_box_shape(inputrec, state->box_rel, state->box);
1669 /* Scale the coordinates */
1670 for (n = start; (n < start+homenr); n++)
1672 tmvmul_ur0(pcoupl_mu, state->x[n], state->x[n]);
1676 switch (inputrec->epct)
1678 case (epctISOTROPIC):
1679 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1680 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1681 Side length scales as exp(veta*dt) */
1683 msmul(state->box, exp(state->veta*dt), state->box);
1685 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1686 o If we assume isotropic scaling, and box length scaling
1687 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1688 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1689 determinant of B is L^DIM det(M), and the determinant
1690 of dB/dt is (dL/dT)^DIM det (M). veta will be
1691 (det(dB/dT)/det(B))^(1/3). Then since M =
1692 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1694 msmul(state->box, state->veta, state->boxv);
1704 if ((!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))) && scale_tot)
1706 /* The transposes of the scaling matrices are stored,
1707 * therefore we need to reverse the order in the multiplication.
1709 mmul_ur0(*scale_tot, pcoupl_mu, *scale_tot);
1712 if (DEFORM(*inputrec))
1714 deform(upd, start, homenr, state->x, state->box, scale_tot, inputrec, step);
1717 dump_it_all(fplog, "After update",
1718 state->natoms, state->x, upd->xp, state->v, force);
1721 void update_coords(FILE *fplog,
1723 t_inputrec *inputrec, /* input record and box stuff */
1727 rvec *f, /* forces on home particles */
1731 gmx_ekindata_t *ekind,
1736 t_commrec *cr, /* these shouldn't be here -- need to think about it */
1738 gmx_constr_t constr,
1741 gmx_bool bNH, bPR, bLastStep, bLog = FALSE, bEner = FALSE;
1743 real *imass, *imassin;
1746 int start, homenr, nrend, i, j, d, n, m, g;
1747 int blen0, blen1, iatom, jatom, nshake, nsettle, nconstr, nexpand;
1750 rvec *vcom, *xcom, *vall, *xall, *xin, *vin, *forcein, *fall, *xpall, *xprimein, *xprime;
1753 /* Running the velocity half does nothing except for velocity verlet */
1754 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1755 !EI_VV(inputrec->eI))
1757 gmx_incons("update_coords called for velocity without VV integrator");
1761 homenr = md->homenr;
1762 nrend = start+homenr;
1764 xprime = get_xprime(state, upd);
1766 dt = inputrec->delta_t;
1769 /* We need to update the NMR restraint history when time averaging is used */
1770 if (state->flags & (1<<estDISRE_RM3TAV))
1772 update_disres_history(fcd, &state->hist);
1774 if (state->flags & (1<<estORIRE_DTAV))
1776 update_orires_history(fcd, &state->hist);
1780 bNH = inputrec->etc == etcNOSEHOOVER;
1781 bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1783 if (bDoLR && inputrec->nstcalclr > 1 && !EI_VV(inputrec->eI)) /* get this working with VV? */
1785 /* Store the total force + nstcalclr-1 times the LR force
1786 * in forces_lr, so it can be used in a normal update algorithm
1787 * to produce twin time stepping.
1789 /* is this correct in the new construction? MRS */
1790 combine_forces(inputrec->nstcalclr, constr, inputrec, md, idef, cr,
1791 step, state, bMolPBC,
1792 start, nrend, f, f_lr, nrnb);
1800 /* ############# START The update of velocities and positions ######### */
1802 dump_it_all(fplog, "Before update",
1803 state->natoms, state->x, xprime, state->v, force);
1805 if (inputrec->eI == eiSD2)
1807 check_sd2_work_data_allocation(upd->sd, nrend);
1809 do_update_sd2_Tconsts(upd->sd,
1810 inputrec->opts.ngtc,
1811 inputrec->opts.tau_t,
1812 inputrec->opts.ref_t);
1814 if (inputrec->eI == eiBD)
1816 do_update_bd_Tconsts(dt, inputrec->bd_fric,
1817 inputrec->opts.ngtc, inputrec->opts.ref_t,
1821 nth = gmx_omp_nthreads_get(emntUpdate);
1823 #pragma omp parallel for num_threads(nth) schedule(static) private(alpha)
1824 for (th = 0; th < nth; th++)
1826 int start_th, end_th;
1828 start_th = start + ((nrend-start)* th )/nth;
1829 end_th = start + ((nrend-start)*(th+1))/nth;
1831 switch (inputrec->eI)
1834 if (ekind->cosacc.cos_accel == 0)
1836 do_update_md(start_th, end_th, dt,
1837 ekind->tcstat, state->nosehoover_vxi,
1838 ekind->bNEMD, ekind->grpstat, inputrec->opts.acc,
1839 inputrec->opts.nFreeze,
1840 md->invmass, md->ptype,
1841 md->cFREEZE, md->cACC, md->cTC,
1842 state->x, xprime, state->v, force, M,
1847 do_update_visc(start_th, end_th, dt,
1848 ekind->tcstat, state->nosehoover_vxi,
1849 md->invmass, md->ptype,
1850 md->cTC, state->x, xprime, state->v, force, M,
1852 ekind->cosacc.cos_accel,
1858 do_update_sd1(upd->sd,
1859 start_th, end_th, dt,
1860 inputrec->opts.acc, inputrec->opts.nFreeze,
1861 md->invmass, md->ptype,
1862 md->cFREEZE, md->cACC, md->cTC,
1863 state->x, xprime, state->v, force,
1864 inputrec->opts.ngtc, inputrec->opts.tau_t, inputrec->opts.ref_t,
1865 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1868 /* The SD update is done in 2 parts, because an extra constraint step
1871 do_update_sd2(upd->sd,
1872 bInitStep, start_th, end_th,
1873 inputrec->opts.acc, inputrec->opts.nFreeze,
1874 md->invmass, md->ptype,
1875 md->cFREEZE, md->cACC, md->cTC,
1876 state->x, xprime, state->v, force, state->sd_X,
1877 inputrec->opts.tau_t,
1878 TRUE, step, inputrec->ld_seed,
1879 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1882 do_update_bd(start_th, end_th, dt,
1883 inputrec->opts.nFreeze, md->invmass, md->ptype,
1884 md->cFREEZE, md->cTC,
1885 state->x, xprime, state->v, force,
1888 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1892 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
1897 do_update_vv_vel(start_th, end_th, dt,
1898 inputrec->opts.acc, inputrec->opts.nFreeze,
1899 md->invmass, md->ptype,
1900 md->cFREEZE, md->cACC,
1902 (bNH || bPR), state->veta, alpha);
1905 do_update_vv_pos(start_th, end_th, dt,
1906 inputrec->opts.nFreeze,
1907 md->ptype, md->cFREEZE,
1908 state->x, xprime, state->v,
1909 (bNH || bPR), state->veta);
1914 gmx_fatal(FARGS, "Don't know how to update coordinates");
1922 void correct_ekin(FILE *log, int start, int end, rvec v[], rvec vcm, real mass[],
1923 real tmass, tensor ekin)
1926 * This is a debugging routine. It should not be called for production code
1928 * The kinetic energy should calculated according to:
1929 * Ekin = 1/2 m (v-vcm)^2
1930 * However the correction is not always applied, since vcm may not be
1931 * known in time and we compute
1932 * Ekin' = 1/2 m v^2 instead
1933 * This can be corrected afterwards by computing
1934 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
1936 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
1943 /* Local particles */
1946 /* Processor dependent part. */
1948 for (i = start; (i < end); i++)
1952 for (j = 0; (j < DIM); j++)
1958 svmul(1/tmass, vcm, vcm);
1959 svmul(0.5, vcm, hvcm);
1961 for (j = 0; (j < DIM); j++)
1963 for (k = 0; (k < DIM); k++)
1965 dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
1968 pr_rvecs(log, 0, "dekin", dekin, DIM);
1969 pr_rvecs(log, 0, " ekin", ekin, DIM);
1970 fprintf(log, "dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
1971 trace(dekin), trace(ekin), vcm[XX], vcm[YY], vcm[ZZ]);
1972 fprintf(log, "mv = (%8.4f %8.4f %8.4f)\n",
1973 mv[XX], mv[YY], mv[ZZ]);
1976 extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_int64_t step, const t_commrec *cr,
1977 t_mdatoms *md, t_state *state, gmx_update_t upd, gmx_constr_t constr)
1981 real rate = (ir->delta_t)/ir->opts.tau_t[0];
1983 if (ir->etc == etcANDERSEN && constr != NULL)
1985 gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
1988 /* proceed with andersen if 1) it's fixed probability per
1989 particle andersen or 2) it's massive andersen and it's tau_t/dt */
1990 if ((ir->etc == etcANDERSEN) || do_per_step(step, (int)(1.0/rate)))
1992 andersen_tcoupl(ir, step, cr, md, state, rate,
1993 upd->sd->randomize_group, upd->sd->boltzfac);