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.
43 #include "types/commrec.h"
46 #include "gromacs/math/units.h"
48 #include "gromacs/math/vec.h"
50 #include "gromacs/random/random.h"
59 #include "gmx_omp_nthreads.h"
61 #include "gromacs/fileio/confio.h"
62 #include "gromacs/pbcutil/mshift.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/pulling/pull.h"
65 #include "gromacs/timing/wallcycle.h"
66 #include "gromacs/utility/futil.h"
67 #include "gromacs/utility/gmxomp.h"
68 #include "gromacs/utility/smalloc.h"
70 /*For debugging, start at v(-dt/2) for velolcity verlet -- uncomment next line */
71 /*#define STARTFROMDT2*/
95 gmx_sd_sigma_t *sdsig;
98 /* andersen temperature control stuff */
99 gmx_bool *randomize_group;
103 typedef struct gmx_update
106 /* xprime for constraint algorithms */
110 /* Variables for the deform algorithm */
111 gmx_int64_t deformref_step;
112 matrix deformref_box;
116 static void do_update_md(int start, int nrend, double dt,
117 t_grp_tcstat *tcstat,
119 gmx_bool bNEMD, t_grp_acc *gstat, rvec accel[],
122 unsigned short ptype[], unsigned short cFREEZE[],
123 unsigned short cACC[], unsigned short cTC[],
124 rvec x[], rvec xprime[], rvec v[],
126 gmx_bool bNH, gmx_bool bPR)
129 int gf = 0, ga = 0, gt = 0;
131 real vn, vv, va, vb, vnrel;
137 /* Update with coupling to extended ensembles, used for
138 * Nose-Hoover and Parrinello-Rahman coupling
139 * Nose-Hoover uses the reversible leap-frog integrator from
140 * Holian et al. Phys Rev E 52(3) : 2338, 1995
142 for (n = start; n < nrend; n++)
157 lg = tcstat[gt].lambda;
162 rvec_sub(v[n], gstat[ga].u, vrel);
164 for (d = 0; d < DIM; d++)
166 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
168 vnrel = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
169 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
170 /* do not scale the mean velocities u */
171 vn = gstat[ga].u[d] + accel[ga][d]*dt + vnrel;
173 xprime[n][d] = x[n][d]+vn*dt;
178 xprime[n][d] = x[n][d];
183 else if (cFREEZE != NULL ||
184 nFreeze[0][XX] || nFreeze[0][YY] || nFreeze[0][ZZ] ||
187 /* Update with Berendsen/v-rescale coupling and freeze or NEMD */
188 for (n = start; n < nrend; n++)
190 w_dt = invmass[n]*dt;
203 lg = tcstat[gt].lambda;
205 for (d = 0; d < DIM; d++)
208 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
210 vv = lg*vn + f[n][d]*w_dt;
212 /* do not scale the mean velocities u */
214 va = vv + accel[ga][d]*dt;
215 vb = va + (1.0-lg)*u;
217 xprime[n][d] = x[n][d]+vb*dt;
222 xprime[n][d] = x[n][d];
229 /* Plain update with Berendsen/v-rescale coupling */
230 for (n = start; n < nrend; n++)
232 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
234 w_dt = invmass[n]*dt;
239 lg = tcstat[gt].lambda;
241 for (d = 0; d < DIM; d++)
243 vn = lg*v[n][d] + f[n][d]*w_dt;
245 xprime[n][d] = x[n][d] + vn*dt;
250 for (d = 0; d < DIM; d++)
253 xprime[n][d] = x[n][d];
260 static void do_update_vv_vel(int start, int nrend, double dt,
261 rvec accel[], ivec nFreeze[], real invmass[],
262 unsigned short ptype[], unsigned short cFREEZE[],
263 unsigned short cACC[], rvec v[], rvec f[],
264 gmx_bool bExtended, real veta, real alpha)
269 real u, vn, vv, va, vb, vnrel;
275 g = 0.25*dt*veta*alpha;
277 mv2 = series_sinhx(g);
284 for (n = start; n < nrend; n++)
286 w_dt = invmass[n]*dt;
296 for (d = 0; d < DIM; d++)
298 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
300 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
308 } /* do_update_vv_vel */
310 static void do_update_vv_pos(int start, int nrend, double dt,
312 unsigned short ptype[], unsigned short cFREEZE[],
313 rvec x[], rvec xprime[], rvec v[],
314 gmx_bool bExtended, real veta)
321 /* Would it make more sense if Parrinello-Rahman was put here? */
326 mr2 = series_sinhx(g);
334 for (n = start; n < nrend; n++)
342 for (d = 0; d < DIM; d++)
344 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
346 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
350 xprime[n][d] = x[n][d];
354 } /* do_update_vv_pos */
356 static void do_update_visc(int start, int nrend, double dt,
357 t_grp_tcstat *tcstat,
360 unsigned short ptype[], unsigned short cTC[],
361 rvec x[], rvec xprime[], rvec v[],
362 rvec f[], matrix M, matrix box, real
363 cos_accel, real vcos,
364 gmx_bool bNH, gmx_bool bPR)
369 real lg, vxi = 0, vv;
374 fac = 2*M_PI/(box[ZZ][ZZ]);
378 /* Update with coupling to extended ensembles, used for
379 * Nose-Hoover and Parrinello-Rahman coupling
381 for (n = start; n < nrend; n++)
388 lg = tcstat[gt].lambda;
389 cosz = cos(fac*x[n][ZZ]);
391 copy_rvec(v[n], vrel);
399 for (d = 0; d < DIM; d++)
403 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
405 vn = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
406 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
409 vn += vc + dt*cosz*cos_accel;
412 xprime[n][d] = x[n][d]+vn*dt;
416 xprime[n][d] = x[n][d];
423 /* Classic version of update, used with berendsen coupling */
424 for (n = start; n < nrend; n++)
426 w_dt = invmass[n]*dt;
431 lg = tcstat[gt].lambda;
432 cosz = cos(fac*x[n][ZZ]);
434 for (d = 0; d < DIM; d++)
438 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
443 /* Do not scale the cosine velocity profile */
444 vv = vc + lg*(vn - vc + f[n][d]*w_dt);
445 /* Add the cosine accelaration profile */
446 vv += dt*cosz*cos_accel;
450 vv = lg*(vn + f[n][d]*w_dt);
453 xprime[n][d] = x[n][d]+vv*dt;
458 xprime[n][d] = x[n][d];
465 static gmx_stochd_t *init_stochd(t_inputrec *ir)
474 ngtc = ir->opts.ngtc;
478 snew(sd->bd_rf, ngtc);
480 else if (EI_SD(ir->eI))
483 snew(sd->sdsig, ngtc);
486 for (n = 0; n < ngtc; n++)
488 if (ir->opts.tau_t[n] > 0)
490 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
491 sdc[n].eph = exp(sdc[n].gdt/2);
492 sdc[n].emh = exp(-sdc[n].gdt/2);
493 sdc[n].em = exp(-sdc[n].gdt);
497 /* No friction and noise on this group */
503 if (sdc[n].gdt >= 0.05)
505 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
506 - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
507 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
508 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
513 /* Seventh order expansions for small y */
514 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
515 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))));
516 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
520 fprintf(debug, "SD const tc-grp %d: b %g c %g d %g\n",
521 n, sdc[n].b, sdc[n].c, sdc[n].d);
525 else if (ETC_ANDERSEN(ir->etc))
534 snew(sd->randomize_group, ngtc);
535 snew(sd->boltzfac, ngtc);
537 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
538 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
540 for (n = 0; n < ngtc; n++)
542 reft = max(0.0, opts->ref_t[n]);
543 if ((opts->tau_t[n] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
545 sd->randomize_group[n] = TRUE;
546 sd->boltzfac[n] = BOLTZ*opts->ref_t[n];
550 sd->randomize_group[n] = FALSE;
557 gmx_update_t init_update(t_inputrec *ir)
563 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
565 upd->sd = init_stochd(ir);
574 static void do_update_sd1(gmx_stochd_t *sd,
575 int start, int nrend, double dt,
576 rvec accel[], ivec nFreeze[],
577 real invmass[], unsigned short ptype[],
578 unsigned short cFREEZE[], unsigned short cACC[],
579 unsigned short cTC[],
580 rvec x[], rvec xprime[], rvec v[], rvec f[],
581 int ngtc, real ref_t[],
583 gmx_bool bFirstHalfConstr,
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));
605 for (n = start; n < nrend; n++)
608 int ng = gatindex ? gatindex[n] : n;
610 ism = sqrt(invmass[n]);
624 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
626 for (d = 0; d < DIM; d++)
628 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
632 sd_V = ism*sig[gt].V*rnd[d];
633 vn = v[n][d] + (invmass[n]*f[n][d] + accel[ga][d])*dt;
634 v[n][d] = vn*sdc[gt].em + sd_V;
635 /* Here we include half of the friction+noise
636 * update of v into the integration of x.
638 xprime[n][d] = x[n][d] + 0.5*(vn + v[n][d])*dt;
643 xprime[n][d] = x[n][d];
650 /* We do have constraints */
651 if (bFirstHalfConstr)
653 /* First update without friction and noise */
656 for (n = start; n < nrend; n++)
673 for (d = 0; d < DIM; d++)
675 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
677 v[n][d] = v[n][d] + (im*f[n][d] + accel[ga][d])*dt;
678 xprime[n][d] = x[n][d] + v[n][d]*dt;
683 xprime[n][d] = x[n][d];
690 /* Update friction and noise only */
691 for (n = start; n < nrend; n++)
694 int ng = gatindex ? gatindex[n] : n;
696 ism = sqrt(invmass[n]);
710 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
712 for (d = 0; d < DIM; d++)
714 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
718 sd_V = ism*sig[gt].V*rnd[d];
720 v[n][d] = vn*sdc[gt].em + sd_V;
721 /* Add the friction and noise contribution only */
722 xprime[n][d] = xprime[n][d] + 0.5*(v[n][d] - vn)*dt;
730 static void check_sd2_work_data_allocation(gmx_stochd_t *sd, int nrend)
732 if (nrend > sd->sd_V_nalloc)
734 sd->sd_V_nalloc = over_alloc_dd(nrend);
735 srenew(sd->sd_V, sd->sd_V_nalloc);
739 static void do_update_sd2_Tconsts(gmx_stochd_t *sd,
744 /* This is separated from the update below, because it is single threaded */
753 for (gt = 0; gt < ngtc; gt++)
755 kT = BOLTZ*ref_t[gt];
756 /* The mass is encounted for later, since this differs per atom */
757 sig[gt].V = sqrt(kT*(1-sdc[gt].em));
758 sig[gt].X = sqrt(kT*sqr(tau_t[gt])*sdc[gt].c);
759 sig[gt].Yv = sqrt(kT*sdc[gt].b/sdc[gt].c);
760 sig[gt].Yx = sqrt(kT*sqr(tau_t[gt])*sdc[gt].b/(1-sdc[gt].em));
764 static void do_update_sd2(gmx_stochd_t *sd,
766 int start, int nrend,
767 rvec accel[], ivec nFreeze[],
768 real invmass[], unsigned short ptype[],
769 unsigned short cFREEZE[], unsigned short cACC[],
770 unsigned short cTC[],
771 rvec x[], rvec xprime[], rvec v[], rvec f[],
774 gmx_bool bFirstHalf, gmx_int64_t step, int seed,
779 /* The random part of the velocity update, generated in the first
780 * half of the update, needs to be remembered for the second half.
784 int gf = 0, ga = 0, gt = 0;
785 real vn = 0, Vmh, Xmh;
793 for (n = start; n < nrend; n++)
795 real rnd[6], rndi[3];
796 ng = gatindex ? gatindex[n] : n;
797 ism = sqrt(invmass[n]);
811 gmx_rng_cycle_6gaussian_table(step*2+(bFirstHalf ? 1 : 2), ng, seed, RND_SEED_UPDATE, rnd);
814 gmx_rng_cycle_3gaussian_table(step*2, ng, seed, RND_SEED_UPDATE, rndi);
816 for (d = 0; d < DIM; d++)
822 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
828 sd_X[n][d] = ism*sig[gt].X*rndi[d];
830 Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
831 + ism*sig[gt].Yv*rnd[d*2];
832 sd_V[n][d] = ism*sig[gt].V*rnd[d*2+1];
834 v[n][d] = vn*sdc[gt].em
835 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
836 + sd_V[n][d] - sdc[gt].em*Vmh;
838 xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
842 /* Correct the velocities for the constraints.
843 * This operation introduces some inaccuracy,
844 * since the velocity is determined from differences in coordinates.
847 (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
849 Xmh = sd_V[n][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
850 + ism*sig[gt].Yx*rnd[d*2];
851 sd_X[n][d] = ism*sig[gt].X*rnd[d*2+1];
853 xprime[n][d] += sd_X[n][d] - Xmh;
862 xprime[n][d] = x[n][d];
869 static void do_update_bd_Tconsts(double dt, real friction_coefficient,
870 int ngtc, const real ref_t[],
873 /* This is separated from the update below, because it is single threaded */
876 if (friction_coefficient != 0)
878 for (gt = 0; gt < ngtc; gt++)
880 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]/(friction_coefficient*dt));
885 for (gt = 0; gt < ngtc; gt++)
887 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]);
892 static void do_update_bd(int start, int nrend, double dt,
894 real invmass[], unsigned short ptype[],
895 unsigned short cFREEZE[], unsigned short cTC[],
896 rvec x[], rvec xprime[], rvec v[],
897 rvec f[], real friction_coefficient,
898 real *rf, gmx_int64_t step, int seed,
901 /* note -- these appear to be full step velocities . . . */
907 if (friction_coefficient != 0)
909 invfr = 1.0/friction_coefficient;
912 for (n = start; (n < nrend); n++)
915 int ng = gatindex ? gatindex[n] : n;
925 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
926 for (d = 0; (d < DIM); d++)
928 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
930 if (friction_coefficient != 0)
932 vn = invfr*f[n][d] + rf[gt]*rnd[d];
936 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
937 vn = 0.5*invmass[n]*f[n][d]*dt
938 + sqrt(0.5*invmass[n])*rf[gt]*rnd[d];
942 xprime[n][d] = x[n][d]+vn*dt;
947 xprime[n][d] = x[n][d];
953 static void dump_it_all(FILE gmx_unused *fp, const char gmx_unused *title,
954 int gmx_unused natoms, rvec gmx_unused x[], rvec gmx_unused xp[],
955 rvec gmx_unused v[], rvec gmx_unused f[])
960 fprintf(fp, "%s\n", title);
961 pr_rvecs(fp, 0, "x", x, natoms);
962 pr_rvecs(fp, 0, "xp", xp, natoms);
963 pr_rvecs(fp, 0, "v", v, natoms);
964 pr_rvecs(fp, 0, "f", f, natoms);
969 static void calc_ke_part_normal(rvec v[], t_grpopts *opts, t_mdatoms *md,
970 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel,
971 gmx_bool bSaveEkinOld)
974 t_grp_tcstat *tcstat = ekind->tcstat;
975 t_grp_acc *grpstat = ekind->grpstat;
978 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
979 an option, but not supported now. Additionally, if we are doing iterations.
980 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
981 bSavEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't copy over the ekinh_old.
982 If FALSE, we overrwrite it.
985 /* group velocities are calculated in update_ekindata and
986 * accumulated in acumulate_groups.
987 * Now the partial global and groups ekin.
989 for (g = 0; (g < opts->ngtc); g++)
994 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
998 clear_mat(tcstat[g].ekinf);
1002 clear_mat(tcstat[g].ekinh);
1006 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
1009 ekind->dekindl_old = ekind->dekindl;
1011 nthread = gmx_omp_nthreads_get(emntUpdate);
1013 #pragma omp parallel for num_threads(nthread) schedule(static)
1014 for (thread = 0; thread < nthread; thread++)
1016 int start_t, end_t, n;
1024 start_t = ((thread+0)*md->homenr)/nthread;
1025 end_t = ((thread+1)*md->homenr)/nthread;
1027 ekin_sum = ekind->ekin_work[thread];
1028 dekindl_sum = ekind->dekindl_work[thread];
1030 for (gt = 0; gt < opts->ngtc; gt++)
1032 clear_mat(ekin_sum[gt]);
1038 for (n = start_t; n < end_t; n++)
1048 hm = 0.5*md->massT[n];
1050 for (d = 0; (d < DIM); d++)
1052 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
1054 for (d = 0; (d < DIM); d++)
1056 for (m = 0; (m < DIM); m++)
1058 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1059 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
1062 if (md->nMassPerturbed && md->bPerturbed[n])
1065 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1071 for (thread = 0; thread < nthread; thread++)
1073 for (g = 0; g < opts->ngtc; g++)
1077 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
1082 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1087 ekind->dekindl += *ekind->dekindl_work[thread];
1090 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1093 static void calc_ke_part_visc(matrix box, rvec x[], rvec v[],
1094 t_grpopts *opts, t_mdatoms *md,
1095 gmx_ekindata_t *ekind,
1096 t_nrnb *nrnb, gmx_bool bEkinAveVel)
1098 int start = 0, homenr = md->homenr;
1099 int g, d, n, m, gt = 0;
1102 t_grp_tcstat *tcstat = ekind->tcstat;
1103 t_cos_acc *cosacc = &(ekind->cosacc);
1108 for (g = 0; g < opts->ngtc; g++)
1110 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1111 clear_mat(ekind->tcstat[g].ekinh);
1113 ekind->dekindl_old = ekind->dekindl;
1115 fac = 2*M_PI/box[ZZ][ZZ];
1118 for (n = start; n < start+homenr; n++)
1124 hm = 0.5*md->massT[n];
1126 /* Note that the times of x and v differ by half a step */
1127 /* MRS -- would have to be changed for VV */
1128 cosz = cos(fac*x[n][ZZ]);
1129 /* Calculate the amplitude of the new velocity profile */
1130 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1132 copy_rvec(v[n], v_corrt);
1133 /* Subtract the profile for the kinetic energy */
1134 v_corrt[XX] -= cosz*cosacc->vcos;
1135 for (d = 0; (d < DIM); d++)
1137 for (m = 0; (m < DIM); m++)
1139 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1142 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1146 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1150 if (md->nPerturbed && md->bPerturbed[n])
1152 /* The minus sign here might be confusing.
1153 * The kinetic contribution from dH/dl doesn't come from
1154 * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
1155 * where p are the momenta. The difference is only a minus sign.
1157 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1160 ekind->dekindl = dekindl;
1161 cosacc->mvcos = mvcos;
1163 inc_nrnb(nrnb, eNR_EKIN, homenr);
1166 void calc_ke_part(t_state *state, t_grpopts *opts, t_mdatoms *md,
1167 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
1169 if (ekind->cosacc.cos_accel == 0)
1171 calc_ke_part_normal(state->v, opts, md, ekind, nrnb, bEkinAveVel, bSaveEkinOld);
1175 calc_ke_part_visc(state->box, state->x, state->v, opts, md, ekind, nrnb, bEkinAveVel);
1179 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1181 ekinstate->ekin_n = ir->opts.ngtc;
1182 snew(ekinstate->ekinh, ekinstate->ekin_n);
1183 snew(ekinstate->ekinf, ekinstate->ekin_n);
1184 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1185 snew(ekinstate->ekinscalef_nhc, ekinstate->ekin_n);
1186 snew(ekinstate->ekinscaleh_nhc, ekinstate->ekin_n);
1187 snew(ekinstate->vscale_nhc, ekinstate->ekin_n);
1188 ekinstate->dekindl = 0;
1189 ekinstate->mvcos = 0;
1192 void update_ekinstate(ekinstate_t *ekinstate, gmx_ekindata_t *ekind)
1196 for (i = 0; i < ekinstate->ekin_n; i++)
1198 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1199 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1200 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1201 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1202 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1203 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1206 copy_mat(ekind->ekin, ekinstate->ekin_total);
1207 ekinstate->dekindl = ekind->dekindl;
1208 ekinstate->mvcos = ekind->cosacc.mvcos;
1212 void restore_ekinstate_from_state(t_commrec *cr,
1213 gmx_ekindata_t *ekind, ekinstate_t *ekinstate)
1219 for (i = 0; i < ekinstate->ekin_n; i++)
1221 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1222 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1223 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1224 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1225 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1226 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1229 copy_mat(ekinstate->ekin_total, ekind->ekin);
1231 ekind->dekindl = ekinstate->dekindl;
1232 ekind->cosacc.mvcos = ekinstate->mvcos;
1233 n = ekinstate->ekin_n;
1238 gmx_bcast(sizeof(n), &n, cr);
1239 for (i = 0; i < n; i++)
1241 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1242 ekind->tcstat[i].ekinh[0], cr);
1243 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1244 ekind->tcstat[i].ekinf[0], cr);
1245 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1246 ekind->tcstat[i].ekinh_old[0], cr);
1248 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1249 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1250 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1251 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1252 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1253 &(ekind->tcstat[i].vscale_nhc), cr);
1255 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1256 ekind->ekin[0], cr);
1258 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1259 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1263 void set_deform_reference_box(gmx_update_t upd, gmx_int64_t step, matrix box)
1265 upd->deformref_step = step;
1266 copy_mat(box, upd->deformref_box);
1269 static void deform(gmx_update_t upd,
1270 int start, int homenr, rvec x[], matrix box, matrix *scale_tot,
1271 const t_inputrec *ir, gmx_int64_t step)
1273 matrix bnew, invbox, mu;
1277 elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1278 copy_mat(box, bnew);
1279 for (i = 0; i < DIM; i++)
1281 for (j = 0; j < DIM; j++)
1283 if (ir->deform[i][j] != 0)
1286 upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1290 /* We correct the off-diagonal elements,
1291 * which can grow indefinitely during shearing,
1292 * so the shifts do not get messed up.
1294 for (i = 1; i < DIM; i++)
1296 for (j = i-1; j >= 0; j--)
1298 while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1300 rvec_dec(bnew[i], bnew[j]);
1302 while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1304 rvec_inc(bnew[i], bnew[j]);
1308 m_inv_ur0(box, invbox);
1309 copy_mat(bnew, box);
1310 mmul_ur0(box, invbox, mu);
1312 for (i = start; i < start+homenr; i++)
1314 x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1315 x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1316 x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1320 /* The transposes of the scaling matrices are stored,
1321 * so we need to do matrix multiplication in the inverse order.
1323 mmul_ur0(*scale_tot, mu, *scale_tot);
1327 void update_tcouple(gmx_int64_t step,
1328 t_inputrec *inputrec,
1330 gmx_ekindata_t *ekind,
1335 gmx_bool bTCouple = FALSE;
1337 int i, start, end, homenr, offset;
1339 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1340 if (inputrec->etc != etcNO &&
1341 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec)))
1343 /* We should only couple after a step where energies were determined (for leapfrog versions)
1344 or the step energies are determined, for velocity verlet versions */
1346 if (EI_VV(inputrec->eI))
1354 bTCouple = (inputrec->nsttcouple == 1 ||
1355 do_per_step(step+inputrec->nsttcouple-offset,
1356 inputrec->nsttcouple));
1361 dttc = inputrec->nsttcouple*inputrec->delta_t;
1363 switch (inputrec->etc)
1368 berendsen_tcoupl(inputrec, ekind, dttc);
1371 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1372 state->nosehoover_xi, state->nosehoover_vxi, MassQ);
1375 vrescale_tcoupl(inputrec, step, ekind, dttc,
1376 state->therm_integral);
1379 /* rescale in place here */
1380 if (EI_VV(inputrec->eI))
1382 rescale_velocities(ekind, md, 0, md->homenr, state->v);
1387 /* Set the T scaling lambda to 1 to have no scaling */
1388 for (i = 0; (i < inputrec->opts.ngtc); i++)
1390 ekind->tcstat[i].lambda = 1.0;
1395 void update_pcouple(FILE *fplog,
1397 t_inputrec *inputrec,
1403 gmx_bool bPCouple = FALSE;
1407 /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1408 if (inputrec->epc != epcNO && (!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))))
1410 /* We should only couple after a step where energies were determined */
1411 bPCouple = (inputrec->nstpcouple == 1 ||
1412 do_per_step(step+inputrec->nstpcouple-1,
1413 inputrec->nstpcouple));
1416 clear_mat(pcoupl_mu);
1417 for (i = 0; i < DIM; i++)
1419 pcoupl_mu[i][i] = 1.0;
1426 dtpc = inputrec->nstpcouple*inputrec->delta_t;
1428 switch (inputrec->epc)
1430 /* We can always pcoupl, even if we did not sum the energies
1431 * the previous step, since state->pres_prev is only updated
1432 * when the energies have been summed.
1436 case (epcBERENDSEN):
1439 berendsen_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev, state->box,
1443 case (epcPARRINELLORAHMAN):
1444 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1445 state->box, state->box_rel, state->boxv,
1446 M, pcoupl_mu, bInitStep);
1454 static rvec *get_xprime(const t_state *state, gmx_update_t upd)
1456 if (state->nalloc > upd->xp_nalloc)
1458 upd->xp_nalloc = state->nalloc;
1459 srenew(upd->xp, upd->xp_nalloc);
1465 static void combine_forces(gmx_update_t upd,
1467 gmx_constr_t constr,
1468 t_inputrec *ir, t_mdatoms *md, t_idef *idef,
1471 t_state *state, gmx_bool bMolPBC,
1472 int start, int nrend,
1473 rvec f[], rvec f_lr[],
1474 tensor *vir_lr_constr,
1479 /* f contains the short-range forces + the long range forces
1480 * which are stored separately in f_lr.
1483 if (constr != NULL && vir_lr_constr != NULL &&
1484 !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1486 /* We need to constrain the LR forces separately,
1487 * because due to the different pre-factor for the SR and LR
1488 * forces in the update algorithm, we have to correct
1489 * the constraint virial for the nstcalclr-1 extra f_lr.
1490 * Constrain only the additional LR part of the force.
1492 /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1497 xp = get_xprime(state, upd);
1499 fac = (nstcalclr - 1)*ir->delta_t*ir->delta_t;
1501 for (i = 0; i < md->homenr; i++)
1503 if (md->cFREEZE != NULL)
1505 gf = md->cFREEZE[i];
1507 for (d = 0; d < DIM; d++)
1509 if ((md->ptype[i] != eptVSite) &&
1510 (md->ptype[i] != eptShell) &&
1511 !ir->opts.nFreeze[gf][d])
1513 xp[i][d] = state->x[i][d] + fac*f_lr[i][d]*md->invmass[i];
1517 constrain(NULL, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, 1.0, md,
1518 state->x, xp, xp, bMolPBC, state->box, state->lambda[efptBONDED], NULL,
1519 NULL, vir_lr_constr, nrnb, econqCoord, ir->epc == epcMTTK, state->veta, state->veta);
1522 /* Add nstcalclr-1 times the LR force to the sum of both forces
1523 * and store the result in forces_lr.
1525 for (i = start; i < nrend; i++)
1527 for (d = 0; d < DIM; d++)
1529 f_lr[i][d] = f[i][d] + (nstcalclr - 1)*f_lr[i][d];
1534 void update_constraints(FILE *fplog,
1536 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1537 t_inputrec *inputrec, /* input record and box stuff */
1538 gmx_ekindata_t *ekind,
1543 rvec force[], /* forces on home particles */
1548 gmx_wallcycle_t wcycle,
1550 gmx_constr_t constr,
1551 gmx_bool bFirstHalf,
1555 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1558 int start, homenr, nrend, i, n, m, g, d;
1560 rvec *vbuf, *xprime = NULL;
1567 if (bFirstHalf && !EI_VV(inputrec->eI))
1572 /* for now, SD update is here -- though it really seems like it
1573 should be reformulated as a velocity verlet method, since it has two parts */
1576 homenr = md->homenr;
1577 nrend = start+homenr;
1579 dt = inputrec->delta_t;
1584 * APPLY CONSTRAINTS:
1587 * When doing PR pressure coupling we have to constrain the
1588 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1589 * it is enough to do this once though, since the relative velocities
1590 * after this will be normal to the bond vector
1595 /* clear out constraints before applying */
1596 clear_mat(vir_part);
1598 xprime = get_xprime(state, upd);
1600 bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1601 bLog = (do_per_step(step, inputrec->nstlog) || bLastStep || (step < 0));
1602 bEner = (do_per_step(step, inputrec->nstenergy) || bLastStep);
1603 /* Constrain the coordinates xprime */
1604 wallcycle_start(wcycle, ewcCONSTR);
1605 if (EI_VV(inputrec->eI) && bFirstHalf)
1607 constrain(NULL, bLog, bEner, constr, idef,
1608 inputrec, ekind, cr, step, 1, 1.0, md,
1609 state->x, state->v, state->v,
1610 bMolPBC, state->box,
1611 state->lambda[efptBONDED], dvdlambda,
1612 NULL, bCalcVir ? &vir_con : NULL, nrnb, econqVeloc,
1613 inputrec->epc == epcMTTK, state->veta, vetanew);
1617 constrain(NULL, bLog, bEner, constr, idef,
1618 inputrec, ekind, cr, step, 1, 1.0, md,
1619 state->x, xprime, NULL,
1620 bMolPBC, state->box,
1621 state->lambda[efptBONDED], dvdlambda,
1622 state->v, bCalcVir ? &vir_con : NULL, nrnb, econqCoord,
1623 inputrec->epc == epcMTTK, state->veta, state->veta);
1625 wallcycle_stop(wcycle, ewcCONSTR);
1629 dump_it_all(fplog, "After Shake",
1630 state->natoms, state->x, xprime, state->v, force);
1634 if (inputrec->eI == eiSD2)
1636 /* A correction factor eph is needed for the SD constraint force */
1637 /* Here we can, unfortunately, not have proper corrections
1638 * for different friction constants, so we use the first one.
1640 for (i = 0; i < DIM; i++)
1642 for (m = 0; m < DIM; m++)
1644 vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1650 m_add(vir_part, vir_con, vir_part);
1654 pr_rvecs(debug, 0, "constraint virial", vir_part, DIM);
1661 if (inputrec->eI == eiSD1 && bDoConstr && !bFirstHalf)
1663 xprime = get_xprime(state, upd);
1665 nth = gmx_omp_nthreads_get(emntUpdate);
1667 #pragma omp parallel for num_threads(nth) schedule(static)
1669 for (th = 0; th < nth; th++)
1671 int start_th, end_th;
1673 start_th = start + ((nrend-start)* th )/nth;
1674 end_th = start + ((nrend-start)*(th+1))/nth;
1676 /* The second part of the SD integration */
1677 do_update_sd1(upd->sd,
1678 start_th, end_th, dt,
1679 inputrec->opts.acc, inputrec->opts.nFreeze,
1680 md->invmass, md->ptype,
1681 md->cFREEZE, md->cACC, md->cTC,
1682 state->x, xprime, state->v, force,
1683 inputrec->opts.ngtc, inputrec->opts.ref_t,
1685 step, inputrec->ld_seed,
1686 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1688 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1692 /* Constrain the coordinates xprime for half a time step */
1693 wallcycle_start(wcycle, ewcCONSTR);
1695 constrain(NULL, bLog, bEner, constr, idef,
1696 inputrec, NULL, cr, step, 1, 0.5, md,
1697 state->x, xprime, NULL,
1698 bMolPBC, state->box,
1699 state->lambda[efptBONDED], dvdlambda,
1700 state->v, NULL, nrnb, econqCoord, FALSE, 0, 0);
1702 wallcycle_stop(wcycle, ewcCONSTR);
1706 if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1708 xprime = get_xprime(state, upd);
1710 nth = gmx_omp_nthreads_get(emntUpdate);
1712 #pragma omp parallel for num_threads(nth) schedule(static)
1713 for (th = 0; th < nth; th++)
1715 int start_th, end_th;
1717 start_th = start + ((nrend-start)* th )/nth;
1718 end_th = start + ((nrend-start)*(th+1))/nth;
1720 /* The second part of the SD integration */
1721 do_update_sd2(upd->sd,
1722 FALSE, start_th, end_th,
1723 inputrec->opts.acc, inputrec->opts.nFreeze,
1724 md->invmass, md->ptype,
1725 md->cFREEZE, md->cACC, md->cTC,
1726 state->x, xprime, state->v, force, state->sd_X,
1727 inputrec->opts.tau_t,
1728 FALSE, step, inputrec->ld_seed,
1729 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1731 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1735 /* Constrain the coordinates xprime */
1736 wallcycle_start(wcycle, ewcCONSTR);
1737 constrain(NULL, bLog, bEner, constr, idef,
1738 inputrec, NULL, cr, step, 1, 1.0, md,
1739 state->x, xprime, NULL,
1740 bMolPBC, state->box,
1741 state->lambda[efptBONDED], dvdlambda,
1742 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
1743 wallcycle_stop(wcycle, ewcCONSTR);
1748 /* We must always unshift after updating coordinates; if we did not shake
1749 x was shifted in do_force */
1751 if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1753 if (graph && (graph->nnodes > 0))
1755 unshift_x(graph, state->box, state->x, upd->xp);
1756 if (TRICLINIC(state->box))
1758 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1762 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1767 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntUpdate)) schedule(static)
1768 for (i = start; i < nrend; i++)
1770 copy_rvec(upd->xp[i], state->x[i]);
1774 dump_it_all(fplog, "After unshift",
1775 state->natoms, state->x, upd->xp, state->v, force);
1777 /* ############# END the update of velocities and positions ######### */
1780 void update_box(FILE *fplog,
1782 t_inputrec *inputrec, /* input record and box stuff */
1785 rvec force[], /* forces on home particles */
1791 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE;
1794 int start, homenr, nrend, i, n, m, g;
1798 homenr = md->homenr;
1799 nrend = start+homenr;
1802 (inputrec->etc == etcNOSEHOOVER) ||
1803 (inputrec->epc == epcPARRINELLORAHMAN) ||
1804 (inputrec->epc == epcMTTK);
1806 dt = inputrec->delta_t;
1810 /* now update boxes */
1811 switch (inputrec->epc)
1815 case (epcBERENDSEN):
1816 berendsen_pscale(inputrec, pcoupl_mu, state->box, state->box_rel,
1817 start, homenr, state->x, md->cFREEZE, nrnb);
1819 case (epcPARRINELLORAHMAN):
1820 /* The box velocities were updated in do_pr_pcoupl in the update
1821 * iteration, but we dont change the box vectors until we get here
1822 * since we need to be able to shift/unshift above.
1824 for (i = 0; i < DIM; i++)
1826 for (m = 0; m <= i; m++)
1828 state->box[i][m] += dt*state->boxv[i][m];
1831 preserve_box_shape(inputrec, state->box_rel, state->box);
1833 /* Scale the coordinates */
1834 for (n = start; (n < start+homenr); n++)
1836 tmvmul_ur0(pcoupl_mu, state->x[n], state->x[n]);
1840 switch (inputrec->epct)
1842 case (epctISOTROPIC):
1843 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1844 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1845 Side length scales as exp(veta*dt) */
1847 msmul(state->box, exp(state->veta*dt), state->box);
1849 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1850 o If we assume isotropic scaling, and box length scaling
1851 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1852 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1853 determinant of B is L^DIM det(M), and the determinant
1854 of dB/dt is (dL/dT)^DIM det (M). veta will be
1855 (det(dB/dT)/det(B))^(1/3). Then since M =
1856 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1858 msmul(state->box, state->veta, state->boxv);
1868 if ((!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))) && scale_tot)
1870 /* The transposes of the scaling matrices are stored,
1871 * therefore we need to reverse the order in the multiplication.
1873 mmul_ur0(*scale_tot, pcoupl_mu, *scale_tot);
1876 if (DEFORM(*inputrec))
1878 deform(upd, start, homenr, state->x, state->box, scale_tot, inputrec, step);
1881 dump_it_all(fplog, "After update",
1882 state->natoms, state->x, upd->xp, state->v, force);
1885 void update_coords(FILE *fplog,
1887 t_inputrec *inputrec, /* input record and box stuff */
1891 rvec *f, /* forces on home particles */
1894 tensor *vir_lr_constr,
1896 gmx_ekindata_t *ekind,
1901 t_commrec *cr, /* these shouldn't be here -- need to think about it */
1903 gmx_constr_t constr,
1906 gmx_bool bNH, bPR, bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1908 real *imass, *imassin;
1911 int start, homenr, nrend, i, j, d, n, m, g;
1912 int blen0, blen1, iatom, jatom, nshake, nsettle, nconstr, nexpand;
1915 rvec *vcom, *xcom, *vall, *xall, *xin, *vin, *forcein, *fall, *xpall, *xprimein, *xprime;
1918 bDoConstr = (NULL != constr);
1920 /* Running the velocity half does nothing except for velocity verlet */
1921 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1922 !EI_VV(inputrec->eI))
1924 gmx_incons("update_coords called for velocity without VV integrator");
1928 homenr = md->homenr;
1929 nrend = start+homenr;
1931 xprime = get_xprime(state, upd);
1933 dt = inputrec->delta_t;
1936 /* We need to update the NMR restraint history when time averaging is used */
1937 if (state->flags & (1<<estDISRE_RM3TAV))
1939 update_disres_history(fcd, &state->hist);
1941 if (state->flags & (1<<estORIRE_DTAV))
1943 update_orires_history(fcd, &state->hist);
1947 bNH = inputrec->etc == etcNOSEHOOVER;
1948 bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1950 if (bDoLR && inputrec->nstcalclr > 1 && !EI_VV(inputrec->eI)) /* get this working with VV? */
1952 /* Store the total force + nstcalclr-1 times the LR force
1953 * in forces_lr, so it can be used in a normal update algorithm
1954 * to produce twin time stepping.
1956 /* is this correct in the new construction? MRS */
1958 inputrec->nstcalclr, constr, inputrec, md, idef, cr,
1959 step, state, bMolPBC,
1960 start, nrend, f, f_lr, vir_lr_constr, nrnb);
1968 /* ############# START The update of velocities and positions ######### */
1970 dump_it_all(fplog, "Before update",
1971 state->natoms, state->x, xprime, state->v, force);
1973 if (inputrec->eI == eiSD2)
1975 check_sd2_work_data_allocation(upd->sd, nrend);
1977 do_update_sd2_Tconsts(upd->sd,
1978 inputrec->opts.ngtc,
1979 inputrec->opts.tau_t,
1980 inputrec->opts.ref_t);
1982 if (inputrec->eI == eiBD)
1984 do_update_bd_Tconsts(dt, inputrec->bd_fric,
1985 inputrec->opts.ngtc, inputrec->opts.ref_t,
1989 nth = gmx_omp_nthreads_get(emntUpdate);
1991 #pragma omp parallel for num_threads(nth) schedule(static) private(alpha)
1992 for (th = 0; th < nth; th++)
1994 int start_th, end_th;
1996 start_th = start + ((nrend-start)* th )/nth;
1997 end_th = start + ((nrend-start)*(th+1))/nth;
1999 switch (inputrec->eI)
2002 if (ekind->cosacc.cos_accel == 0)
2004 do_update_md(start_th, end_th, dt,
2005 ekind->tcstat, state->nosehoover_vxi,
2006 ekind->bNEMD, ekind->grpstat, inputrec->opts.acc,
2007 inputrec->opts.nFreeze,
2008 md->invmass, md->ptype,
2009 md->cFREEZE, md->cACC, md->cTC,
2010 state->x, xprime, state->v, force, M,
2015 do_update_visc(start_th, end_th, dt,
2016 ekind->tcstat, state->nosehoover_vxi,
2017 md->invmass, md->ptype,
2018 md->cTC, state->x, xprime, state->v, force, M,
2020 ekind->cosacc.cos_accel,
2026 /* With constraints, the SD1 update is done in 2 parts */
2027 do_update_sd1(upd->sd,
2028 start_th, end_th, dt,
2029 inputrec->opts.acc, inputrec->opts.nFreeze,
2030 md->invmass, md->ptype,
2031 md->cFREEZE, md->cACC, md->cTC,
2032 state->x, xprime, state->v, force,
2033 inputrec->opts.ngtc, inputrec->opts.ref_t,
2035 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2038 /* The SD2 update is always done in 2 parts,
2039 * because an extra constraint step is needed
2041 do_update_sd2(upd->sd,
2042 bInitStep, start_th, end_th,
2043 inputrec->opts.acc, inputrec->opts.nFreeze,
2044 md->invmass, md->ptype,
2045 md->cFREEZE, md->cACC, md->cTC,
2046 state->x, xprime, state->v, force, state->sd_X,
2047 inputrec->opts.tau_t,
2048 TRUE, step, inputrec->ld_seed,
2049 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2052 do_update_bd(start_th, end_th, dt,
2053 inputrec->opts.nFreeze, md->invmass, md->ptype,
2054 md->cFREEZE, md->cTC,
2055 state->x, xprime, state->v, force,
2058 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2062 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
2067 do_update_vv_vel(start_th, end_th, dt,
2068 inputrec->opts.acc, inputrec->opts.nFreeze,
2069 md->invmass, md->ptype,
2070 md->cFREEZE, md->cACC,
2072 (bNH || bPR), state->veta, alpha);
2075 do_update_vv_pos(start_th, end_th, dt,
2076 inputrec->opts.nFreeze,
2077 md->ptype, md->cFREEZE,
2078 state->x, xprime, state->v,
2079 (bNH || bPR), state->veta);
2084 gmx_fatal(FARGS, "Don't know how to update coordinates");
2092 void correct_ekin(FILE *log, int start, int end, rvec v[], rvec vcm, real mass[],
2093 real tmass, tensor ekin)
2096 * This is a debugging routine. It should not be called for production code
2098 * The kinetic energy should calculated according to:
2099 * Ekin = 1/2 m (v-vcm)^2
2100 * However the correction is not always applied, since vcm may not be
2101 * known in time and we compute
2102 * Ekin' = 1/2 m v^2 instead
2103 * This can be corrected afterwards by computing
2104 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
2106 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
2113 /* Local particles */
2116 /* Processor dependent part. */
2118 for (i = start; (i < end); i++)
2122 for (j = 0; (j < DIM); j++)
2128 svmul(1/tmass, vcm, vcm);
2129 svmul(0.5, vcm, hvcm);
2131 for (j = 0; (j < DIM); j++)
2133 for (k = 0; (k < DIM); k++)
2135 dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
2138 pr_rvecs(log, 0, "dekin", dekin, DIM);
2139 pr_rvecs(log, 0, " ekin", ekin, DIM);
2140 fprintf(log, "dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
2141 trace(dekin), trace(ekin), vcm[XX], vcm[YY], vcm[ZZ]);
2142 fprintf(log, "mv = (%8.4f %8.4f %8.4f)\n",
2143 mv[XX], mv[YY], mv[ZZ]);
2146 extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_int64_t step, const t_commrec *cr,
2147 t_mdatoms *md, t_state *state, gmx_update_t upd, gmx_constr_t constr)
2151 real rate = (ir->delta_t)/ir->opts.tau_t[0];
2153 if (ir->etc == etcANDERSEN && constr != NULL)
2155 gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
2158 /* proceed with andersen if 1) it's fixed probability per
2159 particle andersen or 2) it's massive andersen and it's tau_t/dt */
2160 if ((ir->etc == etcANDERSEN) || do_per_step(step, (int)(1.0/rate)))
2162 andersen_tcoupl(ir, step, cr, md, state, rate,
2163 upd->sd->randomize_group, upd->sd->boltzfac);