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"
48 #include "gromacs/math/units.h"
50 #include "gromacs/math/vec.h"
52 #include "gromacs/random/random.h"
61 #include "gmx_omp_nthreads.h"
63 #include "gromacs/fileio/confio.h"
64 #include "gromacs/pbcutil/mshift.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/pulling/pull.h"
67 #include "gromacs/timing/wallcycle.h"
68 #include "gromacs/utility/futil.h"
69 #include "gromacs/utility/gmxomp.h"
70 #include "gromacs/utility/smalloc.h"
72 /*For debugging, start at v(-dt/2) for velolcity verlet -- uncomment next line */
73 /*#define STARTFROMDT2*/
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 /* Variables for the deform algorithm */
113 gmx_int64_t deformref_step;
114 matrix deformref_box;
118 static void do_update_md(int start, int nrend, double dt,
119 t_grp_tcstat *tcstat,
121 gmx_bool bNEMD, t_grp_acc *gstat, rvec accel[],
124 unsigned short ptype[], unsigned short cFREEZE[],
125 unsigned short cACC[], unsigned short cTC[],
126 rvec x[], rvec xprime[], rvec v[],
128 gmx_bool bNH, gmx_bool bPR)
131 int gf = 0, ga = 0, gt = 0;
133 real vn, vv, va, vb, vnrel;
139 /* Update with coupling to extended ensembles, used for
140 * Nose-Hoover and Parrinello-Rahman coupling
141 * Nose-Hoover uses the reversible leap-frog integrator from
142 * Holian et al. Phys Rev E 52(3) : 2338, 1995
144 for (n = start; n < nrend; n++)
159 lg = tcstat[gt].lambda;
164 rvec_sub(v[n], gstat[ga].u, vrel);
166 for (d = 0; d < DIM; d++)
168 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
170 vnrel = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
171 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
172 /* do not scale the mean velocities u */
173 vn = gstat[ga].u[d] + accel[ga][d]*dt + vnrel;
175 xprime[n][d] = x[n][d]+vn*dt;
180 xprime[n][d] = x[n][d];
185 else if (cFREEZE != NULL ||
186 nFreeze[0][XX] || nFreeze[0][YY] || nFreeze[0][ZZ] ||
189 /* Update with Berendsen/v-rescale coupling and freeze or NEMD */
190 for (n = start; n < nrend; n++)
192 w_dt = invmass[n]*dt;
205 lg = tcstat[gt].lambda;
207 for (d = 0; d < DIM; d++)
210 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
212 vv = lg*vn + f[n][d]*w_dt;
214 /* do not scale the mean velocities u */
216 va = vv + accel[ga][d]*dt;
217 vb = va + (1.0-lg)*u;
219 xprime[n][d] = x[n][d]+vb*dt;
224 xprime[n][d] = x[n][d];
231 /* Plain update with Berendsen/v-rescale coupling */
232 for (n = start; n < nrend; n++)
234 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
236 w_dt = invmass[n]*dt;
241 lg = tcstat[gt].lambda;
243 for (d = 0; d < DIM; d++)
245 vn = lg*v[n][d] + f[n][d]*w_dt;
247 xprime[n][d] = x[n][d] + vn*dt;
252 for (d = 0; d < DIM; d++)
255 xprime[n][d] = x[n][d];
262 static void do_update_vv_vel(int start, int nrend, double dt,
263 rvec accel[], ivec nFreeze[], real invmass[],
264 unsigned short ptype[], unsigned short cFREEZE[],
265 unsigned short cACC[], rvec v[], rvec f[],
266 gmx_bool bExtended, real veta, real alpha)
271 real u, vn, vv, va, vb, vnrel;
277 g = 0.25*dt*veta*alpha;
279 mv2 = series_sinhx(g);
286 for (n = start; n < nrend; n++)
288 w_dt = invmass[n]*dt;
298 for (d = 0; d < DIM; d++)
300 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
302 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
310 } /* do_update_vv_vel */
312 static void do_update_vv_pos(int start, int nrend, double dt,
314 unsigned short ptype[], unsigned short cFREEZE[],
315 rvec x[], rvec xprime[], rvec v[],
316 gmx_bool bExtended, real veta)
323 /* Would it make more sense if Parrinello-Rahman was put here? */
328 mr2 = series_sinhx(g);
336 for (n = start; n < nrend; n++)
344 for (d = 0; d < DIM; d++)
346 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
348 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
352 xprime[n][d] = x[n][d];
356 } /* do_update_vv_pos */
358 static void do_update_visc(int start, int nrend, double dt,
359 t_grp_tcstat *tcstat,
362 unsigned short ptype[], unsigned short cTC[],
363 rvec x[], rvec xprime[], rvec v[],
364 rvec f[], matrix M, matrix box, real
365 cos_accel, real vcos,
366 gmx_bool bNH, gmx_bool bPR)
371 real lg, vxi = 0, vv;
376 fac = 2*M_PI/(box[ZZ][ZZ]);
380 /* Update with coupling to extended ensembles, used for
381 * Nose-Hoover and Parrinello-Rahman coupling
383 for (n = start; n < nrend; n++)
390 lg = tcstat[gt].lambda;
391 cosz = cos(fac*x[n][ZZ]);
393 copy_rvec(v[n], vrel);
401 for (d = 0; d < DIM; d++)
405 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
407 vn = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
408 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
411 vn += vc + dt*cosz*cos_accel;
414 xprime[n][d] = x[n][d]+vn*dt;
418 xprime[n][d] = x[n][d];
425 /* Classic version of update, used with berendsen coupling */
426 for (n = start; n < nrend; n++)
428 w_dt = invmass[n]*dt;
433 lg = tcstat[gt].lambda;
434 cosz = cos(fac*x[n][ZZ]);
436 for (d = 0; d < DIM; d++)
440 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
445 /* Do not scale the cosine velocity profile */
446 vv = vc + lg*(vn - vc + f[n][d]*w_dt);
447 /* Add the cosine accelaration profile */
448 vv += dt*cosz*cos_accel;
452 vv = lg*(vn + f[n][d]*w_dt);
455 xprime[n][d] = x[n][d]+vv*dt;
460 xprime[n][d] = x[n][d];
467 static gmx_stochd_t *init_stochd(t_inputrec *ir)
476 ngtc = ir->opts.ngtc;
480 snew(sd->bd_rf, ngtc);
482 else if (EI_SD(ir->eI))
485 snew(sd->sdsig, ngtc);
488 for (n = 0; n < ngtc; n++)
490 if (ir->opts.tau_t[n] > 0)
492 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
493 sdc[n].eph = exp(sdc[n].gdt/2);
494 sdc[n].emh = exp(-sdc[n].gdt/2);
495 sdc[n].em = exp(-sdc[n].gdt);
499 /* No friction and noise on this group */
505 if (sdc[n].gdt >= 0.05)
507 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
508 - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
509 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
510 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
515 /* Seventh order expansions for small y */
516 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
517 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))));
518 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
522 fprintf(debug, "SD const tc-grp %d: b %g c %g d %g\n",
523 n, sdc[n].b, sdc[n].c, sdc[n].d);
527 else if (ETC_ANDERSEN(ir->etc))
536 snew(sd->randomize_group, ngtc);
537 snew(sd->boltzfac, ngtc);
539 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
540 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
542 for (n = 0; n < ngtc; n++)
544 reft = max(0.0, opts->ref_t[n]);
545 if ((opts->tau_t[n] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
547 sd->randomize_group[n] = TRUE;
548 sd->boltzfac[n] = BOLTZ*opts->ref_t[n];
552 sd->randomize_group[n] = FALSE;
559 gmx_update_t init_update(t_inputrec *ir)
565 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
567 upd->sd = init_stochd(ir);
576 static void do_update_sd1(gmx_stochd_t *sd,
577 int start, int nrend, double dt,
578 rvec accel[], ivec nFreeze[],
579 real invmass[], unsigned short ptype[],
580 unsigned short cFREEZE[], unsigned short cACC[],
581 unsigned short cTC[],
582 rvec x[], rvec xprime[], rvec v[], rvec f[],
583 int ngtc, real tau_t[], real ref_t[],
584 gmx_int64_t step, int seed, int* gatindex)
589 int gf = 0, ga = 0, gt = 0;
596 for (n = 0; n < ngtc; n++)
599 /* The mass is encounted for later, since this differs per atom */
600 sig[n].V = sqrt(kT*(1 - sdc[n].em*sdc[n].em));
603 for (n = start; n < nrend; n++)
606 int ng = gatindex ? gatindex[n] : n;
607 ism = sqrt(invmass[n]);
621 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
622 for (d = 0; d < DIM; d++)
624 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
626 sd_V = ism*sig[gt].V*rnd[d];
628 v[n][d] = v[n][d]*sdc[gt].em
629 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
632 xprime[n][d] = x[n][d] + v[n][d]*dt;
637 xprime[n][d] = x[n][d];
643 static void check_sd2_work_data_allocation(gmx_stochd_t *sd, int nrend)
645 if (nrend > sd->sd_V_nalloc)
647 sd->sd_V_nalloc = over_alloc_dd(nrend);
648 srenew(sd->sd_V, sd->sd_V_nalloc);
652 static void do_update_sd2_Tconsts(gmx_stochd_t *sd,
657 /* This is separated from the update below, because it is single threaded */
666 for (gt = 0; gt < ngtc; gt++)
668 kT = BOLTZ*ref_t[gt];
669 /* The mass is encounted for later, since this differs per atom */
670 sig[gt].V = sqrt(kT*(1-sdc[gt].em));
671 sig[gt].X = sqrt(kT*sqr(tau_t[gt])*sdc[gt].c);
672 sig[gt].Yv = sqrt(kT*sdc[gt].b/sdc[gt].c);
673 sig[gt].Yx = sqrt(kT*sqr(tau_t[gt])*sdc[gt].b/(1-sdc[gt].em));
677 static void do_update_sd2(gmx_stochd_t *sd,
679 int start, int nrend,
680 rvec accel[], ivec nFreeze[],
681 real invmass[], unsigned short ptype[],
682 unsigned short cFREEZE[], unsigned short cACC[],
683 unsigned short cTC[],
684 rvec x[], rvec xprime[], rvec v[], rvec f[],
687 gmx_bool bFirstHalf, gmx_int64_t step, int seed,
692 /* The random part of the velocity update, generated in the first
693 * half of the update, needs to be remembered for the second half.
697 int gf = 0, ga = 0, gt = 0;
698 real vn = 0, Vmh, Xmh;
706 for (n = start; n < nrend; n++)
708 real rnd[6], rndi[3];
709 ng = gatindex ? gatindex[n] : n;
710 ism = sqrt(invmass[n]);
724 gmx_rng_cycle_6gaussian_table(step*2+(bFirstHalf ? 1 : 2), ng, seed, RND_SEED_UPDATE, rnd);
727 gmx_rng_cycle_3gaussian_table(step*2, ng, seed, RND_SEED_UPDATE, rndi);
729 for (d = 0; d < DIM; d++)
735 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
741 sd_X[n][d] = ism*sig[gt].X*rndi[d];
743 Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
744 + ism*sig[gt].Yv*rnd[d*2];
745 sd_V[n][d] = ism*sig[gt].V*rnd[d*2+1];
747 v[n][d] = vn*sdc[gt].em
748 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
749 + sd_V[n][d] - sdc[gt].em*Vmh;
751 xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
755 /* Correct the velocities for the constraints.
756 * This operation introduces some inaccuracy,
757 * since the velocity is determined from differences in coordinates.
760 (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
762 Xmh = sd_V[n][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
763 + ism*sig[gt].Yx*rnd[d*2];
764 sd_X[n][d] = ism*sig[gt].X*rnd[d*2+1];
766 xprime[n][d] += sd_X[n][d] - Xmh;
775 xprime[n][d] = x[n][d];
782 static void do_update_bd_Tconsts(double dt, real friction_coefficient,
783 int ngtc, const real ref_t[],
786 /* This is separated from the update below, because it is single threaded */
789 if (friction_coefficient != 0)
791 for (gt = 0; gt < ngtc; gt++)
793 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]/(friction_coefficient*dt));
798 for (gt = 0; gt < ngtc; gt++)
800 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]);
805 static void do_update_bd(int start, int nrend, double dt,
807 real invmass[], unsigned short ptype[],
808 unsigned short cFREEZE[], unsigned short cTC[],
809 rvec x[], rvec xprime[], rvec v[],
810 rvec f[], real friction_coefficient,
811 real *rf, gmx_int64_t step, int seed,
814 /* note -- these appear to be full step velocities . . . */
820 if (friction_coefficient != 0)
822 invfr = 1.0/friction_coefficient;
825 for (n = start; (n < nrend); n++)
828 int ng = gatindex ? gatindex[n] : n;
838 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
839 for (d = 0; (d < DIM); d++)
841 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
843 if (friction_coefficient != 0)
845 vn = invfr*f[n][d] + rf[gt]*rnd[d];
849 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
850 vn = 0.5*invmass[n]*f[n][d]*dt
851 + sqrt(0.5*invmass[n])*rf[gt]*rnd[d];
855 xprime[n][d] = x[n][d]+vn*dt;
860 xprime[n][d] = x[n][d];
866 static void dump_it_all(FILE gmx_unused *fp, const char gmx_unused *title,
867 int gmx_unused natoms, rvec gmx_unused x[], rvec gmx_unused xp[],
868 rvec gmx_unused v[], rvec gmx_unused f[])
873 fprintf(fp, "%s\n", title);
874 pr_rvecs(fp, 0, "x", x, natoms);
875 pr_rvecs(fp, 0, "xp", xp, natoms);
876 pr_rvecs(fp, 0, "v", v, natoms);
877 pr_rvecs(fp, 0, "f", f, natoms);
882 static void calc_ke_part_normal(rvec v[], t_grpopts *opts, t_mdatoms *md,
883 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel,
884 gmx_bool bSaveEkinOld)
887 t_grp_tcstat *tcstat = ekind->tcstat;
888 t_grp_acc *grpstat = ekind->grpstat;
891 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
892 an option, but not supported now. Additionally, if we are doing iterations.
893 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
894 bSavEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't copy over the ekinh_old.
895 If FALSE, we overrwrite it.
898 /* group velocities are calculated in update_ekindata and
899 * accumulated in acumulate_groups.
900 * Now the partial global and groups ekin.
902 for (g = 0; (g < opts->ngtc); g++)
907 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
911 clear_mat(tcstat[g].ekinf);
915 clear_mat(tcstat[g].ekinh);
919 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
922 ekind->dekindl_old = ekind->dekindl;
924 nthread = gmx_omp_nthreads_get(emntUpdate);
926 #pragma omp parallel for num_threads(nthread) schedule(static)
927 for (thread = 0; thread < nthread; thread++)
929 int start_t, end_t, n;
937 start_t = ((thread+0)*md->homenr)/nthread;
938 end_t = ((thread+1)*md->homenr)/nthread;
940 ekin_sum = ekind->ekin_work[thread];
941 dekindl_sum = ekind->dekindl_work[thread];
943 for (gt = 0; gt < opts->ngtc; gt++)
945 clear_mat(ekin_sum[gt]);
951 for (n = start_t; n < end_t; n++)
961 hm = 0.5*md->massT[n];
963 for (d = 0; (d < DIM); d++)
965 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
967 for (d = 0; (d < DIM); d++)
969 for (m = 0; (m < DIM); m++)
971 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
972 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
975 if (md->nMassPerturbed && md->bPerturbed[n])
978 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
984 for (thread = 0; thread < nthread; thread++)
986 for (g = 0; g < opts->ngtc; g++)
990 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
995 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1000 ekind->dekindl += *ekind->dekindl_work[thread];
1003 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1006 static void calc_ke_part_visc(matrix box, rvec x[], rvec v[],
1007 t_grpopts *opts, t_mdatoms *md,
1008 gmx_ekindata_t *ekind,
1009 t_nrnb *nrnb, gmx_bool bEkinAveVel)
1011 int start = 0, homenr = md->homenr;
1012 int g, d, n, m, gt = 0;
1015 t_grp_tcstat *tcstat = ekind->tcstat;
1016 t_cos_acc *cosacc = &(ekind->cosacc);
1021 for (g = 0; g < opts->ngtc; g++)
1023 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1024 clear_mat(ekind->tcstat[g].ekinh);
1026 ekind->dekindl_old = ekind->dekindl;
1028 fac = 2*M_PI/box[ZZ][ZZ];
1031 for (n = start; n < start+homenr; n++)
1037 hm = 0.5*md->massT[n];
1039 /* Note that the times of x and v differ by half a step */
1040 /* MRS -- would have to be changed for VV */
1041 cosz = cos(fac*x[n][ZZ]);
1042 /* Calculate the amplitude of the new velocity profile */
1043 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1045 copy_rvec(v[n], v_corrt);
1046 /* Subtract the profile for the kinetic energy */
1047 v_corrt[XX] -= cosz*cosacc->vcos;
1048 for (d = 0; (d < DIM); d++)
1050 for (m = 0; (m < DIM); m++)
1052 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1055 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1059 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1063 if (md->nPerturbed && md->bPerturbed[n])
1065 dekindl += 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1068 ekind->dekindl = dekindl;
1069 cosacc->mvcos = mvcos;
1071 inc_nrnb(nrnb, eNR_EKIN, homenr);
1074 void calc_ke_part(t_state *state, t_grpopts *opts, t_mdatoms *md,
1075 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
1077 if (ekind->cosacc.cos_accel == 0)
1079 calc_ke_part_normal(state->v, opts, md, ekind, nrnb, bEkinAveVel, bSaveEkinOld);
1083 calc_ke_part_visc(state->box, state->x, state->v, opts, md, ekind, nrnb, bEkinAveVel);
1087 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1089 ekinstate->ekin_n = ir->opts.ngtc;
1090 snew(ekinstate->ekinh, ekinstate->ekin_n);
1091 snew(ekinstate->ekinf, ekinstate->ekin_n);
1092 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1093 snew(ekinstate->ekinscalef_nhc, ekinstate->ekin_n);
1094 snew(ekinstate->ekinscaleh_nhc, ekinstate->ekin_n);
1095 snew(ekinstate->vscale_nhc, ekinstate->ekin_n);
1096 ekinstate->dekindl = 0;
1097 ekinstate->mvcos = 0;
1100 void update_ekinstate(ekinstate_t *ekinstate, gmx_ekindata_t *ekind)
1104 for (i = 0; i < ekinstate->ekin_n; i++)
1106 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1107 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1108 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1109 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1110 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1111 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1114 copy_mat(ekind->ekin, ekinstate->ekin_total);
1115 ekinstate->dekindl = ekind->dekindl;
1116 ekinstate->mvcos = ekind->cosacc.mvcos;
1120 void restore_ekinstate_from_state(t_commrec *cr,
1121 gmx_ekindata_t *ekind, ekinstate_t *ekinstate)
1127 for (i = 0; i < ekinstate->ekin_n; i++)
1129 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1130 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1131 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1132 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1133 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1134 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1137 copy_mat(ekinstate->ekin_total, ekind->ekin);
1139 ekind->dekindl = ekinstate->dekindl;
1140 ekind->cosacc.mvcos = ekinstate->mvcos;
1141 n = ekinstate->ekin_n;
1146 gmx_bcast(sizeof(n), &n, cr);
1147 for (i = 0; i < n; i++)
1149 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1150 ekind->tcstat[i].ekinh[0], cr);
1151 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1152 ekind->tcstat[i].ekinf[0], cr);
1153 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1154 ekind->tcstat[i].ekinh_old[0], cr);
1156 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1157 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1158 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1159 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1160 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1161 &(ekind->tcstat[i].vscale_nhc), cr);
1163 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1164 ekind->ekin[0], cr);
1166 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1167 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1171 void set_deform_reference_box(gmx_update_t upd, gmx_int64_t step, matrix box)
1173 upd->deformref_step = step;
1174 copy_mat(box, upd->deformref_box);
1177 static void deform(gmx_update_t upd,
1178 int start, int homenr, rvec x[], matrix box, matrix *scale_tot,
1179 const t_inputrec *ir, gmx_int64_t step)
1181 matrix bnew, invbox, mu;
1185 elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1186 copy_mat(box, bnew);
1187 for (i = 0; i < DIM; i++)
1189 for (j = 0; j < DIM; j++)
1191 if (ir->deform[i][j] != 0)
1194 upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1198 /* We correct the off-diagonal elements,
1199 * which can grow indefinitely during shearing,
1200 * so the shifts do not get messed up.
1202 for (i = 1; i < DIM; i++)
1204 for (j = i-1; j >= 0; j--)
1206 while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1208 rvec_dec(bnew[i], bnew[j]);
1210 while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1212 rvec_inc(bnew[i], bnew[j]);
1216 m_inv_ur0(box, invbox);
1217 copy_mat(bnew, box);
1218 mmul_ur0(box, invbox, mu);
1220 for (i = start; i < start+homenr; i++)
1222 x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1223 x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1224 x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1228 /* The transposes of the scaling matrices are stored,
1229 * so we need to do matrix multiplication in the inverse order.
1231 mmul_ur0(*scale_tot, mu, *scale_tot);
1235 void update_tcouple(gmx_int64_t step,
1236 t_inputrec *inputrec,
1238 gmx_ekindata_t *ekind,
1243 gmx_bool bTCouple = FALSE;
1245 int i, start, end, homenr, offset;
1247 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1248 if (inputrec->etc != etcNO &&
1249 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec)))
1251 /* We should only couple after a step where energies were determined (for leapfrog versions)
1252 or the step energies are determined, for velocity verlet versions */
1254 if (EI_VV(inputrec->eI))
1262 bTCouple = (inputrec->nsttcouple == 1 ||
1263 do_per_step(step+inputrec->nsttcouple-offset,
1264 inputrec->nsttcouple));
1269 dttc = inputrec->nsttcouple*inputrec->delta_t;
1271 switch (inputrec->etc)
1276 berendsen_tcoupl(inputrec, ekind, dttc);
1279 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1280 state->nosehoover_xi, state->nosehoover_vxi, MassQ);
1283 vrescale_tcoupl(inputrec, step, ekind, dttc,
1284 state->therm_integral);
1287 /* rescale in place here */
1288 if (EI_VV(inputrec->eI))
1290 rescale_velocities(ekind, md, 0, md->homenr, state->v);
1295 /* Set the T scaling lambda to 1 to have no scaling */
1296 for (i = 0; (i < inputrec->opts.ngtc); i++)
1298 ekind->tcstat[i].lambda = 1.0;
1303 void update_pcouple(FILE *fplog,
1305 t_inputrec *inputrec,
1311 gmx_bool bPCouple = FALSE;
1315 /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1316 if (inputrec->epc != epcNO && (!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))))
1318 /* We should only couple after a step where energies were determined */
1319 bPCouple = (inputrec->nstpcouple == 1 ||
1320 do_per_step(step+inputrec->nstpcouple-1,
1321 inputrec->nstpcouple));
1324 clear_mat(pcoupl_mu);
1325 for (i = 0; i < DIM; i++)
1327 pcoupl_mu[i][i] = 1.0;
1334 dtpc = inputrec->nstpcouple*inputrec->delta_t;
1336 switch (inputrec->epc)
1338 /* We can always pcoupl, even if we did not sum the energies
1339 * the previous step, since state->pres_prev is only updated
1340 * when the energies have been summed.
1344 case (epcBERENDSEN):
1347 berendsen_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev, state->box,
1351 case (epcPARRINELLORAHMAN):
1352 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1353 state->box, state->box_rel, state->boxv,
1354 M, pcoupl_mu, bInitStep);
1362 static rvec *get_xprime(const t_state *state, gmx_update_t upd)
1364 if (state->nalloc > upd->xp_nalloc)
1366 upd->xp_nalloc = state->nalloc;
1367 srenew(upd->xp, upd->xp_nalloc);
1373 static void combine_forces(gmx_update_t upd,
1375 gmx_constr_t constr,
1376 t_inputrec *ir, t_mdatoms *md, t_idef *idef,
1379 t_state *state, gmx_bool bMolPBC,
1380 int start, int nrend,
1381 rvec f[], rvec f_lr[],
1382 tensor *vir_lr_constr,
1387 /* f contains the short-range forces + the long range forces
1388 * which are stored separately in f_lr.
1391 if (constr != NULL && vir_lr_constr != NULL &&
1392 !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1394 /* We need to constrain the LR forces separately,
1395 * because due to the different pre-factor for the SR and LR
1396 * forces in the update algorithm, we have to correct
1397 * the constraint virial for the nstcalclr-1 extra f_lr.
1398 * Constrain only the additional LR part of the force.
1400 /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1405 xp = get_xprime(state, upd);
1407 fac = (nstcalclr - 1)*ir->delta_t*ir->delta_t;
1409 for (i = 0; i < md->homenr; i++)
1411 if (md->cFREEZE != NULL)
1413 gf = md->cFREEZE[i];
1415 for (d = 0; d < DIM; d++)
1417 if ((md->ptype[i] != eptVSite) &&
1418 (md->ptype[i] != eptShell) &&
1419 !ir->opts.nFreeze[gf][d])
1421 xp[i][d] = state->x[i][d] + fac*f_lr[i][d]*md->invmass[i];
1425 constrain(NULL, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
1426 state->x, xp, xp, bMolPBC, state->box, state->lambda[efptBONDED], NULL,
1427 NULL, vir_lr_constr, nrnb, econqCoord, ir->epc == epcMTTK, state->veta, state->veta);
1430 /* Add nstcalclr-1 times the LR force to the sum of both forces
1431 * and store the result in forces_lr.
1433 for (i = start; i < nrend; i++)
1435 for (d = 0; d < DIM; d++)
1437 f_lr[i][d] = f[i][d] + (nstcalclr - 1)*f_lr[i][d];
1442 void update_constraints(FILE *fplog,
1444 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1445 t_inputrec *inputrec, /* input record and box stuff */
1446 gmx_ekindata_t *ekind,
1451 rvec force[], /* forces on home particles */
1456 gmx_wallcycle_t wcycle,
1458 gmx_constr_t constr,
1459 gmx_bool bFirstHalf,
1463 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1466 int start, homenr, nrend, i, n, m, g, d;
1468 rvec *vbuf, *xprime = NULL;
1475 if (bFirstHalf && !EI_VV(inputrec->eI))
1480 /* for now, SD update is here -- though it really seems like it
1481 should be reformulated as a velocity verlet method, since it has two parts */
1484 homenr = md->homenr;
1485 nrend = start+homenr;
1487 dt = inputrec->delta_t;
1492 * APPLY CONSTRAINTS:
1495 * When doing PR pressure coupling we have to constrain the
1496 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1497 * it is enough to do this once though, since the relative velocities
1498 * after this will be normal to the bond vector
1503 /* clear out constraints before applying */
1504 clear_mat(vir_part);
1506 xprime = get_xprime(state, upd);
1508 bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1509 bLog = (do_per_step(step, inputrec->nstlog) || bLastStep || (step < 0));
1510 bEner = (do_per_step(step, inputrec->nstenergy) || bLastStep);
1511 /* Constrain the coordinates xprime */
1512 wallcycle_start(wcycle, ewcCONSTR);
1513 if (EI_VV(inputrec->eI) && bFirstHalf)
1515 constrain(NULL, bLog, bEner, constr, idef,
1516 inputrec, ekind, cr, step, 1, md,
1517 state->x, state->v, state->v,
1518 bMolPBC, state->box,
1519 state->lambda[efptBONDED], dvdlambda,
1520 NULL, bCalcVir ? &vir_con : NULL, nrnb, econqVeloc,
1521 inputrec->epc == epcMTTK, state->veta, vetanew);
1525 constrain(NULL, bLog, bEner, constr, idef,
1526 inputrec, ekind, cr, step, 1, md,
1527 state->x, xprime, NULL,
1528 bMolPBC, state->box,
1529 state->lambda[efptBONDED], dvdlambda,
1530 state->v, bCalcVir ? &vir_con : NULL, nrnb, econqCoord,
1531 inputrec->epc == epcMTTK, state->veta, state->veta);
1533 wallcycle_stop(wcycle, ewcCONSTR);
1537 dump_it_all(fplog, "After Shake",
1538 state->natoms, state->x, xprime, state->v, force);
1542 if (inputrec->eI == eiSD2)
1544 /* A correction factor eph is needed for the SD constraint force */
1545 /* Here we can, unfortunately, not have proper corrections
1546 * for different friction constants, so we use the first one.
1548 for (i = 0; i < DIM; i++)
1550 for (m = 0; m < DIM; m++)
1552 vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1558 m_add(vir_part, vir_con, vir_part);
1562 pr_rvecs(debug, 0, "constraint virial", vir_part, DIM);
1568 if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1570 xprime = get_xprime(state, upd);
1572 nth = gmx_omp_nthreads_get(emntUpdate);
1574 #pragma omp parallel for num_threads(nth) schedule(static)
1575 for (th = 0; th < nth; th++)
1577 int start_th, end_th;
1579 start_th = start + ((nrend-start)* th )/nth;
1580 end_th = start + ((nrend-start)*(th+1))/nth;
1582 /* The second part of the SD integration */
1583 do_update_sd2(upd->sd,
1584 FALSE, start_th, end_th,
1585 inputrec->opts.acc, inputrec->opts.nFreeze,
1586 md->invmass, md->ptype,
1587 md->cFREEZE, md->cACC, md->cTC,
1588 state->x, xprime, state->v, force, state->sd_X,
1589 inputrec->opts.tau_t,
1590 FALSE, step, inputrec->ld_seed,
1591 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1593 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1597 /* Constrain the coordinates xprime */
1598 wallcycle_start(wcycle, ewcCONSTR);
1599 constrain(NULL, bLog, bEner, constr, idef,
1600 inputrec, NULL, cr, step, 1, md,
1601 state->x, xprime, NULL,
1602 bMolPBC, state->box,
1603 state->lambda[efptBONDED], dvdlambda,
1604 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
1605 wallcycle_stop(wcycle, ewcCONSTR);
1609 /* We must always unshift after updating coordinates; if we did not shake
1610 x was shifted in do_force */
1612 if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1614 if (graph && (graph->nnodes > 0))
1616 unshift_x(graph, state->box, state->x, upd->xp);
1617 if (TRICLINIC(state->box))
1619 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1623 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1628 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntUpdate)) schedule(static)
1629 for (i = start; i < nrend; i++)
1631 copy_rvec(upd->xp[i], state->x[i]);
1635 dump_it_all(fplog, "After unshift",
1636 state->natoms, state->x, upd->xp, state->v, force);
1638 /* ############# END the update of velocities and positions ######### */
1641 void update_box(FILE *fplog,
1643 t_inputrec *inputrec, /* input record and box stuff */
1646 rvec force[], /* forces on home particles */
1652 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE;
1655 int start, homenr, nrend, i, n, m, g;
1659 homenr = md->homenr;
1660 nrend = start+homenr;
1663 (inputrec->etc == etcNOSEHOOVER) ||
1664 (inputrec->epc == epcPARRINELLORAHMAN) ||
1665 (inputrec->epc == epcMTTK);
1667 dt = inputrec->delta_t;
1671 /* now update boxes */
1672 switch (inputrec->epc)
1676 case (epcBERENDSEN):
1677 berendsen_pscale(inputrec, pcoupl_mu, state->box, state->box_rel,
1678 start, homenr, state->x, md->cFREEZE, nrnb);
1680 case (epcPARRINELLORAHMAN):
1681 /* The box velocities were updated in do_pr_pcoupl in the update
1682 * iteration, but we dont change the box vectors until we get here
1683 * since we need to be able to shift/unshift above.
1685 for (i = 0; i < DIM; i++)
1687 for (m = 0; m <= i; m++)
1689 state->box[i][m] += dt*state->boxv[i][m];
1692 preserve_box_shape(inputrec, state->box_rel, state->box);
1694 /* Scale the coordinates */
1695 for (n = start; (n < start+homenr); n++)
1697 tmvmul_ur0(pcoupl_mu, state->x[n], state->x[n]);
1701 switch (inputrec->epct)
1703 case (epctISOTROPIC):
1704 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1705 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1706 Side length scales as exp(veta*dt) */
1708 msmul(state->box, exp(state->veta*dt), state->box);
1710 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1711 o If we assume isotropic scaling, and box length scaling
1712 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1713 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1714 determinant of B is L^DIM det(M), and the determinant
1715 of dB/dt is (dL/dT)^DIM det (M). veta will be
1716 (det(dB/dT)/det(B))^(1/3). Then since M =
1717 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1719 msmul(state->box, state->veta, state->boxv);
1729 if ((!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))) && scale_tot)
1731 /* The transposes of the scaling matrices are stored,
1732 * therefore we need to reverse the order in the multiplication.
1734 mmul_ur0(*scale_tot, pcoupl_mu, *scale_tot);
1737 if (DEFORM(*inputrec))
1739 deform(upd, start, homenr, state->x, state->box, scale_tot, inputrec, step);
1742 dump_it_all(fplog, "After update",
1743 state->natoms, state->x, upd->xp, state->v, force);
1746 void update_coords(FILE *fplog,
1748 t_inputrec *inputrec, /* input record and box stuff */
1752 rvec *f, /* forces on home particles */
1755 tensor *vir_lr_constr,
1757 gmx_ekindata_t *ekind,
1762 t_commrec *cr, /* these shouldn't be here -- need to think about it */
1764 gmx_constr_t constr,
1767 gmx_bool bNH, bPR, bLastStep, bLog = FALSE, bEner = FALSE;
1769 real *imass, *imassin;
1772 int start, homenr, nrend, i, j, d, n, m, g;
1773 int blen0, blen1, iatom, jatom, nshake, nsettle, nconstr, nexpand;
1776 rvec *vcom, *xcom, *vall, *xall, *xin, *vin, *forcein, *fall, *xpall, *xprimein, *xprime;
1779 /* Running the velocity half does nothing except for velocity verlet */
1780 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1781 !EI_VV(inputrec->eI))
1783 gmx_incons("update_coords called for velocity without VV integrator");
1787 homenr = md->homenr;
1788 nrend = start+homenr;
1790 xprime = get_xprime(state, upd);
1792 dt = inputrec->delta_t;
1795 /* We need to update the NMR restraint history when time averaging is used */
1796 if (state->flags & (1<<estDISRE_RM3TAV))
1798 update_disres_history(fcd, &state->hist);
1800 if (state->flags & (1<<estORIRE_DTAV))
1802 update_orires_history(fcd, &state->hist);
1806 bNH = inputrec->etc == etcNOSEHOOVER;
1807 bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1809 if (bDoLR && inputrec->nstcalclr > 1 && !EI_VV(inputrec->eI)) /* get this working with VV? */
1811 /* Store the total force + nstcalclr-1 times the LR force
1812 * in forces_lr, so it can be used in a normal update algorithm
1813 * to produce twin time stepping.
1815 /* is this correct in the new construction? MRS */
1817 inputrec->nstcalclr, constr, inputrec, md, idef, cr,
1818 step, state, bMolPBC,
1819 start, nrend, f, f_lr, vir_lr_constr, nrnb);
1827 /* ############# START The update of velocities and positions ######### */
1829 dump_it_all(fplog, "Before update",
1830 state->natoms, state->x, xprime, state->v, force);
1832 if (inputrec->eI == eiSD2)
1834 check_sd2_work_data_allocation(upd->sd, nrend);
1836 do_update_sd2_Tconsts(upd->sd,
1837 inputrec->opts.ngtc,
1838 inputrec->opts.tau_t,
1839 inputrec->opts.ref_t);
1841 if (inputrec->eI == eiBD)
1843 do_update_bd_Tconsts(dt, inputrec->bd_fric,
1844 inputrec->opts.ngtc, inputrec->opts.ref_t,
1848 nth = gmx_omp_nthreads_get(emntUpdate);
1850 #pragma omp parallel for num_threads(nth) schedule(static) private(alpha)
1851 for (th = 0; th < nth; th++)
1853 int start_th, end_th;
1855 start_th = start + ((nrend-start)* th )/nth;
1856 end_th = start + ((nrend-start)*(th+1))/nth;
1858 switch (inputrec->eI)
1861 if (ekind->cosacc.cos_accel == 0)
1863 do_update_md(start_th, end_th, dt,
1864 ekind->tcstat, state->nosehoover_vxi,
1865 ekind->bNEMD, ekind->grpstat, inputrec->opts.acc,
1866 inputrec->opts.nFreeze,
1867 md->invmass, md->ptype,
1868 md->cFREEZE, md->cACC, md->cTC,
1869 state->x, xprime, state->v, force, M,
1874 do_update_visc(start_th, end_th, dt,
1875 ekind->tcstat, state->nosehoover_vxi,
1876 md->invmass, md->ptype,
1877 md->cTC, state->x, xprime, state->v, force, M,
1879 ekind->cosacc.cos_accel,
1885 do_update_sd1(upd->sd,
1886 start_th, end_th, dt,
1887 inputrec->opts.acc, inputrec->opts.nFreeze,
1888 md->invmass, md->ptype,
1889 md->cFREEZE, md->cACC, md->cTC,
1890 state->x, xprime, state->v, force,
1891 inputrec->opts.ngtc, inputrec->opts.tau_t, inputrec->opts.ref_t,
1892 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1895 /* The SD update is done in 2 parts, because an extra constraint step
1898 do_update_sd2(upd->sd,
1899 bInitStep, start_th, end_th,
1900 inputrec->opts.acc, inputrec->opts.nFreeze,
1901 md->invmass, md->ptype,
1902 md->cFREEZE, md->cACC, md->cTC,
1903 state->x, xprime, state->v, force, state->sd_X,
1904 inputrec->opts.tau_t,
1905 TRUE, step, inputrec->ld_seed,
1906 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1909 do_update_bd(start_th, end_th, dt,
1910 inputrec->opts.nFreeze, md->invmass, md->ptype,
1911 md->cFREEZE, md->cTC,
1912 state->x, xprime, state->v, force,
1915 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1919 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
1924 do_update_vv_vel(start_th, end_th, dt,
1925 inputrec->opts.acc, inputrec->opts.nFreeze,
1926 md->invmass, md->ptype,
1927 md->cFREEZE, md->cACC,
1929 (bNH || bPR), state->veta, alpha);
1932 do_update_vv_pos(start_th, end_th, dt,
1933 inputrec->opts.nFreeze,
1934 md->ptype, md->cFREEZE,
1935 state->x, xprime, state->v,
1936 (bNH || bPR), state->veta);
1941 gmx_fatal(FARGS, "Don't know how to update coordinates");
1949 void correct_ekin(FILE *log, int start, int end, rvec v[], rvec vcm, real mass[],
1950 real tmass, tensor ekin)
1953 * This is a debugging routine. It should not be called for production code
1955 * The kinetic energy should calculated according to:
1956 * Ekin = 1/2 m (v-vcm)^2
1957 * However the correction is not always applied, since vcm may not be
1958 * known in time and we compute
1959 * Ekin' = 1/2 m v^2 instead
1960 * This can be corrected afterwards by computing
1961 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
1963 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
1970 /* Local particles */
1973 /* Processor dependent part. */
1975 for (i = start; (i < end); i++)
1979 for (j = 0; (j < DIM); j++)
1985 svmul(1/tmass, vcm, vcm);
1986 svmul(0.5, vcm, hvcm);
1988 for (j = 0; (j < DIM); j++)
1990 for (k = 0; (k < DIM); k++)
1992 dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
1995 pr_rvecs(log, 0, "dekin", dekin, DIM);
1996 pr_rvecs(log, 0, " ekin", ekin, DIM);
1997 fprintf(log, "dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
1998 trace(dekin), trace(ekin), vcm[XX], vcm[YY], vcm[ZZ]);
1999 fprintf(log, "mv = (%8.4f %8.4f %8.4f)\n",
2000 mv[XX], mv[YY], mv[ZZ]);
2003 extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_int64_t step, const t_commrec *cr,
2004 t_mdatoms *md, t_state *state, gmx_update_t upd, gmx_constr_t constr)
2008 real rate = (ir->delta_t)/ir->opts.tau_t[0];
2010 if (ir->etc == etcANDERSEN && constr != NULL)
2012 gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
2015 /* proceed with andersen if 1) it's fixed probability per
2016 particle andersen or 2) it's massive andersen and it's tau_t/dt */
2017 if ((ir->etc == etcANDERSEN) || do_per_step(step, (int)(1.0/rate)))
2019 andersen_tcoupl(ir, step, cr, md, state, rate,
2020 upd->sd->randomize_group, upd->sd->boltzfac);