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,2015,2016, 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.
39 #include "gromacs/legacyheaders/update.h"
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/legacyheaders/constr.h"
48 #include "gromacs/legacyheaders/disre.h"
49 #include "gromacs/legacyheaders/force.h"
50 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
51 #include "gromacs/legacyheaders/macros.h"
52 #include "gromacs/legacyheaders/mdrun.h"
53 #include "gromacs/legacyheaders/names.h"
54 #include "gromacs/legacyheaders/nrnb.h"
55 #include "gromacs/legacyheaders/orires.h"
56 #include "gromacs/legacyheaders/tgroup.h"
57 #include "gromacs/legacyheaders/txtdump.h"
58 #include "gromacs/legacyheaders/typedefs.h"
59 #include "gromacs/legacyheaders/types/commrec.h"
60 #include "gromacs/math/units.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/pbcutil/mshift.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/pulling/pull.h"
65 #include "gromacs/random/random.h"
66 #include "gromacs/timing/wallcycle.h"
67 #include "gromacs/utility/futil.h"
68 #include "gromacs/utility/gmxomp.h"
69 #include "gromacs/utility/smalloc.h"
71 /*For debugging, start at v(-dt/2) for velolcity verlet -- uncomment next line */
72 /*#define STARTFROMDT2*/
96 gmx_sd_sigma_t *sdsig;
99 /* andersen temperature control stuff */
100 gmx_bool *randomize_group;
104 typedef struct gmx_update
107 /* xprime for constraint algorithms */
111 /* Variables for the deform algorithm */
112 gmx_int64_t deformref_step;
113 matrix deformref_box;
117 static void do_update_md(int start, int nrend,
118 double dt, int nstpcouple,
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 - nstpcouple*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)
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)
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,
356 double dt, int nstpcouple,
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++)
401 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
403 vn = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
404 - nstpcouple*iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
407 vn += vc + dt*cosz*cos_accel;
410 xprime[n][d] = x[n][d]+vn*dt;
414 xprime[n][d] = x[n][d];
421 /* Classic version of update, used with berendsen coupling */
422 for (n = start; n < nrend; n++)
424 w_dt = invmass[n]*dt;
429 lg = tcstat[gt].lambda;
430 cosz = cos(fac*x[n][ZZ]);
432 for (d = 0; d < DIM; d++)
436 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
441 /* Do not scale the cosine velocity profile */
442 vv = vc + lg*(vn - vc + f[n][d]*w_dt);
443 /* Add the cosine accelaration profile */
444 vv += dt*cosz*cos_accel;
448 vv = lg*(vn + f[n][d]*w_dt);
451 xprime[n][d] = x[n][d]+vv*dt;
456 xprime[n][d] = x[n][d];
463 static gmx_stochd_t *init_stochd(t_inputrec *ir)
472 ngtc = ir->opts.ngtc;
476 snew(sd->bd_rf, ngtc);
478 else if (EI_SD(ir->eI))
481 snew(sd->sdsig, ngtc);
484 for (n = 0; n < ngtc; n++)
486 if (ir->opts.tau_t[n] > 0)
488 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
489 sdc[n].eph = exp(sdc[n].gdt/2);
490 sdc[n].emh = exp(-sdc[n].gdt/2);
491 sdc[n].em = exp(-sdc[n].gdt);
495 /* No friction and noise on this group */
501 if (sdc[n].gdt >= 0.05)
503 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
504 - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
505 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
506 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
511 /* Seventh order expansions for small y */
512 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
513 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))));
514 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
518 fprintf(debug, "SD const tc-grp %d: b %g c %g d %g\n",
519 n, sdc[n].b, sdc[n].c, sdc[n].d);
523 else if (ETC_ANDERSEN(ir->etc))
532 snew(sd->randomize_group, ngtc);
533 snew(sd->boltzfac, ngtc);
535 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
536 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
538 for (n = 0; n < ngtc; n++)
540 reft = std::max<real>(0, opts->ref_t[n]);
541 if ((opts->tau_t[n] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
543 sd->randomize_group[n] = TRUE;
544 sd->boltzfac[n] = BOLTZ*opts->ref_t[n];
548 sd->randomize_group[n] = FALSE;
555 gmx_update_t init_update(t_inputrec *ir)
561 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
563 upd->sd = init_stochd(ir);
572 static void do_update_sd1(gmx_stochd_t *sd,
573 int start, int nrend, double dt,
574 rvec accel[], ivec nFreeze[],
575 real invmass[], unsigned short ptype[],
576 unsigned short cFREEZE[], unsigned short cACC[],
577 unsigned short cTC[],
578 rvec x[], rvec xprime[], rvec v[], rvec f[],
579 int ngtc, real ref_t[],
581 gmx_bool bFirstHalfConstr,
582 gmx_int64_t step, int seed, int* gatindex)
587 int gf = 0, ga = 0, gt = 0;
594 for (n = 0; n < ngtc; n++)
597 /* The mass is encounted for later, since this differs per atom */
598 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;
608 ism = sqrt(invmass[n]);
622 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
624 for (d = 0; d < DIM; d++)
626 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
630 sd_V = ism*sig[gt].V*rnd[d];
631 vn = v[n][d] + (invmass[n]*f[n][d] + accel[ga][d])*dt;
632 v[n][d] = vn*sdc[gt].em + sd_V;
633 /* Here we include half of the friction+noise
634 * update of v into the integration of x.
636 xprime[n][d] = x[n][d] + 0.5*(vn + v[n][d])*dt;
641 xprime[n][d] = x[n][d];
648 /* We do have constraints */
649 if (bFirstHalfConstr)
651 /* First update without friction and noise */
654 for (n = start; n < nrend; n++)
667 for (d = 0; d < DIM; d++)
669 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
671 v[n][d] = v[n][d] + (im*f[n][d] + accel[ga][d])*dt;
672 xprime[n][d] = x[n][d] + v[n][d]*dt;
677 xprime[n][d] = x[n][d];
684 /* Update friction and noise only */
685 for (n = start; n < nrend; n++)
688 int ng = gatindex ? gatindex[n] : n;
690 ism = sqrt(invmass[n]);
700 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
702 for (d = 0; d < DIM; d++)
704 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
708 sd_V = ism*sig[gt].V*rnd[d];
710 v[n][d] = vn*sdc[gt].em + sd_V;
711 /* Add the friction and noise contribution only */
712 xprime[n][d] = xprime[n][d] + 0.5*(v[n][d] - vn)*dt;
720 static void check_sd2_work_data_allocation(gmx_stochd_t *sd, int nrend)
722 if (nrend > sd->sd_V_nalloc)
724 sd->sd_V_nalloc = over_alloc_dd(nrend);
725 srenew(sd->sd_V, sd->sd_V_nalloc);
729 static void do_update_sd2_Tconsts(gmx_stochd_t *sd,
734 /* This is separated from the update below, because it is single threaded */
743 for (gt = 0; gt < ngtc; gt++)
745 kT = BOLTZ*ref_t[gt];
746 /* The mass is encounted for later, since this differs per atom */
747 sig[gt].V = sqrt(kT*(1-sdc[gt].em));
748 sig[gt].X = sqrt(kT*sqr(tau_t[gt])*sdc[gt].c);
749 sig[gt].Yv = sqrt(kT*sdc[gt].b/sdc[gt].c);
750 sig[gt].Yx = sqrt(kT*sqr(tau_t[gt])*sdc[gt].b/(1-sdc[gt].em));
754 static void do_update_sd2(gmx_stochd_t *sd,
756 int start, int nrend,
757 rvec accel[], ivec nFreeze[],
758 real invmass[], unsigned short ptype[],
759 unsigned short cFREEZE[], unsigned short cACC[],
760 unsigned short cTC[],
761 rvec x[], rvec xprime[], rvec v[], rvec f[],
764 gmx_bool bFirstHalf, gmx_int64_t step, int seed,
769 /* The random part of the velocity update, generated in the first
770 * half of the update, needs to be remembered for the second half.
773 int gf = 0, ga = 0, gt = 0;
774 real vn = 0, Vmh, Xmh;
782 for (n = start; n < nrend; n++)
784 real rnd[6], rndi[3];
785 ng = gatindex ? gatindex[n] : n;
786 ism = sqrt(invmass[n]);
800 gmx_rng_cycle_6gaussian_table(step*2+(bFirstHalf ? 1 : 2), ng, seed, RND_SEED_UPDATE, rnd);
803 gmx_rng_cycle_3gaussian_table(step*2, ng, seed, RND_SEED_UPDATE, rndi);
805 for (d = 0; d < DIM; d++)
811 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
817 sd_X[n][d] = ism*sig[gt].X*rndi[d];
819 Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
820 + ism*sig[gt].Yv*rnd[d*2];
821 sd_V[n][d] = ism*sig[gt].V*rnd[d*2+1];
823 v[n][d] = vn*sdc[gt].em
824 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
825 + sd_V[n][d] - sdc[gt].em*Vmh;
827 xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
831 /* Correct the velocities for the constraints.
832 * This operation introduces some inaccuracy,
833 * since the velocity is determined from differences in coordinates.
836 (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
838 Xmh = sd_V[n][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
839 + ism*sig[gt].Yx*rnd[d*2];
840 sd_X[n][d] = ism*sig[gt].X*rnd[d*2+1];
842 xprime[n][d] += sd_X[n][d] - Xmh;
851 xprime[n][d] = x[n][d];
858 static void do_update_bd_Tconsts(double dt, real friction_coefficient,
859 int ngtc, const real ref_t[],
862 /* This is separated from the update below, because it is single threaded */
865 if (friction_coefficient != 0)
867 for (gt = 0; gt < ngtc; gt++)
869 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]/(friction_coefficient*dt));
874 for (gt = 0; gt < ngtc; gt++)
876 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]);
881 static void do_update_bd(int start, int nrend, double dt,
883 real invmass[], unsigned short ptype[],
884 unsigned short cFREEZE[], unsigned short cTC[],
885 rvec x[], rvec xprime[], rvec v[],
886 rvec f[], real friction_coefficient,
887 real *rf, gmx_int64_t step, int seed,
890 /* note -- these appear to be full step velocities . . . */
896 if (friction_coefficient != 0)
898 invfr = 1.0/friction_coefficient;
901 for (n = start; (n < nrend); n++)
904 int ng = gatindex ? gatindex[n] : n;
914 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
915 for (d = 0; (d < DIM); d++)
917 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
919 if (friction_coefficient != 0)
921 vn = invfr*f[n][d] + rf[gt]*rnd[d];
925 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
926 vn = 0.5*invmass[n]*f[n][d]*dt
927 + sqrt(0.5*invmass[n])*rf[gt]*rnd[d];
931 xprime[n][d] = x[n][d]+vn*dt;
936 xprime[n][d] = x[n][d];
942 static void dump_it_all(FILE gmx_unused *fp, const char gmx_unused *title,
943 int gmx_unused natoms, rvec gmx_unused x[], rvec gmx_unused xp[],
944 rvec gmx_unused v[], rvec gmx_unused f[])
949 fprintf(fp, "%s\n", title);
950 pr_rvecs(fp, 0, "x", x, natoms);
951 pr_rvecs(fp, 0, "xp", xp, natoms);
952 pr_rvecs(fp, 0, "v", v, natoms);
953 pr_rvecs(fp, 0, "f", f, natoms);
958 static void calc_ke_part_normal(rvec v[], t_grpopts *opts, t_mdatoms *md,
959 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel)
962 t_grp_tcstat *tcstat = ekind->tcstat;
963 t_grp_acc *grpstat = ekind->grpstat;
966 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
967 an option, but not supported now.
968 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
971 /* group velocities are calculated in update_ekindata and
972 * accumulated in acumulate_groups.
973 * Now the partial global and groups ekin.
975 for (g = 0; (g < opts->ngtc); g++)
977 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
980 clear_mat(tcstat[g].ekinf);
981 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
985 clear_mat(tcstat[g].ekinh);
988 ekind->dekindl_old = ekind->dekindl;
989 nthread = gmx_omp_nthreads_get(emntUpdate);
991 #pragma omp parallel for num_threads(nthread) schedule(static)
992 for (thread = 0; thread < nthread; thread++)
994 int start_t, end_t, n;
1002 start_t = ((thread+0)*md->homenr)/nthread;
1003 end_t = ((thread+1)*md->homenr)/nthread;
1005 ekin_sum = ekind->ekin_work[thread];
1006 dekindl_sum = ekind->dekindl_work[thread];
1008 for (gt = 0; gt < opts->ngtc; gt++)
1010 clear_mat(ekin_sum[gt]);
1016 for (n = start_t; n < end_t; n++)
1026 hm = 0.5*md->massT[n];
1028 for (d = 0; (d < DIM); d++)
1030 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
1032 for (d = 0; (d < DIM); d++)
1034 for (m = 0; (m < DIM); m++)
1036 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1037 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
1040 if (md->nMassPerturbed && md->bPerturbed[n])
1043 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1049 for (thread = 0; thread < nthread; thread++)
1051 for (g = 0; g < opts->ngtc; g++)
1055 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
1060 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1065 ekind->dekindl += *ekind->dekindl_work[thread];
1068 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1071 static void calc_ke_part_visc(matrix box, rvec x[], rvec v[],
1072 t_grpopts *opts, t_mdatoms *md,
1073 gmx_ekindata_t *ekind,
1074 t_nrnb *nrnb, gmx_bool bEkinAveVel)
1076 int start = 0, homenr = md->homenr;
1077 int g, d, n, m, gt = 0;
1080 t_grp_tcstat *tcstat = ekind->tcstat;
1081 t_cos_acc *cosacc = &(ekind->cosacc);
1086 for (g = 0; g < opts->ngtc; g++)
1088 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1089 clear_mat(ekind->tcstat[g].ekinh);
1091 ekind->dekindl_old = ekind->dekindl;
1093 fac = 2*M_PI/box[ZZ][ZZ];
1096 for (n = start; n < start+homenr; n++)
1102 hm = 0.5*md->massT[n];
1104 /* Note that the times of x and v differ by half a step */
1105 /* MRS -- would have to be changed for VV */
1106 cosz = cos(fac*x[n][ZZ]);
1107 /* Calculate the amplitude of the new velocity profile */
1108 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1110 copy_rvec(v[n], v_corrt);
1111 /* Subtract the profile for the kinetic energy */
1112 v_corrt[XX] -= cosz*cosacc->vcos;
1113 for (d = 0; (d < DIM); d++)
1115 for (m = 0; (m < DIM); m++)
1117 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1120 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1124 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1128 if (md->nPerturbed && md->bPerturbed[n])
1130 /* The minus sign here might be confusing.
1131 * The kinetic contribution from dH/dl doesn't come from
1132 * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
1133 * where p are the momenta. The difference is only a minus sign.
1135 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1138 ekind->dekindl = dekindl;
1139 cosacc->mvcos = mvcos;
1141 inc_nrnb(nrnb, eNR_EKIN, homenr);
1144 void calc_ke_part(t_state *state, t_grpopts *opts, t_mdatoms *md,
1145 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel)
1147 if (ekind->cosacc.cos_accel == 0)
1149 calc_ke_part_normal(state->v, opts, md, ekind, nrnb, bEkinAveVel);
1153 calc_ke_part_visc(state->box, state->x, state->v, opts, md, ekind, nrnb, bEkinAveVel);
1157 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1159 ekinstate->ekin_n = ir->opts.ngtc;
1160 snew(ekinstate->ekinh, ekinstate->ekin_n);
1161 snew(ekinstate->ekinf, ekinstate->ekin_n);
1162 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1163 snew(ekinstate->ekinscalef_nhc, ekinstate->ekin_n);
1164 snew(ekinstate->ekinscaleh_nhc, ekinstate->ekin_n);
1165 snew(ekinstate->vscale_nhc, ekinstate->ekin_n);
1166 ekinstate->dekindl = 0;
1167 ekinstate->mvcos = 0;
1170 void update_ekinstate(ekinstate_t *ekinstate, gmx_ekindata_t *ekind)
1174 for (i = 0; i < ekinstate->ekin_n; i++)
1176 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1177 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1178 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1179 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1180 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1181 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1184 copy_mat(ekind->ekin, ekinstate->ekin_total);
1185 ekinstate->dekindl = ekind->dekindl;
1186 ekinstate->mvcos = ekind->cosacc.mvcos;
1190 void restore_ekinstate_from_state(t_commrec *cr,
1191 gmx_ekindata_t *ekind, ekinstate_t *ekinstate)
1197 for (i = 0; i < ekinstate->ekin_n; i++)
1199 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1200 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1201 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1202 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1203 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1204 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1207 copy_mat(ekinstate->ekin_total, ekind->ekin);
1209 ekind->dekindl = ekinstate->dekindl;
1210 ekind->cosacc.mvcos = ekinstate->mvcos;
1211 n = ekinstate->ekin_n;
1216 gmx_bcast(sizeof(n), &n, cr);
1217 for (i = 0; i < n; i++)
1219 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1220 ekind->tcstat[i].ekinh[0], cr);
1221 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1222 ekind->tcstat[i].ekinf[0], cr);
1223 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1224 ekind->tcstat[i].ekinh_old[0], cr);
1226 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1227 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1228 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1229 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1230 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1231 &(ekind->tcstat[i].vscale_nhc), cr);
1233 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1234 ekind->ekin[0], cr);
1236 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1237 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1241 void set_deform_reference_box(gmx_update_t upd, gmx_int64_t step, matrix box)
1243 upd->deformref_step = step;
1244 copy_mat(box, upd->deformref_box);
1247 static void deform(gmx_update_t upd,
1248 int start, int homenr, rvec x[], matrix box,
1249 const t_inputrec *ir, gmx_int64_t step)
1251 matrix bnew, invbox, mu;
1255 elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1256 copy_mat(box, bnew);
1257 for (i = 0; i < DIM; i++)
1259 for (j = 0; j < DIM; j++)
1261 if (ir->deform[i][j] != 0)
1264 upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1268 /* We correct the off-diagonal elements,
1269 * which can grow indefinitely during shearing,
1270 * so the shifts do not get messed up.
1272 for (i = 1; i < DIM; i++)
1274 for (j = i-1; j >= 0; j--)
1276 while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1278 rvec_dec(bnew[i], bnew[j]);
1280 while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1282 rvec_inc(bnew[i], bnew[j]);
1286 m_inv_ur0(box, invbox);
1287 copy_mat(bnew, box);
1288 mmul_ur0(box, invbox, mu);
1290 for (i = start; i < start+homenr; i++)
1292 x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1293 x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1294 x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1298 void update_tcouple(gmx_int64_t step,
1299 t_inputrec *inputrec,
1301 gmx_ekindata_t *ekind,
1306 gmx_bool bTCouple = FALSE;
1310 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1311 if (inputrec->etc != etcNO &&
1312 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec)))
1314 /* We should only couple after a step where energies were determined (for leapfrog versions)
1315 or the step energies are determined, for velocity verlet versions */
1317 if (EI_VV(inputrec->eI))
1325 bTCouple = (inputrec->nsttcouple == 1 ||
1326 do_per_step(step+inputrec->nsttcouple-offset,
1327 inputrec->nsttcouple));
1332 dttc = inputrec->nsttcouple*inputrec->delta_t;
1334 switch (inputrec->etc)
1339 berendsen_tcoupl(inputrec, ekind, dttc);
1342 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1343 state->nosehoover_xi, state->nosehoover_vxi, MassQ);
1346 vrescale_tcoupl(inputrec, step, ekind, dttc,
1347 state->therm_integral);
1350 /* rescale in place here */
1351 if (EI_VV(inputrec->eI))
1353 rescale_velocities(ekind, md, 0, md->homenr, state->v);
1358 /* Set the T scaling lambda to 1 to have no scaling */
1359 for (i = 0; (i < inputrec->opts.ngtc); i++)
1361 ekind->tcstat[i].lambda = 1.0;
1366 void update_pcouple(FILE *fplog,
1368 t_inputrec *inputrec,
1374 gmx_bool bPCouple = FALSE;
1378 /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1379 if (inputrec->epc != epcNO && (!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))))
1381 /* We should only couple after a step where energies were determined */
1382 bPCouple = (inputrec->nstpcouple == 1 ||
1383 do_per_step(step+inputrec->nstpcouple-1,
1384 inputrec->nstpcouple));
1387 clear_mat(pcoupl_mu);
1388 for (i = 0; i < DIM; i++)
1390 pcoupl_mu[i][i] = 1.0;
1397 dtpc = inputrec->nstpcouple*inputrec->delta_t;
1399 switch (inputrec->epc)
1401 /* We can always pcoupl, even if we did not sum the energies
1402 * the previous step, since state->pres_prev is only updated
1403 * when the energies have been summed.
1407 case (epcBERENDSEN):
1410 berendsen_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev, state->box,
1414 case (epcPARRINELLORAHMAN):
1415 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1416 state->box, state->box_rel, state->boxv,
1417 M, pcoupl_mu, bInitStep);
1425 static rvec *get_xprime(const t_state *state, gmx_update_t upd)
1427 if (state->nalloc > upd->xp_nalloc)
1429 upd->xp_nalloc = state->nalloc;
1430 srenew(upd->xp, upd->xp_nalloc);
1436 static void combine_forces(gmx_update_t upd,
1438 gmx_constr_t constr,
1439 t_inputrec *ir, t_mdatoms *md, t_idef *idef,
1442 t_state *state, gmx_bool bMolPBC,
1443 int start, int nrend,
1444 rvec f[], rvec f_lr[],
1445 tensor *vir_lr_constr,
1450 /* f contains the short-range forces + the long range forces
1451 * which are stored separately in f_lr.
1454 if (constr != NULL && vir_lr_constr != NULL &&
1455 !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1457 /* We need to constrain the LR forces separately,
1458 * because due to the different pre-factor for the SR and LR
1459 * forces in the update algorithm, we have to correct
1460 * the constraint virial for the nstcalclr-1 extra f_lr.
1461 * Constrain only the additional LR part of the force.
1463 /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1468 xp = get_xprime(state, upd);
1470 fac = (nstcalclr - 1)*ir->delta_t*ir->delta_t;
1472 for (i = 0; i < md->homenr; i++)
1474 if (md->cFREEZE != NULL)
1476 gf = md->cFREEZE[i];
1478 for (d = 0; d < DIM; d++)
1480 if ((md->ptype[i] != eptVSite) &&
1481 (md->ptype[i] != eptShell) &&
1482 !ir->opts.nFreeze[gf][d])
1484 xp[i][d] = state->x[i][d] + fac*f_lr[i][d]*md->invmass[i];
1488 xp[i][d] = state->x[i][d];
1492 constrain(NULL, FALSE, FALSE, constr, idef, ir, cr, step, 0, 1.0, md,
1493 state->x, xp, xp, bMolPBC, state->box, state->lambda[efptBONDED], NULL,
1494 NULL, vir_lr_constr, nrnb, econqForce);
1497 /* Add nstcalclr-1 times the LR force to the sum of both forces
1498 * and store the result in forces_lr.
1500 for (i = start; i < nrend; i++)
1502 for (d = 0; d < DIM; d++)
1504 f_lr[i][d] = f[i][d] + (nstcalclr - 1)*f_lr[i][d];
1509 void update_constraints(FILE *fplog,
1511 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1512 t_inputrec *inputrec, /* input record and box stuff */
1517 rvec force[], /* forces on home particles */
1522 gmx_wallcycle_t wcycle,
1524 gmx_constr_t constr,
1525 gmx_bool bFirstHalf,
1528 gmx_bool bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1530 int start, homenr, nrend, i, m;
1532 rvec *xprime = NULL;
1539 if (bFirstHalf && !EI_VV(inputrec->eI))
1544 /* for now, SD update is here -- though it really seems like it
1545 should be reformulated as a velocity verlet method, since it has two parts */
1548 homenr = md->homenr;
1549 nrend = start+homenr;
1551 dt = inputrec->delta_t;
1555 * APPLY CONSTRAINTS:
1558 * When doing PR pressure coupling we have to constrain the
1559 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1560 * it is enough to do this once though, since the relative velocities
1561 * after this will be normal to the bond vector
1566 /* clear out constraints before applying */
1567 clear_mat(vir_part);
1569 xprime = get_xprime(state, upd);
1571 bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1572 bLog = (do_per_step(step, inputrec->nstlog) || bLastStep || (step < 0));
1573 bEner = (do_per_step(step, inputrec->nstenergy) || bLastStep);
1574 /* Constrain the coordinates xprime */
1575 wallcycle_start(wcycle, ewcCONSTR);
1576 if (EI_VV(inputrec->eI) && bFirstHalf)
1578 constrain(NULL, bLog, bEner, constr, idef,
1579 inputrec, cr, step, 1, 1.0, md,
1580 state->x, state->v, state->v,
1581 bMolPBC, state->box,
1582 state->lambda[efptBONDED], dvdlambda,
1583 NULL, bCalcVir ? &vir_con : NULL, nrnb, econqVeloc);
1587 constrain(NULL, bLog, bEner, constr, idef,
1588 inputrec, cr, step, 1, 1.0, md,
1589 state->x, xprime, NULL,
1590 bMolPBC, state->box,
1591 state->lambda[efptBONDED], dvdlambda,
1592 state->v, bCalcVir ? &vir_con : NULL, nrnb, econqCoord);
1594 wallcycle_stop(wcycle, ewcCONSTR);
1598 dump_it_all(fplog, "After Shake",
1599 state->natoms, state->x, xprime, state->v, force);
1603 if (inputrec->eI == eiSD2)
1605 /* A correction factor eph is needed for the SD constraint force */
1606 /* Here we can, unfortunately, not have proper corrections
1607 * for different friction constants, so we use the first one.
1609 for (i = 0; i < DIM; i++)
1611 for (m = 0; m < DIM; m++)
1613 vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1619 m_add(vir_part, vir_con, vir_part);
1623 pr_rvecs(debug, 0, "constraint virial", vir_part, DIM);
1630 if (inputrec->eI == eiSD1 && bDoConstr && !bFirstHalf)
1632 wallcycle_start(wcycle, ewcUPDATE);
1633 xprime = get_xprime(state, upd);
1635 nth = gmx_omp_nthreads_get(emntUpdate);
1637 #pragma omp parallel for num_threads(nth) schedule(static)
1639 for (th = 0; th < nth; th++)
1641 int start_th, end_th;
1643 start_th = start + ((nrend-start)* th )/nth;
1644 end_th = start + ((nrend-start)*(th+1))/nth;
1646 /* The second part of the SD integration */
1647 do_update_sd1(upd->sd,
1648 start_th, end_th, dt,
1649 inputrec->opts.acc, inputrec->opts.nFreeze,
1650 md->invmass, md->ptype,
1651 md->cFREEZE, md->cACC, md->cTC,
1652 state->x, xprime, state->v, force,
1653 inputrec->opts.ngtc, inputrec->opts.ref_t,
1655 step, inputrec->ld_seed,
1656 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1658 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1659 wallcycle_stop(wcycle, ewcUPDATE);
1663 /* Constrain the coordinates xprime for half a time step */
1664 wallcycle_start(wcycle, ewcCONSTR);
1666 constrain(NULL, bLog, bEner, constr, idef,
1667 inputrec, cr, step, 1, 0.5, md,
1668 state->x, xprime, NULL,
1669 bMolPBC, state->box,
1670 state->lambda[efptBONDED], dvdlambda,
1671 state->v, NULL, nrnb, econqCoord);
1673 wallcycle_stop(wcycle, ewcCONSTR);
1677 if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1679 wallcycle_start(wcycle, ewcUPDATE);
1680 xprime = get_xprime(state, upd);
1682 nth = gmx_omp_nthreads_get(emntUpdate);
1684 #pragma omp parallel for num_threads(nth) schedule(static)
1685 for (th = 0; th < nth; th++)
1687 int start_th, end_th;
1689 start_th = start + ((nrend-start)* th )/nth;
1690 end_th = start + ((nrend-start)*(th+1))/nth;
1692 /* The second part of the SD integration */
1693 do_update_sd2(upd->sd,
1694 FALSE, start_th, end_th,
1695 inputrec->opts.acc, inputrec->opts.nFreeze,
1696 md->invmass, md->ptype,
1697 md->cFREEZE, md->cACC, md->cTC,
1698 state->x, xprime, state->v, force, state->sd_X,
1699 inputrec->opts.tau_t,
1700 FALSE, step, inputrec->ld_seed,
1701 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1703 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1704 wallcycle_stop(wcycle, ewcUPDATE);
1708 /* Constrain the coordinates xprime */
1709 wallcycle_start(wcycle, ewcCONSTR);
1710 constrain(NULL, bLog, bEner, constr, idef,
1711 inputrec, cr, step, 1, 1.0, md,
1712 state->x, xprime, NULL,
1713 bMolPBC, state->box,
1714 state->lambda[efptBONDED], dvdlambda,
1715 NULL, NULL, nrnb, econqCoord);
1716 wallcycle_stop(wcycle, ewcCONSTR);
1721 /* We must always unshift after updating coordinates; if we did not shake
1722 x was shifted in do_force */
1724 if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1726 /* NOTE This part of the update actually does not belong with
1727 * the constraints, since we also call it without constraints.
1728 * But currently we always integrate to a temporary buffer and
1729 * then copy the results back here.
1731 wallcycle_start_nocount(wcycle, ewcUPDATE);
1733 if (md->cFREEZE != NULL && constr != NULL)
1735 /* If we have atoms that are frozen along some, but not all
1736 * dimensions, the constraints will have moved them also along
1737 * the frozen dimensions. To freeze such degrees of freedom
1738 * we copy them back here to later copy them forward. It would
1739 * be more elegant and slightly more efficient to copies zero
1740 * times instead of twice, but the graph case below prevents this.
1742 const ivec *nFreeze = inputrec->opts.nFreeze;
1743 bool partialFreezeAndConstraints = false;
1744 for (int g = 0; g < inputrec->opts.ngfrz; g++)
1746 int numFreezeDim = nFreeze[g][XX] + nFreeze[g][YY] + nFreeze[g][ZZ];
1747 if (numFreezeDim > 0 && numFreezeDim < 3)
1749 partialFreezeAndConstraints = true;
1752 if (partialFreezeAndConstraints)
1754 for (int i = start; i < nrend; i++)
1756 int g = md->cFREEZE[i];
1758 for (int d = 0; d < DIM; d++)
1762 upd->xp[i][d] = state->x[i][d];
1769 if (graph && (graph->nnodes > 0))
1771 unshift_x(graph, state->box, state->x, upd->xp);
1772 if (TRICLINIC(state->box))
1774 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1778 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1783 #ifndef __clang_analyzer__
1784 // cppcheck-suppress unreadVariable
1785 nth = gmx_omp_nthreads_get(emntUpdate);
1787 #pragma omp parallel for num_threads(nth) schedule(static)
1788 for (i = start; i < nrend; i++)
1790 copy_rvec(upd->xp[i], state->x[i]);
1793 wallcycle_stop(wcycle, ewcUPDATE);
1795 dump_it_all(fplog, "After unshift",
1796 state->natoms, state->x, upd->xp, state->v, force);
1798 /* ############# END the update of velocities and positions ######### */
1801 void update_box(FILE *fplog,
1803 t_inputrec *inputrec, /* input record and box stuff */
1806 rvec force[], /* forces on home particles */
1812 int start, homenr, i, n, m;
1815 homenr = md->homenr;
1817 dt = inputrec->delta_t;
1821 /* now update boxes */
1822 switch (inputrec->epc)
1826 case (epcBERENDSEN):
1827 /* We should only scale after a step where the pressure (kinetic
1828 * energy and virial) was calculated. This happens after the
1829 * coordinate update, whereas the current routine is called before
1830 * that, so we scale when step % nstpcouple = 1 instead of 0.
1832 if (inputrec->nstpcouple == 1 || (step % inputrec->nstpcouple == 1))
1834 berendsen_pscale(inputrec, pcoupl_mu, state->box, state->box_rel,
1835 start, homenr, state->x, md->cFREEZE, nrnb);
1838 case (epcPARRINELLORAHMAN):
1839 /* The box velocities were updated in do_pr_pcoupl in the update
1840 * iteration, but we dont change the box vectors until we get here
1841 * since we need to be able to shift/unshift above.
1843 for (i = 0; i < DIM; i++)
1845 for (m = 0; m <= i; m++)
1847 state->box[i][m] += dt*state->boxv[i][m];
1850 preserve_box_shape(inputrec, state->box_rel, state->box);
1852 /* Scale the coordinates */
1853 for (n = start; (n < start+homenr); n++)
1855 tmvmul_ur0(pcoupl_mu, state->x[n], state->x[n]);
1859 switch (inputrec->epct)
1861 case (epctISOTROPIC):
1862 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1863 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1864 Side length scales as exp(veta*dt) */
1866 msmul(state->box, exp(state->veta*dt), state->box);
1868 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1869 o If we assume isotropic scaling, and box length scaling
1870 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1871 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1872 determinant of B is L^DIM det(M), and the determinant
1873 of dB/dt is (dL/dT)^DIM det (M). veta will be
1874 (det(dB/dT)/det(B))^(1/3). Then since M =
1875 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1877 msmul(state->box, state->veta, state->boxv);
1887 if (DEFORM(*inputrec))
1889 deform(upd, start, homenr, state->x, state->box, inputrec, step);
1892 dump_it_all(fplog, "After update",
1893 state->natoms, state->x, upd->xp, state->v, force);
1896 void update_coords(FILE *fplog,
1898 t_inputrec *inputrec, /* input record and box stuff */
1902 rvec *f, /* forces on home particles */
1905 tensor *vir_lr_constr,
1907 gmx_ekindata_t *ekind,
1912 t_commrec *cr, /* these shouldn't be here -- need to think about it */
1914 gmx_constr_t constr,
1917 gmx_bool bNH, bPR, bDoConstr = FALSE;
1920 int start, homenr, nrend;
1924 bDoConstr = (NULL != constr);
1926 /* Running the velocity half does nothing except for velocity verlet */
1927 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1928 !EI_VV(inputrec->eI))
1930 gmx_incons("update_coords called for velocity without VV integrator");
1934 homenr = md->homenr;
1935 nrend = start+homenr;
1937 xprime = get_xprime(state, upd);
1939 dt = inputrec->delta_t;
1941 /* We need to update the NMR restraint history when time averaging is used */
1942 if (state->flags & (1<<estDISRE_RM3TAV))
1944 update_disres_history(fcd, &state->hist);
1946 if (state->flags & (1<<estORIRE_DTAV))
1948 update_orires_history(fcd, &state->hist);
1952 bNH = inputrec->etc == etcNOSEHOOVER;
1953 bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1955 if (bDoLR && inputrec->nstcalclr > 1 && !EI_VV(inputrec->eI)) /* get this working with VV? */
1957 /* Store the total force + nstcalclr-1 times the LR force
1958 * in forces_lr, so it can be used in a normal update algorithm
1959 * to produce twin time stepping.
1961 /* is this correct in the new construction? MRS */
1963 inputrec->nstcalclr, constr, inputrec, md, idef, cr,
1964 step, state, bMolPBC,
1965 start, nrend, f, f_lr, vir_lr_constr, nrnb);
1973 /* ############# START The update of velocities and positions ######### */
1975 dump_it_all(fplog, "Before update",
1976 state->natoms, state->x, xprime, state->v, force);
1978 if (inputrec->eI == eiSD2)
1980 check_sd2_work_data_allocation(upd->sd, nrend);
1982 do_update_sd2_Tconsts(upd->sd,
1983 inputrec->opts.ngtc,
1984 inputrec->opts.tau_t,
1985 inputrec->opts.ref_t);
1987 if (inputrec->eI == eiBD)
1989 do_update_bd_Tconsts(dt, inputrec->bd_fric,
1990 inputrec->opts.ngtc, inputrec->opts.ref_t,
1994 nth = gmx_omp_nthreads_get(emntUpdate);
1996 #pragma omp parallel for num_threads(nth) schedule(static) private(alpha)
1997 for (th = 0; th < nth; th++)
1999 int start_th, end_th;
2001 start_th = start + ((nrend-start)* th )/nth;
2002 end_th = start + ((nrend-start)*(th+1))/nth;
2004 switch (inputrec->eI)
2007 if (ekind->cosacc.cos_accel == 0)
2009 do_update_md(start_th, end_th,
2010 dt, inputrec->nstpcouple,
2011 ekind->tcstat, state->nosehoover_vxi,
2012 ekind->bNEMD, ekind->grpstat, inputrec->opts.acc,
2013 inputrec->opts.nFreeze,
2014 md->invmass, md->ptype,
2015 md->cFREEZE, md->cACC, md->cTC,
2016 state->x, xprime, state->v, force, M,
2021 do_update_visc(start_th, end_th,
2022 dt, inputrec->nstpcouple,
2023 ekind->tcstat, state->nosehoover_vxi,
2024 md->invmass, md->ptype,
2025 md->cTC, state->x, xprime, state->v, force, M,
2027 ekind->cosacc.cos_accel,
2033 /* With constraints, the SD1 update is done in 2 parts */
2034 do_update_sd1(upd->sd,
2035 start_th, end_th, dt,
2036 inputrec->opts.acc, inputrec->opts.nFreeze,
2037 md->invmass, md->ptype,
2038 md->cFREEZE, md->cACC, md->cTC,
2039 state->x, xprime, state->v, force,
2040 inputrec->opts.ngtc, inputrec->opts.ref_t,
2042 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2045 /* The SD2 update is always done in 2 parts,
2046 * because an extra constraint step is needed
2048 do_update_sd2(upd->sd,
2049 bInitStep, start_th, end_th,
2050 inputrec->opts.acc, inputrec->opts.nFreeze,
2051 md->invmass, md->ptype,
2052 md->cFREEZE, md->cACC, md->cTC,
2053 state->x, xprime, state->v, force, state->sd_X,
2054 inputrec->opts.tau_t,
2055 TRUE, step, inputrec->ld_seed,
2056 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2059 do_update_bd(start_th, end_th, dt,
2060 inputrec->opts.nFreeze, md->invmass, md->ptype,
2061 md->cFREEZE, md->cTC,
2062 state->x, xprime, state->v, force,
2065 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2069 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
2074 do_update_vv_vel(start_th, end_th, dt,
2075 inputrec->opts.acc, inputrec->opts.nFreeze,
2076 md->invmass, md->ptype,
2077 md->cFREEZE, md->cACC,
2079 (bNH || bPR), state->veta, alpha);
2082 do_update_vv_pos(start_th, end_th, dt,
2083 inputrec->opts.nFreeze,
2084 md->ptype, md->cFREEZE,
2085 state->x, xprime, state->v,
2086 (bNH || bPR), state->veta);
2091 gmx_fatal(FARGS, "Don't know how to update coordinates");
2099 void correct_ekin(FILE *log, int start, int end, rvec v[], rvec vcm, real mass[],
2100 real tmass, tensor ekin)
2103 * This is a debugging routine. It should not be called for production code
2105 * The kinetic energy should calculated according to:
2106 * Ekin = 1/2 m (v-vcm)^2
2107 * However the correction is not always applied, since vcm may not be
2108 * known in time and we compute
2109 * Ekin' = 1/2 m v^2 instead
2110 * This can be corrected afterwards by computing
2111 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
2113 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
2120 /* Local particles */
2123 /* Processor dependent part. */
2125 for (i = start; (i < end); i++)
2129 for (j = 0; (j < DIM); j++)
2135 svmul(1/tmass, vcm, vcm);
2136 svmul(0.5, vcm, hvcm);
2138 for (j = 0; (j < DIM); j++)
2140 for (k = 0; (k < DIM); k++)
2142 dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
2145 pr_rvecs(log, 0, "dekin", dekin, DIM);
2146 pr_rvecs(log, 0, " ekin", ekin, DIM);
2147 fprintf(log, "dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
2148 trace(dekin), trace(ekin), vcm[XX], vcm[YY], vcm[ZZ]);
2149 fprintf(log, "mv = (%8.4f %8.4f %8.4f)\n",
2150 mv[XX], mv[YY], mv[ZZ]);
2153 extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_int64_t step, const t_commrec *cr,
2154 t_mdatoms *md, t_state *state, gmx_update_t upd, gmx_constr_t constr)
2157 real rate = (ir->delta_t)/ir->opts.tau_t[0];
2159 if (ir->etc == etcANDERSEN && constr != NULL)
2161 /* Currently, Andersen thermostat does not support constrained
2162 systems. Functionality exists in the andersen_tcoupl
2163 function in GROMACS 4.5.7 to allow this combination. That
2164 code could be ported to the current random-number
2165 generation approach, but has not yet been done because of
2166 lack of time and resources. */
2167 gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
2170 /* proceed with andersen if 1) it's fixed probability per
2171 particle andersen or 2) it's massive andersen and it's tau_t/dt */
2172 if ((ir->etc == etcANDERSEN) || do_per_step(step, (int)(1.0/rate)))
2174 andersen_tcoupl(ir, step, cr, md, state, rate,
2175 upd->sd->randomize_group, upd->sd->boltzfac);