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.
42 #include "gromacs/legacyheaders/types/commrec.h"
43 #include "gromacs/legacyheaders/typedefs.h"
44 #include "gromacs/legacyheaders/nrnb.h"
45 #include "gromacs/math/units.h"
46 #include "gromacs/legacyheaders/macros.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/legacyheaders/update.h"
49 #include "gromacs/random/random.h"
50 #include "gromacs/legacyheaders/tgroup.h"
51 #include "gromacs/legacyheaders/force.h"
52 #include "gromacs/legacyheaders/names.h"
53 #include "gromacs/legacyheaders/txtdump.h"
54 #include "gromacs/legacyheaders/mdrun.h"
55 #include "gromacs/legacyheaders/constr.h"
56 #include "gromacs/legacyheaders/disre.h"
57 #include "gromacs/legacyheaders/orires.h"
58 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
60 #include "gromacs/fileio/confio.h"
61 #include "gromacs/pbcutil/mshift.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/pulling/pull.h"
64 #include "gromacs/timing/wallcycle.h"
65 #include "gromacs/utility/futil.h"
66 #include "gromacs/utility/gmxomp.h"
67 #include "gromacs/utility/smalloc.h"
69 /*For debugging, start at v(-dt/2) for velolcity verlet -- uncomment next line */
70 /*#define STARTFROMDT2*/
94 gmx_sd_sigma_t *sdsig;
97 /* andersen temperature control stuff */
98 gmx_bool *randomize_group;
102 typedef struct gmx_update
105 /* xprime for constraint algorithms */
109 /* Variables for the deform algorithm */
110 gmx_int64_t deformref_step;
111 matrix deformref_box;
115 static void do_update_md(int start, int nrend, double dt,
116 t_grp_tcstat *tcstat,
118 gmx_bool bNEMD, t_grp_acc *gstat, rvec accel[],
121 unsigned short ptype[], unsigned short cFREEZE[],
122 unsigned short cACC[], unsigned short cTC[],
123 rvec x[], rvec xprime[], rvec v[],
125 gmx_bool bNH, gmx_bool bPR)
128 int gf = 0, ga = 0, gt = 0;
130 real vn, vv, va, vb, vnrel;
136 /* Update with coupling to extended ensembles, used for
137 * Nose-Hoover and Parrinello-Rahman coupling
138 * Nose-Hoover uses the reversible leap-frog integrator from
139 * Holian et al. Phys Rev E 52(3) : 2338, 1995
141 for (n = start; n < nrend; n++)
156 lg = tcstat[gt].lambda;
161 rvec_sub(v[n], gstat[ga].u, vrel);
163 for (d = 0; d < DIM; d++)
165 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
167 vnrel = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
168 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
169 /* do not scale the mean velocities u */
170 vn = gstat[ga].u[d] + accel[ga][d]*dt + vnrel;
172 xprime[n][d] = x[n][d]+vn*dt;
177 xprime[n][d] = x[n][d];
182 else if (cFREEZE != NULL ||
183 nFreeze[0][XX] || nFreeze[0][YY] || nFreeze[0][ZZ] ||
186 /* Update with Berendsen/v-rescale coupling and freeze or NEMD */
187 for (n = start; n < nrend; n++)
189 w_dt = invmass[n]*dt;
202 lg = tcstat[gt].lambda;
204 for (d = 0; d < DIM; d++)
207 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
209 vv = lg*vn + f[n][d]*w_dt;
211 /* do not scale the mean velocities u */
213 va = vv + accel[ga][d]*dt;
214 vb = va + (1.0-lg)*u;
216 xprime[n][d] = x[n][d]+vb*dt;
221 xprime[n][d] = x[n][d];
228 /* Plain update with Berendsen/v-rescale coupling */
229 for (n = start; n < nrend; n++)
231 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
233 w_dt = invmass[n]*dt;
238 lg = tcstat[gt].lambda;
240 for (d = 0; d < DIM; d++)
242 vn = lg*v[n][d] + f[n][d]*w_dt;
244 xprime[n][d] = x[n][d] + vn*dt;
249 for (d = 0; d < DIM; d++)
252 xprime[n][d] = x[n][d];
259 static void do_update_vv_vel(int start, int nrend, double dt,
260 rvec accel[], ivec nFreeze[], real invmass[],
261 unsigned short ptype[], unsigned short cFREEZE[],
262 unsigned short cACC[], rvec v[], rvec f[],
263 gmx_bool bExtended, real veta, real alpha)
268 real u, vn, vv, va, vb, vnrel;
274 g = 0.25*dt*veta*alpha;
276 mv2 = series_sinhx(g);
283 for (n = start; n < nrend; n++)
285 w_dt = invmass[n]*dt;
295 for (d = 0; d < DIM; d++)
297 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
299 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
307 } /* do_update_vv_vel */
309 static void do_update_vv_pos(int start, int nrend, double dt,
311 unsigned short ptype[], unsigned short cFREEZE[],
312 rvec x[], rvec xprime[], rvec v[],
313 gmx_bool bExtended, real veta)
320 /* Would it make more sense if Parrinello-Rahman was put here? */
325 mr2 = series_sinhx(g);
333 for (n = start; n < nrend; n++)
341 for (d = 0; d < DIM; d++)
343 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
345 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
349 xprime[n][d] = x[n][d];
353 } /* do_update_vv_pos */
355 static void do_update_visc(int start, int nrend, double dt,
356 t_grp_tcstat *tcstat,
359 unsigned short ptype[], unsigned short cTC[],
360 rvec x[], rvec xprime[], rvec v[],
361 rvec f[], matrix M, matrix box, real
362 cos_accel, real vcos,
363 gmx_bool bNH, gmx_bool bPR)
368 real lg, vxi = 0, vv;
373 fac = 2*M_PI/(box[ZZ][ZZ]);
377 /* Update with coupling to extended ensembles, used for
378 * Nose-Hoover and Parrinello-Rahman coupling
380 for (n = start; n < nrend; n++)
387 lg = tcstat[gt].lambda;
388 cosz = cos(fac*x[n][ZZ]);
390 copy_rvec(v[n], vrel);
398 for (d = 0; d < DIM; d++)
402 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
404 vn = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
405 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
408 vn += vc + dt*cosz*cos_accel;
411 xprime[n][d] = x[n][d]+vn*dt;
415 xprime[n][d] = x[n][d];
422 /* Classic version of update, used with berendsen coupling */
423 for (n = start; n < nrend; n++)
425 w_dt = invmass[n]*dt;
430 lg = tcstat[gt].lambda;
431 cosz = cos(fac*x[n][ZZ]);
433 for (d = 0; d < DIM; d++)
437 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
442 /* Do not scale the cosine velocity profile */
443 vv = vc + lg*(vn - vc + f[n][d]*w_dt);
444 /* Add the cosine accelaration profile */
445 vv += dt*cosz*cos_accel;
449 vv = lg*(vn + f[n][d]*w_dt);
452 xprime[n][d] = x[n][d]+vv*dt;
457 xprime[n][d] = x[n][d];
464 static gmx_stochd_t *init_stochd(t_inputrec *ir)
473 ngtc = ir->opts.ngtc;
477 snew(sd->bd_rf, ngtc);
479 else if (EI_SD(ir->eI))
482 snew(sd->sdsig, ngtc);
485 for (n = 0; n < ngtc; n++)
487 if (ir->opts.tau_t[n] > 0)
489 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
490 sdc[n].eph = exp(sdc[n].gdt/2);
491 sdc[n].emh = exp(-sdc[n].gdt/2);
492 sdc[n].em = exp(-sdc[n].gdt);
496 /* No friction and noise on this group */
502 if (sdc[n].gdt >= 0.05)
504 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
505 - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
506 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
507 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
512 /* Seventh order expansions for small y */
513 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
514 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))));
515 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
519 fprintf(debug, "SD const tc-grp %d: b %g c %g d %g\n",
520 n, sdc[n].b, sdc[n].c, sdc[n].d);
524 else if (ETC_ANDERSEN(ir->etc))
533 snew(sd->randomize_group, ngtc);
534 snew(sd->boltzfac, ngtc);
536 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
537 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
539 for (n = 0; n < ngtc; n++)
541 reft = max(0.0, opts->ref_t[n]);
542 if ((opts->tau_t[n] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
544 sd->randomize_group[n] = TRUE;
545 sd->boltzfac[n] = BOLTZ*opts->ref_t[n];
549 sd->randomize_group[n] = FALSE;
556 gmx_update_t init_update(t_inputrec *ir)
562 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
564 upd->sd = init_stochd(ir);
573 static void do_update_sd1(gmx_stochd_t *sd,
574 int start, int nrend, double dt,
575 rvec accel[], ivec nFreeze[],
576 real invmass[], unsigned short ptype[],
577 unsigned short cFREEZE[], unsigned short cACC[],
578 unsigned short cTC[],
579 rvec x[], rvec xprime[], rvec v[], rvec f[],
580 int ngtc, real ref_t[],
582 gmx_bool bFirstHalfConstr,
583 gmx_int64_t step, int seed, int* gatindex)
588 int gf = 0, ga = 0, gt = 0;
595 for (n = 0; n < ngtc; n++)
598 /* The mass is encounted for later, since this differs per atom */
599 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;
609 ism = sqrt(invmass[n]);
623 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
625 for (d = 0; d < DIM; d++)
627 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
631 sd_V = ism*sig[gt].V*rnd[d];
632 vn = v[n][d] + (invmass[n]*f[n][d] + accel[ga][d])*dt;
633 v[n][d] = vn*sdc[gt].em + sd_V;
634 /* Here we include half of the friction+noise
635 * update of v into the integration of x.
637 xprime[n][d] = x[n][d] + 0.5*(vn + v[n][d])*dt;
642 xprime[n][d] = x[n][d];
649 /* We do have constraints */
650 if (bFirstHalfConstr)
652 /* First update without friction and noise */
655 for (n = start; n < nrend; n++)
672 for (d = 0; d < DIM; d++)
674 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
676 v[n][d] = v[n][d] + (im*f[n][d] + accel[ga][d])*dt;
677 xprime[n][d] = x[n][d] + v[n][d]*dt;
682 xprime[n][d] = x[n][d];
689 /* Update friction and noise only */
690 for (n = start; n < nrend; n++)
693 int ng = gatindex ? gatindex[n] : n;
695 ism = sqrt(invmass[n]);
709 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
711 for (d = 0; d < DIM; d++)
713 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
717 sd_V = ism*sig[gt].V*rnd[d];
719 v[n][d] = vn*sdc[gt].em + sd_V;
720 /* Add the friction and noise contribution only */
721 xprime[n][d] = xprime[n][d] + 0.5*(v[n][d] - vn)*dt;
729 static void check_sd2_work_data_allocation(gmx_stochd_t *sd, int nrend)
731 if (nrend > sd->sd_V_nalloc)
733 sd->sd_V_nalloc = over_alloc_dd(nrend);
734 srenew(sd->sd_V, sd->sd_V_nalloc);
738 static void do_update_sd2_Tconsts(gmx_stochd_t *sd,
743 /* This is separated from the update below, because it is single threaded */
752 for (gt = 0; gt < ngtc; gt++)
754 kT = BOLTZ*ref_t[gt];
755 /* The mass is encounted for later, since this differs per atom */
756 sig[gt].V = sqrt(kT*(1-sdc[gt].em));
757 sig[gt].X = sqrt(kT*sqr(tau_t[gt])*sdc[gt].c);
758 sig[gt].Yv = sqrt(kT*sdc[gt].b/sdc[gt].c);
759 sig[gt].Yx = sqrt(kT*sqr(tau_t[gt])*sdc[gt].b/(1-sdc[gt].em));
763 static void do_update_sd2(gmx_stochd_t *sd,
765 int start, int nrend,
766 rvec accel[], ivec nFreeze[],
767 real invmass[], unsigned short ptype[],
768 unsigned short cFREEZE[], unsigned short cACC[],
769 unsigned short cTC[],
770 rvec x[], rvec xprime[], rvec v[], rvec f[],
773 gmx_bool bFirstHalf, gmx_int64_t step, int seed,
778 /* The random part of the velocity update, generated in the first
779 * half of the update, needs to be remembered for the second half.
783 int gf = 0, ga = 0, gt = 0;
784 real vn = 0, Vmh, Xmh;
792 for (n = start; n < nrend; n++)
794 real rnd[6], rndi[3];
795 ng = gatindex ? gatindex[n] : n;
796 ism = sqrt(invmass[n]);
810 gmx_rng_cycle_6gaussian_table(step*2+(bFirstHalf ? 1 : 2), ng, seed, RND_SEED_UPDATE, rnd);
813 gmx_rng_cycle_3gaussian_table(step*2, ng, seed, RND_SEED_UPDATE, rndi);
815 for (d = 0; d < DIM; d++)
821 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
827 sd_X[n][d] = ism*sig[gt].X*rndi[d];
829 Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
830 + ism*sig[gt].Yv*rnd[d*2];
831 sd_V[n][d] = ism*sig[gt].V*rnd[d*2+1];
833 v[n][d] = vn*sdc[gt].em
834 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
835 + sd_V[n][d] - sdc[gt].em*Vmh;
837 xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
841 /* Correct the velocities for the constraints.
842 * This operation introduces some inaccuracy,
843 * since the velocity is determined from differences in coordinates.
846 (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
848 Xmh = sd_V[n][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
849 + ism*sig[gt].Yx*rnd[d*2];
850 sd_X[n][d] = ism*sig[gt].X*rnd[d*2+1];
852 xprime[n][d] += sd_X[n][d] - Xmh;
861 xprime[n][d] = x[n][d];
868 static void do_update_bd_Tconsts(double dt, real friction_coefficient,
869 int ngtc, const real ref_t[],
872 /* This is separated from the update below, because it is single threaded */
875 if (friction_coefficient != 0)
877 for (gt = 0; gt < ngtc; gt++)
879 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]/(friction_coefficient*dt));
884 for (gt = 0; gt < ngtc; gt++)
886 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]);
891 static void do_update_bd(int start, int nrend, double dt,
893 real invmass[], unsigned short ptype[],
894 unsigned short cFREEZE[], unsigned short cTC[],
895 rvec x[], rvec xprime[], rvec v[],
896 rvec f[], real friction_coefficient,
897 real *rf, gmx_int64_t step, int seed,
900 /* note -- these appear to be full step velocities . . . */
906 if (friction_coefficient != 0)
908 invfr = 1.0/friction_coefficient;
911 for (n = start; (n < nrend); n++)
914 int ng = gatindex ? gatindex[n] : n;
924 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
925 for (d = 0; (d < DIM); d++)
927 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
929 if (friction_coefficient != 0)
931 vn = invfr*f[n][d] + rf[gt]*rnd[d];
935 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
936 vn = 0.5*invmass[n]*f[n][d]*dt
937 + sqrt(0.5*invmass[n])*rf[gt]*rnd[d];
941 xprime[n][d] = x[n][d]+vn*dt;
946 xprime[n][d] = x[n][d];
952 static void dump_it_all(FILE gmx_unused *fp, const char gmx_unused *title,
953 int gmx_unused natoms, rvec gmx_unused x[], rvec gmx_unused xp[],
954 rvec gmx_unused v[], rvec gmx_unused f[])
959 fprintf(fp, "%s\n", title);
960 pr_rvecs(fp, 0, "x", x, natoms);
961 pr_rvecs(fp, 0, "xp", xp, natoms);
962 pr_rvecs(fp, 0, "v", v, natoms);
963 pr_rvecs(fp, 0, "f", f, natoms);
968 static void calc_ke_part_normal(rvec v[], t_grpopts *opts, t_mdatoms *md,
969 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel,
970 gmx_bool bSaveEkinOld)
973 t_grp_tcstat *tcstat = ekind->tcstat;
974 t_grp_acc *grpstat = ekind->grpstat;
977 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
978 an option, but not supported now. Additionally, if we are doing iterations.
979 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
980 bSavEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't copy over the ekinh_old.
981 If FALSE, we overrwrite it.
984 /* group velocities are calculated in update_ekindata and
985 * accumulated in acumulate_groups.
986 * Now the partial global and groups ekin.
988 for (g = 0; (g < opts->ngtc); g++)
993 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
997 clear_mat(tcstat[g].ekinf);
1001 clear_mat(tcstat[g].ekinh);
1005 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
1008 ekind->dekindl_old = ekind->dekindl;
1010 nthread = gmx_omp_nthreads_get(emntUpdate);
1012 #pragma omp parallel for num_threads(nthread) schedule(static)
1013 for (thread = 0; thread < nthread; thread++)
1015 int start_t, end_t, n;
1023 start_t = ((thread+0)*md->homenr)/nthread;
1024 end_t = ((thread+1)*md->homenr)/nthread;
1026 ekin_sum = ekind->ekin_work[thread];
1027 dekindl_sum = ekind->dekindl_work[thread];
1029 for (gt = 0; gt < opts->ngtc; gt++)
1031 clear_mat(ekin_sum[gt]);
1037 for (n = start_t; n < end_t; n++)
1047 hm = 0.5*md->massT[n];
1049 for (d = 0; (d < DIM); d++)
1051 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
1053 for (d = 0; (d < DIM); d++)
1055 for (m = 0; (m < DIM); m++)
1057 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1058 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
1061 if (md->nMassPerturbed && md->bPerturbed[n])
1064 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1070 for (thread = 0; thread < nthread; thread++)
1072 for (g = 0; g < opts->ngtc; g++)
1076 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
1081 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1086 ekind->dekindl += *ekind->dekindl_work[thread];
1089 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1092 static void calc_ke_part_visc(matrix box, rvec x[], rvec v[],
1093 t_grpopts *opts, t_mdatoms *md,
1094 gmx_ekindata_t *ekind,
1095 t_nrnb *nrnb, gmx_bool bEkinAveVel)
1097 int start = 0, homenr = md->homenr;
1098 int g, d, n, m, gt = 0;
1101 t_grp_tcstat *tcstat = ekind->tcstat;
1102 t_cos_acc *cosacc = &(ekind->cosacc);
1107 for (g = 0; g < opts->ngtc; g++)
1109 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1110 clear_mat(ekind->tcstat[g].ekinh);
1112 ekind->dekindl_old = ekind->dekindl;
1114 fac = 2*M_PI/box[ZZ][ZZ];
1117 for (n = start; n < start+homenr; n++)
1123 hm = 0.5*md->massT[n];
1125 /* Note that the times of x and v differ by half a step */
1126 /* MRS -- would have to be changed for VV */
1127 cosz = cos(fac*x[n][ZZ]);
1128 /* Calculate the amplitude of the new velocity profile */
1129 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1131 copy_rvec(v[n], v_corrt);
1132 /* Subtract the profile for the kinetic energy */
1133 v_corrt[XX] -= cosz*cosacc->vcos;
1134 for (d = 0; (d < DIM); d++)
1136 for (m = 0; (m < DIM); m++)
1138 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1141 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1145 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1149 if (md->nPerturbed && md->bPerturbed[n])
1151 /* The minus sign here might be confusing.
1152 * The kinetic contribution from dH/dl doesn't come from
1153 * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
1154 * where p are the momenta. The difference is only a minus sign.
1156 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1159 ekind->dekindl = dekindl;
1160 cosacc->mvcos = mvcos;
1162 inc_nrnb(nrnb, eNR_EKIN, homenr);
1165 void calc_ke_part(t_state *state, t_grpopts *opts, t_mdatoms *md,
1166 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
1168 if (ekind->cosacc.cos_accel == 0)
1170 calc_ke_part_normal(state->v, opts, md, ekind, nrnb, bEkinAveVel, bSaveEkinOld);
1174 calc_ke_part_visc(state->box, state->x, state->v, opts, md, ekind, nrnb, bEkinAveVel);
1178 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1180 ekinstate->ekin_n = ir->opts.ngtc;
1181 snew(ekinstate->ekinh, ekinstate->ekin_n);
1182 snew(ekinstate->ekinf, ekinstate->ekin_n);
1183 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1184 snew(ekinstate->ekinscalef_nhc, ekinstate->ekin_n);
1185 snew(ekinstate->ekinscaleh_nhc, ekinstate->ekin_n);
1186 snew(ekinstate->vscale_nhc, ekinstate->ekin_n);
1187 ekinstate->dekindl = 0;
1188 ekinstate->mvcos = 0;
1191 void update_ekinstate(ekinstate_t *ekinstate, gmx_ekindata_t *ekind)
1195 for (i = 0; i < ekinstate->ekin_n; i++)
1197 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1198 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1199 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1200 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1201 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1202 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1205 copy_mat(ekind->ekin, ekinstate->ekin_total);
1206 ekinstate->dekindl = ekind->dekindl;
1207 ekinstate->mvcos = ekind->cosacc.mvcos;
1211 void restore_ekinstate_from_state(t_commrec *cr,
1212 gmx_ekindata_t *ekind, ekinstate_t *ekinstate)
1218 for (i = 0; i < ekinstate->ekin_n; i++)
1220 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1221 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1222 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1223 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1224 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1225 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1228 copy_mat(ekinstate->ekin_total, ekind->ekin);
1230 ekind->dekindl = ekinstate->dekindl;
1231 ekind->cosacc.mvcos = ekinstate->mvcos;
1232 n = ekinstate->ekin_n;
1237 gmx_bcast(sizeof(n), &n, cr);
1238 for (i = 0; i < n; i++)
1240 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1241 ekind->tcstat[i].ekinh[0], cr);
1242 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1243 ekind->tcstat[i].ekinf[0], cr);
1244 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1245 ekind->tcstat[i].ekinh_old[0], cr);
1247 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1248 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1249 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1250 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1251 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1252 &(ekind->tcstat[i].vscale_nhc), cr);
1254 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1255 ekind->ekin[0], cr);
1257 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1258 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1262 void set_deform_reference_box(gmx_update_t upd, gmx_int64_t step, matrix box)
1264 upd->deformref_step = step;
1265 copy_mat(box, upd->deformref_box);
1268 static void deform(gmx_update_t upd,
1269 int start, int homenr, rvec x[], matrix box, matrix *scale_tot,
1270 const t_inputrec *ir, gmx_int64_t step)
1272 matrix bnew, invbox, mu;
1276 elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1277 copy_mat(box, bnew);
1278 for (i = 0; i < DIM; i++)
1280 for (j = 0; j < DIM; j++)
1282 if (ir->deform[i][j] != 0)
1285 upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1289 /* We correct the off-diagonal elements,
1290 * which can grow indefinitely during shearing,
1291 * so the shifts do not get messed up.
1293 for (i = 1; i < DIM; i++)
1295 for (j = i-1; j >= 0; j--)
1297 while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1299 rvec_dec(bnew[i], bnew[j]);
1301 while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1303 rvec_inc(bnew[i], bnew[j]);
1307 m_inv_ur0(box, invbox);
1308 copy_mat(bnew, box);
1309 mmul_ur0(box, invbox, mu);
1311 for (i = start; i < start+homenr; i++)
1313 x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1314 x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1315 x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1317 if (scale_tot != NULL)
1319 /* The transposes of the scaling matrices are stored,
1320 * so we need to do matrix multiplication in the inverse order.
1322 mmul_ur0(*scale_tot, mu, *scale_tot);
1326 void update_tcouple(gmx_int64_t step,
1327 t_inputrec *inputrec,
1329 gmx_ekindata_t *ekind,
1334 gmx_bool bTCouple = FALSE;
1336 int i, start, end, homenr, offset;
1338 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1339 if (inputrec->etc != etcNO &&
1340 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec)))
1342 /* We should only couple after a step where energies were determined (for leapfrog versions)
1343 or the step energies are determined, for velocity verlet versions */
1345 if (EI_VV(inputrec->eI))
1353 bTCouple = (inputrec->nsttcouple == 1 ||
1354 do_per_step(step+inputrec->nsttcouple-offset,
1355 inputrec->nsttcouple));
1360 dttc = inputrec->nsttcouple*inputrec->delta_t;
1362 switch (inputrec->etc)
1367 berendsen_tcoupl(inputrec, ekind, dttc);
1370 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1371 state->nosehoover_xi, state->nosehoover_vxi, MassQ);
1374 vrescale_tcoupl(inputrec, step, ekind, dttc,
1375 state->therm_integral);
1378 /* rescale in place here */
1379 if (EI_VV(inputrec->eI))
1381 rescale_velocities(ekind, md, 0, md->homenr, state->v);
1386 /* Set the T scaling lambda to 1 to have no scaling */
1387 for (i = 0; (i < inputrec->opts.ngtc); i++)
1389 ekind->tcstat[i].lambda = 1.0;
1394 void update_pcouple(FILE *fplog,
1396 t_inputrec *inputrec,
1402 gmx_bool bPCouple = FALSE;
1406 /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1407 if (inputrec->epc != epcNO && (!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))))
1409 /* We should only couple after a step where energies were determined */
1410 bPCouple = (inputrec->nstpcouple == 1 ||
1411 do_per_step(step+inputrec->nstpcouple-1,
1412 inputrec->nstpcouple));
1415 clear_mat(pcoupl_mu);
1416 for (i = 0; i < DIM; i++)
1418 pcoupl_mu[i][i] = 1.0;
1425 dtpc = inputrec->nstpcouple*inputrec->delta_t;
1427 switch (inputrec->epc)
1429 /* We can always pcoupl, even if we did not sum the energies
1430 * the previous step, since state->pres_prev is only updated
1431 * when the energies have been summed.
1435 case (epcBERENDSEN):
1438 berendsen_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev, state->box,
1442 case (epcPARRINELLORAHMAN):
1443 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1444 state->box, state->box_rel, state->boxv,
1445 M, pcoupl_mu, bInitStep);
1453 static rvec *get_xprime(const t_state *state, gmx_update_t upd)
1455 if (state->nalloc > upd->xp_nalloc)
1457 upd->xp_nalloc = state->nalloc;
1458 srenew(upd->xp, upd->xp_nalloc);
1464 static void combine_forces(gmx_update_t upd,
1466 gmx_constr_t constr,
1467 t_inputrec *ir, t_mdatoms *md, t_idef *idef,
1470 t_state *state, gmx_bool bMolPBC,
1471 int start, int nrend,
1472 rvec f[], rvec f_lr[],
1473 tensor *vir_lr_constr,
1478 /* f contains the short-range forces + the long range forces
1479 * which are stored separately in f_lr.
1482 if (constr != NULL && vir_lr_constr != NULL &&
1483 !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1485 /* We need to constrain the LR forces separately,
1486 * because due to the different pre-factor for the SR and LR
1487 * forces in the update algorithm, we have to correct
1488 * the constraint virial for the nstcalclr-1 extra f_lr.
1489 * Constrain only the additional LR part of the force.
1491 /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1496 xp = get_xprime(state, upd);
1498 fac = (nstcalclr - 1)*ir->delta_t*ir->delta_t;
1500 for (i = 0; i < md->homenr; i++)
1502 if (md->cFREEZE != NULL)
1504 gf = md->cFREEZE[i];
1506 for (d = 0; d < DIM; d++)
1508 if ((md->ptype[i] != eptVSite) &&
1509 (md->ptype[i] != eptShell) &&
1510 !ir->opts.nFreeze[gf][d])
1512 xp[i][d] = state->x[i][d] + fac*f_lr[i][d]*md->invmass[i];
1516 constrain(NULL, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, 1.0, md,
1517 state->x, xp, xp, bMolPBC, state->box, state->lambda[efptBONDED], NULL,
1518 NULL, vir_lr_constr, nrnb, econqCoord, ir->epc == epcMTTK, state->veta, state->veta);
1521 /* Add nstcalclr-1 times the LR force to the sum of both forces
1522 * and store the result in forces_lr.
1524 for (i = start; i < nrend; i++)
1526 for (d = 0; d < DIM; d++)
1528 f_lr[i][d] = f[i][d] + (nstcalclr - 1)*f_lr[i][d];
1533 void update_constraints(FILE *fplog,
1535 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1536 t_inputrec *inputrec, /* input record and box stuff */
1537 gmx_ekindata_t *ekind,
1542 rvec force[], /* forces on home particles */
1547 gmx_wallcycle_t wcycle,
1549 gmx_constr_t constr,
1550 gmx_bool bFirstHalf,
1554 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1557 int start, homenr, nrend, i, n, m, g, d;
1559 rvec *vbuf, *xprime = NULL;
1566 if (bFirstHalf && !EI_VV(inputrec->eI))
1571 /* for now, SD update is here -- though it really seems like it
1572 should be reformulated as a velocity verlet method, since it has two parts */
1575 homenr = md->homenr;
1576 nrend = start+homenr;
1578 dt = inputrec->delta_t;
1583 * APPLY CONSTRAINTS:
1586 * When doing PR pressure coupling we have to constrain the
1587 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1588 * it is enough to do this once though, since the relative velocities
1589 * after this will be normal to the bond vector
1594 /* clear out constraints before applying */
1595 clear_mat(vir_part);
1597 xprime = get_xprime(state, upd);
1599 bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1600 bLog = (do_per_step(step, inputrec->nstlog) || bLastStep || (step < 0));
1601 bEner = (do_per_step(step, inputrec->nstenergy) || bLastStep);
1602 /* Constrain the coordinates xprime */
1603 wallcycle_start(wcycle, ewcCONSTR);
1604 if (EI_VV(inputrec->eI) && bFirstHalf)
1606 constrain(NULL, bLog, bEner, constr, idef,
1607 inputrec, ekind, cr, step, 1, 1.0, md,
1608 state->x, state->v, state->v,
1609 bMolPBC, state->box,
1610 state->lambda[efptBONDED], dvdlambda,
1611 NULL, bCalcVir ? &vir_con : NULL, nrnb, econqVeloc,
1612 inputrec->epc == epcMTTK, state->veta, vetanew);
1616 constrain(NULL, bLog, bEner, constr, idef,
1617 inputrec, ekind, cr, step, 1, 1.0, md,
1618 state->x, xprime, NULL,
1619 bMolPBC, state->box,
1620 state->lambda[efptBONDED], dvdlambda,
1621 state->v, bCalcVir ? &vir_con : NULL, nrnb, econqCoord,
1622 inputrec->epc == epcMTTK, state->veta, state->veta);
1624 wallcycle_stop(wcycle, ewcCONSTR);
1628 dump_it_all(fplog, "After Shake",
1629 state->natoms, state->x, xprime, state->v, force);
1633 if (inputrec->eI == eiSD2)
1635 /* A correction factor eph is needed for the SD constraint force */
1636 /* Here we can, unfortunately, not have proper corrections
1637 * for different friction constants, so we use the first one.
1639 for (i = 0; i < DIM; i++)
1641 for (m = 0; m < DIM; m++)
1643 vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1649 m_add(vir_part, vir_con, vir_part);
1653 pr_rvecs(debug, 0, "constraint virial", vir_part, DIM);
1660 if (inputrec->eI == eiSD1 && bDoConstr && !bFirstHalf)
1662 wallcycle_start(wcycle, ewcUPDATE);
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);
1689 wallcycle_stop(wcycle, ewcUPDATE);
1693 /* Constrain the coordinates xprime for half a time step */
1694 wallcycle_start(wcycle, ewcCONSTR);
1696 constrain(NULL, bLog, bEner, constr, idef,
1697 inputrec, NULL, cr, step, 1, 0.5, md,
1698 state->x, xprime, NULL,
1699 bMolPBC, state->box,
1700 state->lambda[efptBONDED], dvdlambda,
1701 state->v, NULL, nrnb, econqCoord, FALSE, 0, 0);
1703 wallcycle_stop(wcycle, ewcCONSTR);
1707 if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1709 wallcycle_start(wcycle, ewcUPDATE);
1710 xprime = get_xprime(state, upd);
1712 nth = gmx_omp_nthreads_get(emntUpdate);
1714 #pragma omp parallel for num_threads(nth) schedule(static)
1715 for (th = 0; th < nth; th++)
1717 int start_th, end_th;
1719 start_th = start + ((nrend-start)* th )/nth;
1720 end_th = start + ((nrend-start)*(th+1))/nth;
1722 /* The second part of the SD integration */
1723 do_update_sd2(upd->sd,
1724 FALSE, start_th, end_th,
1725 inputrec->opts.acc, inputrec->opts.nFreeze,
1726 md->invmass, md->ptype,
1727 md->cFREEZE, md->cACC, md->cTC,
1728 state->x, xprime, state->v, force, state->sd_X,
1729 inputrec->opts.tau_t,
1730 FALSE, step, inputrec->ld_seed,
1731 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1733 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1734 wallcycle_stop(wcycle, ewcUPDATE);
1738 /* Constrain the coordinates xprime */
1739 wallcycle_start(wcycle, ewcCONSTR);
1740 constrain(NULL, bLog, bEner, constr, idef,
1741 inputrec, NULL, cr, step, 1, 1.0, md,
1742 state->x, xprime, NULL,
1743 bMolPBC, state->box,
1744 state->lambda[efptBONDED], dvdlambda,
1745 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
1746 wallcycle_stop(wcycle, ewcCONSTR);
1751 /* We must always unshift after updating coordinates; if we did not shake
1752 x was shifted in do_force */
1754 if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1756 if (graph && (graph->nnodes > 0))
1758 unshift_x(graph, state->box, state->x, upd->xp);
1759 if (TRICLINIC(state->box))
1761 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1765 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1770 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntUpdate)) schedule(static)
1771 for (i = start; i < nrend; i++)
1773 copy_rvec(upd->xp[i], state->x[i]);
1777 dump_it_all(fplog, "After unshift",
1778 state->natoms, state->x, upd->xp, state->v, force);
1780 /* ############# END the update of velocities and positions ######### */
1783 void update_box(FILE *fplog,
1785 t_inputrec *inputrec, /* input record and box stuff */
1788 rvec force[], /* forces on home particles */
1794 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE;
1797 int start, homenr, nrend, i, n, m, g;
1801 homenr = md->homenr;
1802 nrend = start+homenr;
1805 (inputrec->etc == etcNOSEHOOVER) ||
1806 (inputrec->epc == epcPARRINELLORAHMAN) ||
1807 (inputrec->epc == epcMTTK);
1809 dt = inputrec->delta_t;
1813 /* now update boxes */
1814 switch (inputrec->epc)
1818 case (epcBERENDSEN):
1819 berendsen_pscale(inputrec, pcoupl_mu, state->box, state->box_rel,
1820 start, homenr, state->x, md->cFREEZE, nrnb);
1822 case (epcPARRINELLORAHMAN):
1823 /* The box velocities were updated in do_pr_pcoupl in the update
1824 * iteration, but we dont change the box vectors until we get here
1825 * since we need to be able to shift/unshift above.
1827 for (i = 0; i < DIM; i++)
1829 for (m = 0; m <= i; m++)
1831 state->box[i][m] += dt*state->boxv[i][m];
1834 preserve_box_shape(inputrec, state->box_rel, state->box);
1836 /* Scale the coordinates */
1837 for (n = start; (n < start+homenr); n++)
1839 tmvmul_ur0(pcoupl_mu, state->x[n], state->x[n]);
1843 switch (inputrec->epct)
1845 case (epctISOTROPIC):
1846 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1847 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1848 Side length scales as exp(veta*dt) */
1850 msmul(state->box, exp(state->veta*dt), state->box);
1852 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1853 o If we assume isotropic scaling, and box length scaling
1854 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1855 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1856 determinant of B is L^DIM det(M), and the determinant
1857 of dB/dt is (dL/dT)^DIM det (M). veta will be
1858 (det(dB/dT)/det(B))^(1/3). Then since M =
1859 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1861 msmul(state->box, state->veta, state->boxv);
1871 if ((!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))) && scale_tot)
1873 /* The transposes of the scaling matrices are stored,
1874 * therefore we need to reverse the order in the multiplication.
1876 mmul_ur0(*scale_tot, pcoupl_mu, *scale_tot);
1879 if (DEFORM(*inputrec))
1881 deform(upd, start, homenr, state->x, state->box, scale_tot, inputrec, step);
1884 dump_it_all(fplog, "After update",
1885 state->natoms, state->x, upd->xp, state->v, force);
1888 void update_coords(FILE *fplog,
1890 t_inputrec *inputrec, /* input record and box stuff */
1894 rvec *f, /* forces on home particles */
1897 tensor *vir_lr_constr,
1899 gmx_ekindata_t *ekind,
1904 t_commrec *cr, /* these shouldn't be here -- need to think about it */
1906 gmx_constr_t constr,
1909 gmx_bool bNH, bPR, bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1911 real *imass, *imassin;
1914 int start, homenr, nrend, i, j, d, n, m, g;
1915 int blen0, blen1, iatom, jatom, nshake, nsettle, nconstr, nexpand;
1918 rvec *vcom, *xcom, *vall, *xall, *xin, *vin, *forcein, *fall, *xpall, *xprimein, *xprime;
1921 bDoConstr = (NULL != constr);
1923 /* Running the velocity half does nothing except for velocity verlet */
1924 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1925 !EI_VV(inputrec->eI))
1927 gmx_incons("update_coords called for velocity without VV integrator");
1931 homenr = md->homenr;
1932 nrend = start+homenr;
1934 xprime = get_xprime(state, upd);
1936 dt = inputrec->delta_t;
1939 /* We need to update the NMR restraint history when time averaging is used */
1940 if (state->flags & (1<<estDISRE_RM3TAV))
1942 update_disres_history(fcd, &state->hist);
1944 if (state->flags & (1<<estORIRE_DTAV))
1946 update_orires_history(fcd, &state->hist);
1950 bNH = inputrec->etc == etcNOSEHOOVER;
1951 bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1953 if (bDoLR && inputrec->nstcalclr > 1 && !EI_VV(inputrec->eI)) /* get this working with VV? */
1955 /* Store the total force + nstcalclr-1 times the LR force
1956 * in forces_lr, so it can be used in a normal update algorithm
1957 * to produce twin time stepping.
1959 /* is this correct in the new construction? MRS */
1961 inputrec->nstcalclr, constr, inputrec, md, idef, cr,
1962 step, state, bMolPBC,
1963 start, nrend, f, f_lr, vir_lr_constr, nrnb);
1971 /* ############# START The update of velocities and positions ######### */
1973 dump_it_all(fplog, "Before update",
1974 state->natoms, state->x, xprime, state->v, force);
1976 if (inputrec->eI == eiSD2)
1978 check_sd2_work_data_allocation(upd->sd, nrend);
1980 do_update_sd2_Tconsts(upd->sd,
1981 inputrec->opts.ngtc,
1982 inputrec->opts.tau_t,
1983 inputrec->opts.ref_t);
1985 if (inputrec->eI == eiBD)
1987 do_update_bd_Tconsts(dt, inputrec->bd_fric,
1988 inputrec->opts.ngtc, inputrec->opts.ref_t,
1992 nth = gmx_omp_nthreads_get(emntUpdate);
1994 #pragma omp parallel for num_threads(nth) schedule(static) private(alpha)
1995 for (th = 0; th < nth; th++)
1997 int start_th, end_th;
1999 start_th = start + ((nrend-start)* th )/nth;
2000 end_th = start + ((nrend-start)*(th+1))/nth;
2002 switch (inputrec->eI)
2005 if (ekind->cosacc.cos_accel == 0)
2007 do_update_md(start_th, end_th, dt,
2008 ekind->tcstat, state->nosehoover_vxi,
2009 ekind->bNEMD, ekind->grpstat, inputrec->opts.acc,
2010 inputrec->opts.nFreeze,
2011 md->invmass, md->ptype,
2012 md->cFREEZE, md->cACC, md->cTC,
2013 state->x, xprime, state->v, force, M,
2018 do_update_visc(start_th, end_th, dt,
2019 ekind->tcstat, state->nosehoover_vxi,
2020 md->invmass, md->ptype,
2021 md->cTC, state->x, xprime, state->v, force, M,
2023 ekind->cosacc.cos_accel,
2029 /* With constraints, the SD1 update is done in 2 parts */
2030 do_update_sd1(upd->sd,
2031 start_th, end_th, dt,
2032 inputrec->opts.acc, inputrec->opts.nFreeze,
2033 md->invmass, md->ptype,
2034 md->cFREEZE, md->cACC, md->cTC,
2035 state->x, xprime, state->v, force,
2036 inputrec->opts.ngtc, inputrec->opts.ref_t,
2038 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2041 /* The SD2 update is always done in 2 parts,
2042 * because an extra constraint step is needed
2044 do_update_sd2(upd->sd,
2045 bInitStep, start_th, end_th,
2046 inputrec->opts.acc, inputrec->opts.nFreeze,
2047 md->invmass, md->ptype,
2048 md->cFREEZE, md->cACC, md->cTC,
2049 state->x, xprime, state->v, force, state->sd_X,
2050 inputrec->opts.tau_t,
2051 TRUE, step, inputrec->ld_seed,
2052 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2055 do_update_bd(start_th, end_th, dt,
2056 inputrec->opts.nFreeze, md->invmass, md->ptype,
2057 md->cFREEZE, md->cTC,
2058 state->x, xprime, state->v, force,
2061 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2065 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
2070 do_update_vv_vel(start_th, end_th, dt,
2071 inputrec->opts.acc, inputrec->opts.nFreeze,
2072 md->invmass, md->ptype,
2073 md->cFREEZE, md->cACC,
2075 (bNH || bPR), state->veta, alpha);
2078 do_update_vv_pos(start_th, end_th, dt,
2079 inputrec->opts.nFreeze,
2080 md->ptype, md->cFREEZE,
2081 state->x, xprime, state->v,
2082 (bNH || bPR), state->veta);
2087 gmx_fatal(FARGS, "Don't know how to update coordinates");
2095 void correct_ekin(FILE *log, int start, int end, rvec v[], rvec vcm, real mass[],
2096 real tmass, tensor ekin)
2099 * This is a debugging routine. It should not be called for production code
2101 * The kinetic energy should calculated according to:
2102 * Ekin = 1/2 m (v-vcm)^2
2103 * However the correction is not always applied, since vcm may not be
2104 * known in time and we compute
2105 * Ekin' = 1/2 m v^2 instead
2106 * This can be corrected afterwards by computing
2107 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
2109 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
2116 /* Local particles */
2119 /* Processor dependent part. */
2121 for (i = start; (i < end); i++)
2125 for (j = 0; (j < DIM); j++)
2131 svmul(1/tmass, vcm, vcm);
2132 svmul(0.5, vcm, hvcm);
2134 for (j = 0; (j < DIM); j++)
2136 for (k = 0; (k < DIM); k++)
2138 dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
2141 pr_rvecs(log, 0, "dekin", dekin, DIM);
2142 pr_rvecs(log, 0, " ekin", ekin, DIM);
2143 fprintf(log, "dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
2144 trace(dekin), trace(ekin), vcm[XX], vcm[YY], vcm[ZZ]);
2145 fprintf(log, "mv = (%8.4f %8.4f %8.4f)\n",
2146 mv[XX], mv[YY], mv[ZZ]);
2149 extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_int64_t step, const t_commrec *cr,
2150 t_mdatoms *md, t_state *state, gmx_update_t upd, gmx_constr_t constr)
2154 real rate = (ir->delta_t)/ir->opts.tau_t[0];
2156 if (ir->etc == etcANDERSEN && constr != NULL)
2158 gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
2161 /* proceed with andersen if 1) it's fixed probability per
2162 particle andersen or 2) it's massive andersen and it's tau_t/dt */
2163 if ((ir->etc == etcANDERSEN) || do_per_step(step, (int)(1.0/rate)))
2165 andersen_tcoupl(ir, step, cr, md, state, rate,
2166 upd->sd->randomize_group, upd->sd->boltzfac);