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.
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, double dt,
118 t_grp_tcstat *tcstat,
120 gmx_bool bNEMD, t_grp_acc *gstat, rvec accel[],
123 unsigned short ptype[], unsigned short cFREEZE[],
124 unsigned short cACC[], unsigned short cTC[],
125 rvec x[], rvec xprime[], rvec v[],
127 gmx_bool bNH, gmx_bool bPR)
130 int gf = 0, ga = 0, gt = 0;
132 real vn, vv, va, vb, vnrel;
138 /* Update with coupling to extended ensembles, used for
139 * Nose-Hoover and Parrinello-Rahman coupling
140 * Nose-Hoover uses the reversible leap-frog integrator from
141 * Holian et al. Phys Rev E 52(3) : 2338, 1995
143 for (n = start; n < nrend; n++)
158 lg = tcstat[gt].lambda;
163 rvec_sub(v[n], gstat[ga].u, vrel);
165 for (d = 0; d < DIM; d++)
167 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
169 vnrel = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
170 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
171 /* do not scale the mean velocities u */
172 vn = gstat[ga].u[d] + accel[ga][d]*dt + vnrel;
174 xprime[n][d] = x[n][d]+vn*dt;
179 xprime[n][d] = x[n][d];
184 else if (cFREEZE != NULL ||
185 nFreeze[0][XX] || nFreeze[0][YY] || nFreeze[0][ZZ] ||
188 /* Update with Berendsen/v-rescale coupling and freeze or NEMD */
189 for (n = start; n < nrend; n++)
191 w_dt = invmass[n]*dt;
204 lg = tcstat[gt].lambda;
206 for (d = 0; d < DIM; d++)
209 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
211 vv = lg*vn + f[n][d]*w_dt;
213 /* do not scale the mean velocities u */
215 va = vv + accel[ga][d]*dt;
216 vb = va + (1.0-lg)*u;
218 xprime[n][d] = x[n][d]+vb*dt;
223 xprime[n][d] = x[n][d];
230 /* Plain update with Berendsen/v-rescale coupling */
231 for (n = start; n < nrend; n++)
233 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
235 w_dt = invmass[n]*dt;
240 lg = tcstat[gt].lambda;
242 for (d = 0; d < DIM; d++)
244 vn = lg*v[n][d] + f[n][d]*w_dt;
246 xprime[n][d] = x[n][d] + vn*dt;
251 for (d = 0; d < DIM; d++)
254 xprime[n][d] = x[n][d];
261 static void do_update_vv_vel(int start, int nrend, double dt,
262 rvec accel[], ivec nFreeze[], real invmass[],
263 unsigned short ptype[], unsigned short cFREEZE[],
264 unsigned short cACC[], rvec v[], rvec f[],
265 gmx_bool bExtended, real veta, real alpha)
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)
319 /* Would it make more sense if Parrinello-Rahman was put here? */
324 mr2 = series_sinhx(g);
332 for (n = start; n < nrend; n++)
340 for (d = 0; d < DIM; d++)
342 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
344 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
348 xprime[n][d] = x[n][d];
352 } /* do_update_vv_pos */
354 static void do_update_visc(int start, int nrend, double dt,
355 t_grp_tcstat *tcstat,
358 unsigned short ptype[], unsigned short cTC[],
359 rvec x[], rvec xprime[], rvec v[],
360 rvec f[], matrix M, matrix box, real
361 cos_accel, real vcos,
362 gmx_bool bNH, gmx_bool bPR)
367 real lg, vxi = 0, vv;
372 fac = 2*M_PI/(box[ZZ][ZZ]);
376 /* Update with coupling to extended ensembles, used for
377 * Nose-Hoover and Parrinello-Rahman coupling
379 for (n = start; n < nrend; n++)
386 lg = tcstat[gt].lambda;
387 cosz = cos(fac*x[n][ZZ]);
389 copy_rvec(v[n], vrel);
397 for (d = 0; d < DIM; d++)
399 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
401 vn = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
402 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
405 vn += vc + dt*cosz*cos_accel;
408 xprime[n][d] = x[n][d]+vn*dt;
412 xprime[n][d] = x[n][d];
419 /* Classic version of update, used with berendsen coupling */
420 for (n = start; n < nrend; n++)
422 w_dt = invmass[n]*dt;
427 lg = tcstat[gt].lambda;
428 cosz = cos(fac*x[n][ZZ]);
430 for (d = 0; d < DIM; d++)
434 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
439 /* Do not scale the cosine velocity profile */
440 vv = vc + lg*(vn - vc + f[n][d]*w_dt);
441 /* Add the cosine accelaration profile */
442 vv += dt*cosz*cos_accel;
446 vv = lg*(vn + f[n][d]*w_dt);
449 xprime[n][d] = x[n][d]+vv*dt;
454 xprime[n][d] = x[n][d];
461 static gmx_stochd_t *init_stochd(t_inputrec *ir)
470 ngtc = ir->opts.ngtc;
474 snew(sd->bd_rf, ngtc);
476 else if (EI_SD(ir->eI))
479 snew(sd->sdsig, ngtc);
482 for (n = 0; n < ngtc; n++)
484 if (ir->opts.tau_t[n] > 0)
486 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
487 sdc[n].eph = exp(sdc[n].gdt/2);
488 sdc[n].emh = exp(-sdc[n].gdt/2);
489 sdc[n].em = exp(-sdc[n].gdt);
493 /* No friction and noise on this group */
499 if (sdc[n].gdt >= 0.05)
501 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
502 - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
503 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
504 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
509 /* Seventh order expansions for small y */
510 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
511 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))));
512 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
516 fprintf(debug, "SD const tc-grp %d: b %g c %g d %g\n",
517 n, sdc[n].b, sdc[n].c, sdc[n].d);
521 else if (ETC_ANDERSEN(ir->etc))
530 snew(sd->randomize_group, ngtc);
531 snew(sd->boltzfac, ngtc);
533 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
534 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
536 for (n = 0; n < ngtc; n++)
538 reft = std::max<real>(0, opts->ref_t[n]);
539 if ((opts->tau_t[n] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
541 sd->randomize_group[n] = TRUE;
542 sd->boltzfac[n] = BOLTZ*opts->ref_t[n];
546 sd->randomize_group[n] = FALSE;
553 gmx_update_t init_update(t_inputrec *ir)
559 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
561 upd->sd = init_stochd(ir);
570 static void do_update_sd1(gmx_stochd_t *sd,
571 int start, int nrend, double dt,
572 rvec accel[], ivec nFreeze[],
573 real invmass[], unsigned short ptype[],
574 unsigned short cFREEZE[], unsigned short cACC[],
575 unsigned short cTC[],
576 rvec x[], rvec xprime[], rvec v[], rvec f[],
577 int ngtc, real ref_t[],
579 gmx_bool bFirstHalfConstr,
580 gmx_int64_t step, int seed, int* gatindex)
585 int gf = 0, ga = 0, gt = 0;
592 for (n = 0; n < ngtc; n++)
595 /* The mass is encounted for later, since this differs per atom */
596 sig[n].V = sqrt(kT*(1 - sdc[n].em*sdc[n].em));
601 for (n = start; n < nrend; n++)
604 int ng = gatindex ? gatindex[n] : n;
606 ism = sqrt(invmass[n]);
620 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
622 for (d = 0; d < DIM; d++)
624 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
628 sd_V = ism*sig[gt].V*rnd[d];
629 vn = v[n][d] + (invmass[n]*f[n][d] + accel[ga][d])*dt;
630 v[n][d] = vn*sdc[gt].em + sd_V;
631 /* Here we include half of the friction+noise
632 * update of v into the integration of x.
634 xprime[n][d] = x[n][d] + 0.5*(vn + v[n][d])*dt;
639 xprime[n][d] = x[n][d];
646 /* We do have constraints */
647 if (bFirstHalfConstr)
649 /* First update without friction and noise */
652 for (n = start; n < nrend; n++)
665 for (d = 0; d < DIM; d++)
667 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
669 v[n][d] = v[n][d] + (im*f[n][d] + accel[ga][d])*dt;
670 xprime[n][d] = x[n][d] + v[n][d]*dt;
675 xprime[n][d] = x[n][d];
682 /* Update friction and noise only */
683 for (n = start; n < nrend; n++)
686 int ng = gatindex ? gatindex[n] : n;
688 ism = sqrt(invmass[n]);
698 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
700 for (d = 0; d < DIM; d++)
702 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
706 sd_V = ism*sig[gt].V*rnd[d];
708 v[n][d] = vn*sdc[gt].em + sd_V;
709 /* Add the friction and noise contribution only */
710 xprime[n][d] = xprime[n][d] + 0.5*(v[n][d] - vn)*dt;
718 static void check_sd2_work_data_allocation(gmx_stochd_t *sd, int nrend)
720 if (nrend > sd->sd_V_nalloc)
722 sd->sd_V_nalloc = over_alloc_dd(nrend);
723 srenew(sd->sd_V, sd->sd_V_nalloc);
727 static void do_update_sd2_Tconsts(gmx_stochd_t *sd,
732 /* This is separated from the update below, because it is single threaded */
741 for (gt = 0; gt < ngtc; gt++)
743 kT = BOLTZ*ref_t[gt];
744 /* The mass is encounted for later, since this differs per atom */
745 sig[gt].V = sqrt(kT*(1-sdc[gt].em));
746 sig[gt].X = sqrt(kT*sqr(tau_t[gt])*sdc[gt].c);
747 sig[gt].Yv = sqrt(kT*sdc[gt].b/sdc[gt].c);
748 sig[gt].Yx = sqrt(kT*sqr(tau_t[gt])*sdc[gt].b/(1-sdc[gt].em));
752 static void do_update_sd2(gmx_stochd_t *sd,
754 int start, int nrend,
755 rvec accel[], ivec nFreeze[],
756 real invmass[], unsigned short ptype[],
757 unsigned short cFREEZE[], unsigned short cACC[],
758 unsigned short cTC[],
759 rvec x[], rvec xprime[], rvec v[], rvec f[],
762 gmx_bool bFirstHalf, gmx_int64_t step, int seed,
767 /* The random part of the velocity update, generated in the first
768 * half of the update, needs to be remembered for the second half.
771 int gf = 0, ga = 0, gt = 0;
772 real vn = 0, Vmh, Xmh;
780 for (n = start; n < nrend; n++)
782 real rnd[6], rndi[3];
783 ng = gatindex ? gatindex[n] : n;
784 ism = sqrt(invmass[n]);
798 gmx_rng_cycle_6gaussian_table(step*2+(bFirstHalf ? 1 : 2), ng, seed, RND_SEED_UPDATE, rnd);
801 gmx_rng_cycle_3gaussian_table(step*2, ng, seed, RND_SEED_UPDATE, rndi);
803 for (d = 0; d < DIM; d++)
809 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
815 sd_X[n][d] = ism*sig[gt].X*rndi[d];
817 Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
818 + ism*sig[gt].Yv*rnd[d*2];
819 sd_V[n][d] = ism*sig[gt].V*rnd[d*2+1];
821 v[n][d] = vn*sdc[gt].em
822 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
823 + sd_V[n][d] - sdc[gt].em*Vmh;
825 xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
829 /* Correct the velocities for the constraints.
830 * This operation introduces some inaccuracy,
831 * since the velocity is determined from differences in coordinates.
834 (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
836 Xmh = sd_V[n][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
837 + ism*sig[gt].Yx*rnd[d*2];
838 sd_X[n][d] = ism*sig[gt].X*rnd[d*2+1];
840 xprime[n][d] += sd_X[n][d] - Xmh;
849 xprime[n][d] = x[n][d];
856 static void do_update_bd_Tconsts(double dt, real friction_coefficient,
857 int ngtc, const real ref_t[],
860 /* This is separated from the update below, because it is single threaded */
863 if (friction_coefficient != 0)
865 for (gt = 0; gt < ngtc; gt++)
867 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]/(friction_coefficient*dt));
872 for (gt = 0; gt < ngtc; gt++)
874 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]);
879 static void do_update_bd(int start, int nrend, double dt,
881 real invmass[], unsigned short ptype[],
882 unsigned short cFREEZE[], unsigned short cTC[],
883 rvec x[], rvec xprime[], rvec v[],
884 rvec f[], real friction_coefficient,
885 real *rf, gmx_int64_t step, int seed,
888 /* note -- these appear to be full step velocities . . . */
894 if (friction_coefficient != 0)
896 invfr = 1.0/friction_coefficient;
899 for (n = start; (n < nrend); n++)
902 int ng = gatindex ? gatindex[n] : n;
912 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
913 for (d = 0; (d < DIM); d++)
915 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
917 if (friction_coefficient != 0)
919 vn = invfr*f[n][d] + rf[gt]*rnd[d];
923 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
924 vn = 0.5*invmass[n]*f[n][d]*dt
925 + sqrt(0.5*invmass[n])*rf[gt]*rnd[d];
929 xprime[n][d] = x[n][d]+vn*dt;
934 xprime[n][d] = x[n][d];
940 static void dump_it_all(FILE gmx_unused *fp, const char gmx_unused *title,
941 int gmx_unused natoms, rvec gmx_unused x[], rvec gmx_unused xp[],
942 rvec gmx_unused v[], rvec gmx_unused f[])
947 fprintf(fp, "%s\n", title);
948 pr_rvecs(fp, 0, "x", x, natoms);
949 pr_rvecs(fp, 0, "xp", xp, natoms);
950 pr_rvecs(fp, 0, "v", v, natoms);
951 pr_rvecs(fp, 0, "f", f, natoms);
956 static void calc_ke_part_normal(rvec v[], t_grpopts *opts, t_mdatoms *md,
957 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel)
960 t_grp_tcstat *tcstat = ekind->tcstat;
961 t_grp_acc *grpstat = ekind->grpstat;
964 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
965 an option, but not supported now.
966 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
969 /* group velocities are calculated in update_ekindata and
970 * accumulated in acumulate_groups.
971 * Now the partial global and groups ekin.
973 for (g = 0; (g < opts->ngtc); g++)
975 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
978 clear_mat(tcstat[g].ekinf);
979 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
983 clear_mat(tcstat[g].ekinh);
986 ekind->dekindl_old = ekind->dekindl;
987 nthread = gmx_omp_nthreads_get(emntUpdate);
989 #pragma omp parallel for num_threads(nthread) schedule(static)
990 for (thread = 0; thread < nthread; thread++)
992 int start_t, end_t, n;
1000 start_t = ((thread+0)*md->homenr)/nthread;
1001 end_t = ((thread+1)*md->homenr)/nthread;
1003 ekin_sum = ekind->ekin_work[thread];
1004 dekindl_sum = ekind->dekindl_work[thread];
1006 for (gt = 0; gt < opts->ngtc; gt++)
1008 clear_mat(ekin_sum[gt]);
1014 for (n = start_t; n < end_t; n++)
1024 hm = 0.5*md->massT[n];
1026 for (d = 0; (d < DIM); d++)
1028 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
1030 for (d = 0; (d < DIM); d++)
1032 for (m = 0; (m < DIM); m++)
1034 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1035 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
1038 if (md->nMassPerturbed && md->bPerturbed[n])
1041 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1047 for (thread = 0; thread < nthread; thread++)
1049 for (g = 0; g < opts->ngtc; g++)
1053 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
1058 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1063 ekind->dekindl += *ekind->dekindl_work[thread];
1066 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1069 static void calc_ke_part_visc(matrix box, rvec x[], rvec v[],
1070 t_grpopts *opts, t_mdatoms *md,
1071 gmx_ekindata_t *ekind,
1072 t_nrnb *nrnb, gmx_bool bEkinAveVel)
1074 int start = 0, homenr = md->homenr;
1075 int g, d, n, m, gt = 0;
1078 t_grp_tcstat *tcstat = ekind->tcstat;
1079 t_cos_acc *cosacc = &(ekind->cosacc);
1084 for (g = 0; g < opts->ngtc; g++)
1086 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1087 clear_mat(ekind->tcstat[g].ekinh);
1089 ekind->dekindl_old = ekind->dekindl;
1091 fac = 2*M_PI/box[ZZ][ZZ];
1094 for (n = start; n < start+homenr; n++)
1100 hm = 0.5*md->massT[n];
1102 /* Note that the times of x and v differ by half a step */
1103 /* MRS -- would have to be changed for VV */
1104 cosz = cos(fac*x[n][ZZ]);
1105 /* Calculate the amplitude of the new velocity profile */
1106 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1108 copy_rvec(v[n], v_corrt);
1109 /* Subtract the profile for the kinetic energy */
1110 v_corrt[XX] -= cosz*cosacc->vcos;
1111 for (d = 0; (d < DIM); d++)
1113 for (m = 0; (m < DIM); m++)
1115 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1118 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1122 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1126 if (md->nPerturbed && md->bPerturbed[n])
1128 /* The minus sign here might be confusing.
1129 * The kinetic contribution from dH/dl doesn't come from
1130 * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
1131 * where p are the momenta. The difference is only a minus sign.
1133 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1136 ekind->dekindl = dekindl;
1137 cosacc->mvcos = mvcos;
1139 inc_nrnb(nrnb, eNR_EKIN, homenr);
1142 void calc_ke_part(t_state *state, t_grpopts *opts, t_mdatoms *md,
1143 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel)
1145 if (ekind->cosacc.cos_accel == 0)
1147 calc_ke_part_normal(state->v, opts, md, ekind, nrnb, bEkinAveVel);
1151 calc_ke_part_visc(state->box, state->x, state->v, opts, md, ekind, nrnb, bEkinAveVel);
1155 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1157 ekinstate->ekin_n = ir->opts.ngtc;
1158 snew(ekinstate->ekinh, ekinstate->ekin_n);
1159 snew(ekinstate->ekinf, ekinstate->ekin_n);
1160 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1161 snew(ekinstate->ekinscalef_nhc, ekinstate->ekin_n);
1162 snew(ekinstate->ekinscaleh_nhc, ekinstate->ekin_n);
1163 snew(ekinstate->vscale_nhc, ekinstate->ekin_n);
1164 ekinstate->dekindl = 0;
1165 ekinstate->mvcos = 0;
1168 void update_ekinstate(ekinstate_t *ekinstate, gmx_ekindata_t *ekind)
1172 for (i = 0; i < ekinstate->ekin_n; i++)
1174 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1175 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1176 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1177 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1178 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1179 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1182 copy_mat(ekind->ekin, ekinstate->ekin_total);
1183 ekinstate->dekindl = ekind->dekindl;
1184 ekinstate->mvcos = ekind->cosacc.mvcos;
1188 void restore_ekinstate_from_state(t_commrec *cr,
1189 gmx_ekindata_t *ekind, ekinstate_t *ekinstate)
1195 for (i = 0; i < ekinstate->ekin_n; i++)
1197 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1198 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1199 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1200 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1201 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1202 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1205 copy_mat(ekinstate->ekin_total, ekind->ekin);
1207 ekind->dekindl = ekinstate->dekindl;
1208 ekind->cosacc.mvcos = ekinstate->mvcos;
1209 n = ekinstate->ekin_n;
1214 gmx_bcast(sizeof(n), &n, cr);
1215 for (i = 0; i < n; i++)
1217 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1218 ekind->tcstat[i].ekinh[0], cr);
1219 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1220 ekind->tcstat[i].ekinf[0], cr);
1221 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1222 ekind->tcstat[i].ekinh_old[0], cr);
1224 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1225 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1226 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1227 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1228 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1229 &(ekind->tcstat[i].vscale_nhc), cr);
1231 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1232 ekind->ekin[0], cr);
1234 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1235 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1239 void set_deform_reference_box(gmx_update_t upd, gmx_int64_t step, matrix box)
1241 upd->deformref_step = step;
1242 copy_mat(box, upd->deformref_box);
1245 static void deform(gmx_update_t upd,
1246 int start, int homenr, rvec x[], matrix box,
1247 const t_inputrec *ir, gmx_int64_t step)
1249 matrix bnew, invbox, mu;
1253 elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1254 copy_mat(box, bnew);
1255 for (i = 0; i < DIM; i++)
1257 for (j = 0; j < DIM; j++)
1259 if (ir->deform[i][j] != 0)
1262 upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1266 /* We correct the off-diagonal elements,
1267 * which can grow indefinitely during shearing,
1268 * so the shifts do not get messed up.
1270 for (i = 1; i < DIM; i++)
1272 for (j = i-1; j >= 0; j--)
1274 while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1276 rvec_dec(bnew[i], bnew[j]);
1278 while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1280 rvec_inc(bnew[i], bnew[j]);
1284 m_inv_ur0(box, invbox);
1285 copy_mat(bnew, box);
1286 mmul_ur0(box, invbox, mu);
1288 for (i = start; i < start+homenr; i++)
1290 x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1291 x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1292 x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1296 void update_tcouple(gmx_int64_t step,
1297 t_inputrec *inputrec,
1299 gmx_ekindata_t *ekind,
1304 gmx_bool bTCouple = FALSE;
1308 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1309 if (inputrec->etc != etcNO &&
1310 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec)))
1312 /* We should only couple after a step where energies were determined (for leapfrog versions)
1313 or the step energies are determined, for velocity verlet versions */
1315 if (EI_VV(inputrec->eI))
1323 bTCouple = (inputrec->nsttcouple == 1 ||
1324 do_per_step(step+inputrec->nsttcouple-offset,
1325 inputrec->nsttcouple));
1330 dttc = inputrec->nsttcouple*inputrec->delta_t;
1332 switch (inputrec->etc)
1337 berendsen_tcoupl(inputrec, ekind, dttc);
1340 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1341 state->nosehoover_xi, state->nosehoover_vxi, MassQ);
1344 vrescale_tcoupl(inputrec, step, ekind, dttc,
1345 state->therm_integral);
1348 /* rescale in place here */
1349 if (EI_VV(inputrec->eI))
1351 rescale_velocities(ekind, md, 0, md->homenr, state->v);
1356 /* Set the T scaling lambda to 1 to have no scaling */
1357 for (i = 0; (i < inputrec->opts.ngtc); i++)
1359 ekind->tcstat[i].lambda = 1.0;
1364 void update_pcouple(FILE *fplog,
1366 t_inputrec *inputrec,
1372 gmx_bool bPCouple = FALSE;
1376 /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1377 if (inputrec->epc != epcNO && (!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))))
1379 /* We should only couple after a step where energies were determined */
1380 bPCouple = (inputrec->nstpcouple == 1 ||
1381 do_per_step(step+inputrec->nstpcouple-1,
1382 inputrec->nstpcouple));
1385 clear_mat(pcoupl_mu);
1386 for (i = 0; i < DIM; i++)
1388 pcoupl_mu[i][i] = 1.0;
1395 dtpc = inputrec->nstpcouple*inputrec->delta_t;
1397 switch (inputrec->epc)
1399 /* We can always pcoupl, even if we did not sum the energies
1400 * the previous step, since state->pres_prev is only updated
1401 * when the energies have been summed.
1405 case (epcBERENDSEN):
1408 berendsen_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev, state->box,
1412 case (epcPARRINELLORAHMAN):
1413 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1414 state->box, state->box_rel, state->boxv,
1415 M, pcoupl_mu, bInitStep);
1423 static rvec *get_xprime(const t_state *state, gmx_update_t upd)
1425 if (state->nalloc > upd->xp_nalloc)
1427 upd->xp_nalloc = state->nalloc;
1428 srenew(upd->xp, upd->xp_nalloc);
1434 static void combine_forces(gmx_update_t upd,
1436 gmx_constr_t constr,
1437 t_inputrec *ir, t_mdatoms *md, t_idef *idef,
1440 t_state *state, gmx_bool bMolPBC,
1441 int start, int nrend,
1442 rvec f[], rvec f_lr[],
1443 tensor *vir_lr_constr,
1448 /* f contains the short-range forces + the long range forces
1449 * which are stored separately in f_lr.
1452 if (constr != NULL && vir_lr_constr != NULL &&
1453 !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1455 /* We need to constrain the LR forces separately,
1456 * because due to the different pre-factor for the SR and LR
1457 * forces in the update algorithm, we have to correct
1458 * the constraint virial for the nstcalclr-1 extra f_lr.
1459 * Constrain only the additional LR part of the force.
1461 /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1466 xp = get_xprime(state, upd);
1468 fac = (nstcalclr - 1)*ir->delta_t*ir->delta_t;
1470 for (i = 0; i < md->homenr; i++)
1472 if (md->cFREEZE != NULL)
1474 gf = md->cFREEZE[i];
1476 for (d = 0; d < DIM; d++)
1478 if ((md->ptype[i] != eptVSite) &&
1479 (md->ptype[i] != eptShell) &&
1480 !ir->opts.nFreeze[gf][d])
1482 xp[i][d] = state->x[i][d] + fac*f_lr[i][d]*md->invmass[i];
1486 xp[i][d] = state->x[i][d];
1490 constrain(NULL, FALSE, FALSE, constr, idef, ir, cr, step, 0, 1.0, md,
1491 state->x, xp, xp, bMolPBC, state->box, state->lambda[efptBONDED], NULL,
1492 NULL, vir_lr_constr, nrnb, econqForce);
1495 /* Add nstcalclr-1 times the LR force to the sum of both forces
1496 * and store the result in forces_lr.
1498 for (i = start; i < nrend; i++)
1500 for (d = 0; d < DIM; d++)
1502 f_lr[i][d] = f[i][d] + (nstcalclr - 1)*f_lr[i][d];
1507 void update_constraints(FILE *fplog,
1509 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1510 t_inputrec *inputrec, /* input record and box stuff */
1515 rvec force[], /* forces on home particles */
1520 gmx_wallcycle_t wcycle,
1522 gmx_constr_t constr,
1523 gmx_bool bFirstHalf,
1526 gmx_bool bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1528 int start, homenr, nrend, i, m;
1530 rvec *xprime = NULL;
1537 if (bFirstHalf && !EI_VV(inputrec->eI))
1542 /* for now, SD update is here -- though it really seems like it
1543 should be reformulated as a velocity verlet method, since it has two parts */
1546 homenr = md->homenr;
1547 nrend = start+homenr;
1549 dt = inputrec->delta_t;
1553 * APPLY CONSTRAINTS:
1556 * When doing PR pressure coupling we have to constrain the
1557 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1558 * it is enough to do this once though, since the relative velocities
1559 * after this will be normal to the bond vector
1564 /* clear out constraints before applying */
1565 clear_mat(vir_part);
1567 xprime = get_xprime(state, upd);
1569 bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1570 bLog = (do_per_step(step, inputrec->nstlog) || bLastStep || (step < 0));
1571 bEner = (do_per_step(step, inputrec->nstenergy) || bLastStep);
1572 /* Constrain the coordinates xprime */
1573 wallcycle_start(wcycle, ewcCONSTR);
1574 if (EI_VV(inputrec->eI) && bFirstHalf)
1576 constrain(NULL, bLog, bEner, constr, idef,
1577 inputrec, cr, step, 1, 1.0, md,
1578 state->x, state->v, state->v,
1579 bMolPBC, state->box,
1580 state->lambda[efptBONDED], dvdlambda,
1581 NULL, bCalcVir ? &vir_con : NULL, nrnb, econqVeloc);
1585 constrain(NULL, bLog, bEner, constr, idef,
1586 inputrec, cr, step, 1, 1.0, md,
1587 state->x, xprime, NULL,
1588 bMolPBC, state->box,
1589 state->lambda[efptBONDED], dvdlambda,
1590 state->v, bCalcVir ? &vir_con : NULL, nrnb, econqCoord);
1592 wallcycle_stop(wcycle, ewcCONSTR);
1596 dump_it_all(fplog, "After Shake",
1597 state->natoms, state->x, xprime, state->v, force);
1601 if (inputrec->eI == eiSD2)
1603 /* A correction factor eph is needed for the SD constraint force */
1604 /* Here we can, unfortunately, not have proper corrections
1605 * for different friction constants, so we use the first one.
1607 for (i = 0; i < DIM; i++)
1609 for (m = 0; m < DIM; m++)
1611 vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1617 m_add(vir_part, vir_con, vir_part);
1621 pr_rvecs(debug, 0, "constraint virial", vir_part, DIM);
1628 if (inputrec->eI == eiSD1 && bDoConstr && !bFirstHalf)
1630 wallcycle_start(wcycle, ewcUPDATE);
1631 xprime = get_xprime(state, upd);
1633 nth = gmx_omp_nthreads_get(emntUpdate);
1635 #pragma omp parallel for num_threads(nth) schedule(static)
1637 for (th = 0; th < nth; th++)
1639 int start_th, end_th;
1641 start_th = start + ((nrend-start)* th )/nth;
1642 end_th = start + ((nrend-start)*(th+1))/nth;
1644 /* The second part of the SD integration */
1645 do_update_sd1(upd->sd,
1646 start_th, end_th, dt,
1647 inputrec->opts.acc, inputrec->opts.nFreeze,
1648 md->invmass, md->ptype,
1649 md->cFREEZE, md->cACC, md->cTC,
1650 state->x, xprime, state->v, force,
1651 inputrec->opts.ngtc, inputrec->opts.ref_t,
1653 step, inputrec->ld_seed,
1654 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1656 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1657 wallcycle_stop(wcycle, ewcUPDATE);
1661 /* Constrain the coordinates xprime for half a time step */
1662 wallcycle_start(wcycle, ewcCONSTR);
1664 constrain(NULL, bLog, bEner, constr, idef,
1665 inputrec, cr, step, 1, 0.5, md,
1666 state->x, xprime, NULL,
1667 bMolPBC, state->box,
1668 state->lambda[efptBONDED], dvdlambda,
1669 state->v, NULL, nrnb, econqCoord);
1671 wallcycle_stop(wcycle, ewcCONSTR);
1675 if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1677 wallcycle_start(wcycle, ewcUPDATE);
1678 xprime = get_xprime(state, upd);
1680 nth = gmx_omp_nthreads_get(emntUpdate);
1682 #pragma omp parallel for num_threads(nth) schedule(static)
1683 for (th = 0; th < nth; th++)
1685 int start_th, end_th;
1687 start_th = start + ((nrend-start)* th )/nth;
1688 end_th = start + ((nrend-start)*(th+1))/nth;
1690 /* The second part of the SD integration */
1691 do_update_sd2(upd->sd,
1692 FALSE, start_th, end_th,
1693 inputrec->opts.acc, inputrec->opts.nFreeze,
1694 md->invmass, md->ptype,
1695 md->cFREEZE, md->cACC, md->cTC,
1696 state->x, xprime, state->v, force, state->sd_X,
1697 inputrec->opts.tau_t,
1698 FALSE, step, inputrec->ld_seed,
1699 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1701 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1702 wallcycle_stop(wcycle, ewcUPDATE);
1706 /* Constrain the coordinates xprime */
1707 wallcycle_start(wcycle, ewcCONSTR);
1708 constrain(NULL, bLog, bEner, constr, idef,
1709 inputrec, cr, step, 1, 1.0, md,
1710 state->x, xprime, NULL,
1711 bMolPBC, state->box,
1712 state->lambda[efptBONDED], dvdlambda,
1713 NULL, NULL, nrnb, econqCoord);
1714 wallcycle_stop(wcycle, ewcCONSTR);
1719 /* We must always unshift after updating coordinates; if we did not shake
1720 x was shifted in do_force */
1722 if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1724 if (graph && (graph->nnodes > 0))
1726 unshift_x(graph, state->box, state->x, upd->xp);
1727 if (TRICLINIC(state->box))
1729 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1733 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1738 #ifndef __clang_analyzer__
1739 // cppcheck-suppress unreadVariable
1740 nth = gmx_omp_nthreads_get(emntUpdate);
1742 #pragma omp parallel for num_threads(nth) schedule(static)
1743 for (i = start; i < nrend; i++)
1745 copy_rvec(upd->xp[i], state->x[i]);
1749 dump_it_all(fplog, "After unshift",
1750 state->natoms, state->x, upd->xp, state->v, force);
1752 /* ############# END the update of velocities and positions ######### */
1755 void update_box(FILE *fplog,
1757 t_inputrec *inputrec, /* input record and box stuff */
1760 rvec force[], /* forces on home particles */
1766 int start, homenr, i, n, m;
1769 homenr = md->homenr;
1771 dt = inputrec->delta_t;
1775 /* now update boxes */
1776 switch (inputrec->epc)
1780 case (epcBERENDSEN):
1781 berendsen_pscale(inputrec, pcoupl_mu, state->box, state->box_rel,
1782 start, homenr, state->x, md->cFREEZE, nrnb);
1784 case (epcPARRINELLORAHMAN):
1785 /* The box velocities were updated in do_pr_pcoupl in the update
1786 * iteration, but we dont change the box vectors until we get here
1787 * since we need to be able to shift/unshift above.
1789 for (i = 0; i < DIM; i++)
1791 for (m = 0; m <= i; m++)
1793 state->box[i][m] += dt*state->boxv[i][m];
1796 preserve_box_shape(inputrec, state->box_rel, state->box);
1798 /* Scale the coordinates */
1799 for (n = start; (n < start+homenr); n++)
1801 tmvmul_ur0(pcoupl_mu, state->x[n], state->x[n]);
1805 switch (inputrec->epct)
1807 case (epctISOTROPIC):
1808 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1809 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1810 Side length scales as exp(veta*dt) */
1812 msmul(state->box, exp(state->veta*dt), state->box);
1814 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1815 o If we assume isotropic scaling, and box length scaling
1816 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1817 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1818 determinant of B is L^DIM det(M), and the determinant
1819 of dB/dt is (dL/dT)^DIM det (M). veta will be
1820 (det(dB/dT)/det(B))^(1/3). Then since M =
1821 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1823 msmul(state->box, state->veta, state->boxv);
1833 if (DEFORM(*inputrec))
1835 deform(upd, start, homenr, state->x, state->box, inputrec, step);
1838 dump_it_all(fplog, "After update",
1839 state->natoms, state->x, upd->xp, state->v, force);
1842 void update_coords(FILE *fplog,
1844 t_inputrec *inputrec, /* input record and box stuff */
1848 rvec *f, /* forces on home particles */
1851 tensor *vir_lr_constr,
1853 gmx_ekindata_t *ekind,
1858 t_commrec *cr, /* these shouldn't be here -- need to think about it */
1860 gmx_constr_t constr,
1863 gmx_bool bNH, bPR, bDoConstr = FALSE;
1866 int start, homenr, nrend;
1870 bDoConstr = (NULL != constr);
1872 /* Running the velocity half does nothing except for velocity verlet */
1873 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1874 !EI_VV(inputrec->eI))
1876 gmx_incons("update_coords called for velocity without VV integrator");
1880 homenr = md->homenr;
1881 nrend = start+homenr;
1883 xprime = get_xprime(state, upd);
1885 dt = inputrec->delta_t;
1887 /* We need to update the NMR restraint history when time averaging is used */
1888 if (state->flags & (1<<estDISRE_RM3TAV))
1890 update_disres_history(fcd, &state->hist);
1892 if (state->flags & (1<<estORIRE_DTAV))
1894 update_orires_history(fcd, &state->hist);
1898 bNH = inputrec->etc == etcNOSEHOOVER;
1899 bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1901 if (bDoLR && inputrec->nstcalclr > 1 && !EI_VV(inputrec->eI)) /* get this working with VV? */
1903 /* Store the total force + nstcalclr-1 times the LR force
1904 * in forces_lr, so it can be used in a normal update algorithm
1905 * to produce twin time stepping.
1907 /* is this correct in the new construction? MRS */
1909 inputrec->nstcalclr, constr, inputrec, md, idef, cr,
1910 step, state, bMolPBC,
1911 start, nrend, f, f_lr, vir_lr_constr, nrnb);
1919 /* ############# START The update of velocities and positions ######### */
1921 dump_it_all(fplog, "Before update",
1922 state->natoms, state->x, xprime, state->v, force);
1924 if (inputrec->eI == eiSD2)
1926 check_sd2_work_data_allocation(upd->sd, nrend);
1928 do_update_sd2_Tconsts(upd->sd,
1929 inputrec->opts.ngtc,
1930 inputrec->opts.tau_t,
1931 inputrec->opts.ref_t);
1933 if (inputrec->eI == eiBD)
1935 do_update_bd_Tconsts(dt, inputrec->bd_fric,
1936 inputrec->opts.ngtc, inputrec->opts.ref_t,
1940 nth = gmx_omp_nthreads_get(emntUpdate);
1942 #pragma omp parallel for num_threads(nth) schedule(static) private(alpha)
1943 for (th = 0; th < nth; th++)
1945 int start_th, end_th;
1947 start_th = start + ((nrend-start)* th )/nth;
1948 end_th = start + ((nrend-start)*(th+1))/nth;
1950 switch (inputrec->eI)
1953 if (ekind->cosacc.cos_accel == 0)
1955 do_update_md(start_th, end_th, dt,
1956 ekind->tcstat, state->nosehoover_vxi,
1957 ekind->bNEMD, ekind->grpstat, inputrec->opts.acc,
1958 inputrec->opts.nFreeze,
1959 md->invmass, md->ptype,
1960 md->cFREEZE, md->cACC, md->cTC,
1961 state->x, xprime, state->v, force, M,
1966 do_update_visc(start_th, end_th, dt,
1967 ekind->tcstat, state->nosehoover_vxi,
1968 md->invmass, md->ptype,
1969 md->cTC, state->x, xprime, state->v, force, M,
1971 ekind->cosacc.cos_accel,
1977 /* With constraints, the SD1 update is done in 2 parts */
1978 do_update_sd1(upd->sd,
1979 start_th, end_th, dt,
1980 inputrec->opts.acc, inputrec->opts.nFreeze,
1981 md->invmass, md->ptype,
1982 md->cFREEZE, md->cACC, md->cTC,
1983 state->x, xprime, state->v, force,
1984 inputrec->opts.ngtc, inputrec->opts.ref_t,
1986 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1989 /* The SD2 update is always done in 2 parts,
1990 * because an extra constraint step is needed
1992 do_update_sd2(upd->sd,
1993 bInitStep, start_th, end_th,
1994 inputrec->opts.acc, inputrec->opts.nFreeze,
1995 md->invmass, md->ptype,
1996 md->cFREEZE, md->cACC, md->cTC,
1997 state->x, xprime, state->v, force, state->sd_X,
1998 inputrec->opts.tau_t,
1999 TRUE, step, inputrec->ld_seed,
2000 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2003 do_update_bd(start_th, end_th, dt,
2004 inputrec->opts.nFreeze, md->invmass, md->ptype,
2005 md->cFREEZE, md->cTC,
2006 state->x, xprime, state->v, force,
2009 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2013 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
2018 do_update_vv_vel(start_th, end_th, dt,
2019 inputrec->opts.acc, inputrec->opts.nFreeze,
2020 md->invmass, md->ptype,
2021 md->cFREEZE, md->cACC,
2023 (bNH || bPR), state->veta, alpha);
2026 do_update_vv_pos(start_th, end_th, dt,
2027 inputrec->opts.nFreeze,
2028 md->ptype, md->cFREEZE,
2029 state->x, xprime, state->v,
2030 (bNH || bPR), state->veta);
2035 gmx_fatal(FARGS, "Don't know how to update coordinates");
2043 void correct_ekin(FILE *log, int start, int end, rvec v[], rvec vcm, real mass[],
2044 real tmass, tensor ekin)
2047 * This is a debugging routine. It should not be called for production code
2049 * The kinetic energy should calculated according to:
2050 * Ekin = 1/2 m (v-vcm)^2
2051 * However the correction is not always applied, since vcm may not be
2052 * known in time and we compute
2053 * Ekin' = 1/2 m v^2 instead
2054 * This can be corrected afterwards by computing
2055 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
2057 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
2064 /* Local particles */
2067 /* Processor dependent part. */
2069 for (i = start; (i < end); i++)
2073 for (j = 0; (j < DIM); j++)
2079 svmul(1/tmass, vcm, vcm);
2080 svmul(0.5, vcm, hvcm);
2082 for (j = 0; (j < DIM); j++)
2084 for (k = 0; (k < DIM); k++)
2086 dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
2089 pr_rvecs(log, 0, "dekin", dekin, DIM);
2090 pr_rvecs(log, 0, " ekin", ekin, DIM);
2091 fprintf(log, "dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
2092 trace(dekin), trace(ekin), vcm[XX], vcm[YY], vcm[ZZ]);
2093 fprintf(log, "mv = (%8.4f %8.4f %8.4f)\n",
2094 mv[XX], mv[YY], mv[ZZ]);
2097 extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_int64_t step, const t_commrec *cr,
2098 t_mdatoms *md, t_state *state, gmx_update_t upd, gmx_constr_t constr)
2101 real rate = (ir->delta_t)/ir->opts.tau_t[0];
2103 if (ir->etc == etcANDERSEN && constr != NULL)
2105 gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
2108 /* proceed with andersen if 1) it's fixed probability per
2109 particle andersen or 2) it's massive andersen and it's tau_t/dt */
2110 if ((ir->etc == etcANDERSEN) || do_per_step(step, (int)(1.0/rate)))
2112 andersen_tcoupl(ir, step, cr, md, state, rate,
2113 upd->sd->randomize_group, upd->sd->boltzfac);