2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "types/commrec.h"
55 #include "gromacs/random/random.h"
65 #include "gmx_omp_nthreads.h"
67 #include "gromacs/fileio/confio.h"
68 #include "gromacs/fileio/futil.h"
69 #include "gromacs/timing/wallcycle.h"
70 #include "gromacs/utility/gmxomp.h"
71 #include "gromacs/pulling/pull.h"
73 /*For debugging, start at v(-dt/2) for velolcity verlet -- uncomment next line */
74 /*#define STARTFROMDT2*/
94 /* The random state for ngaussrand threads.
95 * Normal thermostats need just 1 random number generator,
96 * but SD and BD with OpenMP parallelization need 1 for each thread.
104 gmx_sd_sigma_t *sdsig;
107 /* andersen temperature control stuff */
108 gmx_bool *randomize_group;
112 typedef struct gmx_update
115 /* xprime for constraint algorithms */
119 /* variable size arrays for andersen */
122 gmx_bool randatom_list_init;
124 /* Variables for the deform algorithm */
125 gmx_int64_t deformref_step;
126 matrix deformref_box;
130 static void do_update_md(int start, int nrend, double dt,
131 t_grp_tcstat *tcstat,
133 gmx_bool bNEMD, t_grp_acc *gstat, rvec accel[],
136 unsigned short ptype[], unsigned short cFREEZE[],
137 unsigned short cACC[], unsigned short cTC[],
138 rvec x[], rvec xprime[], rvec v[],
140 gmx_bool bNH, gmx_bool bPR)
143 int gf = 0, ga = 0, gt = 0;
145 real vn, vv, va, vb, vnrel;
151 /* Update with coupling to extended ensembles, used for
152 * Nose-Hoover and Parrinello-Rahman coupling
153 * Nose-Hoover uses the reversible leap-frog integrator from
154 * Holian et al. Phys Rev E 52(3) : 2338, 1995
156 for (n = start; n < nrend; n++)
171 lg = tcstat[gt].lambda;
176 rvec_sub(v[n], gstat[ga].u, vrel);
178 for (d = 0; d < DIM; d++)
180 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
182 vnrel = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
183 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
184 /* do not scale the mean velocities u */
185 vn = gstat[ga].u[d] + accel[ga][d]*dt + vnrel;
187 xprime[n][d] = x[n][d]+vn*dt;
192 xprime[n][d] = x[n][d];
197 else if (cFREEZE != NULL ||
198 nFreeze[0][XX] || nFreeze[0][YY] || nFreeze[0][ZZ] ||
201 /* Update with Berendsen/v-rescale coupling and freeze or NEMD */
202 for (n = start; n < nrend; n++)
204 w_dt = invmass[n]*dt;
217 lg = tcstat[gt].lambda;
219 for (d = 0; d < DIM; d++)
222 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
224 vv = lg*vn + f[n][d]*w_dt;
226 /* do not scale the mean velocities u */
228 va = vv + accel[ga][d]*dt;
229 vb = va + (1.0-lg)*u;
231 xprime[n][d] = x[n][d]+vb*dt;
236 xprime[n][d] = x[n][d];
243 /* Plain update with Berendsen/v-rescale coupling */
244 for (n = start; n < nrend; n++)
246 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
248 w_dt = invmass[n]*dt;
253 lg = tcstat[gt].lambda;
255 for (d = 0; d < DIM; d++)
257 vn = lg*v[n][d] + f[n][d]*w_dt;
259 xprime[n][d] = x[n][d] + vn*dt;
264 for (d = 0; d < DIM; d++)
267 xprime[n][d] = x[n][d];
274 static void do_update_vv_vel(int start, int nrend, double dt,
275 rvec accel[], ivec nFreeze[], real invmass[],
276 unsigned short ptype[], unsigned short cFREEZE[],
277 unsigned short cACC[], rvec v[], rvec f[],
278 gmx_bool bExtended, real veta, real alpha)
283 real u, vn, vv, va, vb, vnrel;
289 g = 0.25*dt*veta*alpha;
291 mv2 = series_sinhx(g);
298 for (n = start; n < nrend; n++)
300 w_dt = invmass[n]*dt;
310 for (d = 0; d < DIM; d++)
312 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
314 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
322 } /* do_update_vv_vel */
324 static void do_update_vv_pos(int start, int nrend, double dt,
326 unsigned short ptype[], unsigned short cFREEZE[],
327 rvec x[], rvec xprime[], rvec v[],
328 gmx_bool bExtended, real veta)
335 /* Would it make more sense if Parrinello-Rahman was put here? */
340 mr2 = series_sinhx(g);
348 for (n = start; n < nrend; n++)
356 for (d = 0; d < DIM; d++)
358 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
360 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
364 xprime[n][d] = x[n][d];
368 } /* do_update_vv_pos */
370 static void do_update_visc(int start, int nrend, double dt,
371 t_grp_tcstat *tcstat,
374 unsigned short ptype[], unsigned short cTC[],
375 rvec x[], rvec xprime[], rvec v[],
376 rvec f[], matrix M, matrix box, real
377 cos_accel, real vcos,
378 gmx_bool bNH, gmx_bool bPR)
383 real lg, vxi = 0, vv;
388 fac = 2*M_PI/(box[ZZ][ZZ]);
392 /* Update with coupling to extended ensembles, used for
393 * Nose-Hoover and Parrinello-Rahman coupling
395 for (n = start; n < nrend; n++)
402 lg = tcstat[gt].lambda;
403 cosz = cos(fac*x[n][ZZ]);
405 copy_rvec(v[n], vrel);
413 for (d = 0; d < DIM; d++)
417 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
419 vn = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
420 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
423 vn += vc + dt*cosz*cos_accel;
426 xprime[n][d] = x[n][d]+vn*dt;
430 xprime[n][d] = x[n][d];
437 /* Classic version of update, used with berendsen coupling */
438 for (n = start; n < nrend; n++)
440 w_dt = invmass[n]*dt;
445 lg = tcstat[gt].lambda;
446 cosz = cos(fac*x[n][ZZ]);
448 for (d = 0; d < DIM; d++)
452 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
457 /* Do not scale the cosine velocity profile */
458 vv = vc + lg*(vn - vc + f[n][d]*w_dt);
459 /* Add the cosine accelaration profile */
460 vv += dt*cosz*cos_accel;
464 vv = lg*(vn + f[n][d]*w_dt);
467 xprime[n][d] = x[n][d]+vv*dt;
472 xprime[n][d] = x[n][d];
479 /* Allocates and initializes sd->gaussrand[i] for i=1, i<sd->ngaussrand,
480 * Using seeds generated from sd->gaussrand[0].
482 static void init_multiple_gaussrand(gmx_stochd_t *sd)
487 ngr = sd->ngaussrand;
490 for (i = 1; i < ngr; i++)
492 seed[i] = gmx_rng_uniform_uint32(sd->gaussrand[0]);
495 if (ngr != gmx_omp_nthreads_get(emntUpdate))
497 gmx_incons("The number of Gaussian number generators should be equal to gmx_omp_nthreads_get(emntUpdate)");
500 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
504 th = gmx_omp_get_thread_num();
507 /* Initialize on each thread to get memory allocated thread-local */
508 sd->gaussrand[th] = gmx_rng_init(seed[th]);
515 static gmx_stochd_t *init_stochd(t_inputrec *ir, int nthreads)
524 /* Initiate random number generator for langevin type dynamics,
525 * for BD, SD or velocity rescaling temperature coupling.
527 if (ir->eI == eiBD || EI_SD(ir->eI))
529 sd->ngaussrand = nthreads;
535 snew(sd->gaussrand, sd->ngaussrand);
537 /* Initialize the first random generator */
538 sd->gaussrand[0] = gmx_rng_init(ir->ld_seed);
540 if (sd->ngaussrand > 1)
542 /* Initialize the rest of the random number generators,
543 * using the first one to generate seeds.
545 init_multiple_gaussrand(sd);
548 ngtc = ir->opts.ngtc;
552 snew(sd->bd_rf, ngtc);
554 else if (EI_SD(ir->eI))
557 snew(sd->sdsig, ngtc);
560 for (n = 0; n < ngtc; n++)
562 if (ir->opts.tau_t[n] > 0)
564 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
565 sdc[n].eph = exp(sdc[n].gdt/2);
566 sdc[n].emh = exp(-sdc[n].gdt/2);
567 sdc[n].em = exp(-sdc[n].gdt);
571 /* No friction and noise on this group */
577 if (sdc[n].gdt >= 0.05)
579 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
580 - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
581 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
582 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
587 /* Seventh order expansions for small y */
588 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
589 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))));
590 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
594 fprintf(debug, "SD const tc-grp %d: b %g c %g d %g\n",
595 n, sdc[n].b, sdc[n].c, sdc[n].d);
599 else if (ETC_ANDERSEN(ir->etc))
608 snew(sd->randomize_group, ngtc);
609 snew(sd->boltzfac, ngtc);
611 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
612 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
614 for (n = 0; n < ngtc; n++)
616 reft = max(0.0, opts->ref_t[n]);
617 if ((opts->tau_t[n] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
619 sd->randomize_group[n] = TRUE;
620 sd->boltzfac[n] = BOLTZ*opts->ref_t[n];
624 sd->randomize_group[n] = FALSE;
631 void get_stochd_state(gmx_update_t upd, t_state *state)
633 /* Note that we only get the state of the first random generator,
634 * even if there are multiple. This avoids repetition.
636 gmx_rng_get_state(upd->sd->gaussrand[0], state->ld_rng, state->ld_rngi);
639 void set_stochd_state(gmx_update_t upd, t_state *state)
646 gmx_rng_set_state(sd->gaussrand[0], state->ld_rng, state->ld_rngi[0]);
648 if (sd->ngaussrand > 1)
650 /* We only end up here with SD or BD with OpenMP.
651 * Destroy and reinitialize the rest of the random number generators,
652 * using seeds generated from the first one.
653 * Although this doesn't recover the previous state,
654 * it at least avoids repetition, which is most important.
655 * Exaclty restoring states with all MPI+OpenMP setups is difficult
656 * and as the integrator is random to start with, doesn't gain us much.
658 for (i = 1; i < sd->ngaussrand; i++)
660 gmx_rng_destroy(sd->gaussrand[i]);
663 init_multiple_gaussrand(sd);
667 gmx_update_t init_update(t_inputrec *ir)
673 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
675 upd->sd = init_stochd(ir, gmx_omp_nthreads_get(emntUpdate));
680 upd->randatom = NULL;
681 upd->randatom_list = NULL;
682 upd->randatom_list_init = FALSE; /* we have not yet cleared the data structure at this point */
687 static void do_update_sd1(gmx_stochd_t *sd,
689 int start, int nrend, double dt,
690 rvec accel[], ivec nFreeze[],
691 real invmass[], unsigned short ptype[],
692 unsigned short cFREEZE[], unsigned short cACC[],
693 unsigned short cTC[],
694 rvec x[], rvec xprime[], rvec v[], rvec f[],
695 int ngtc, real tau_t[], real ref_t[])
700 int gf = 0, ga = 0, gt = 0;
707 for (n = 0; n < ngtc; n++)
710 /* The mass is encounted for later, since this differs per atom */
711 sig[n].V = sqrt(kT*(1 - sdc[n].em*sdc[n].em));
714 for (n = start; n < nrend; n++)
716 ism = sqrt(invmass[n]);
730 for (d = 0; d < DIM; d++)
732 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
734 sd_V = ism*sig[gt].V*gmx_rng_gaussian_table(gaussrand);
736 v[n][d] = v[n][d]*sdc[gt].em
737 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
740 xprime[n][d] = x[n][d] + v[n][d]*dt;
745 xprime[n][d] = x[n][d];
751 static void check_sd2_work_data_allocation(gmx_stochd_t *sd, int nrend)
753 if (nrend > sd->sd_V_nalloc)
755 sd->sd_V_nalloc = over_alloc_dd(nrend);
756 srenew(sd->sd_V, sd->sd_V_nalloc);
760 static void do_update_sd2_Tconsts(gmx_stochd_t *sd,
765 /* This is separated from the update below, because it is single threaded */
774 for (gt = 0; gt < ngtc; gt++)
776 kT = BOLTZ*ref_t[gt];
777 /* The mass is encounted for later, since this differs per atom */
778 sig[gt].V = sqrt(kT*(1-sdc[gt].em));
779 sig[gt].X = sqrt(kT*sqr(tau_t[gt])*sdc[gt].c);
780 sig[gt].Yv = sqrt(kT*sdc[gt].b/sdc[gt].c);
781 sig[gt].Yx = sqrt(kT*sqr(tau_t[gt])*sdc[gt].b/(1-sdc[gt].em));
785 static void do_update_sd2(gmx_stochd_t *sd,
788 int start, int nrend,
789 rvec accel[], ivec nFreeze[],
790 real invmass[], unsigned short ptype[],
791 unsigned short cFREEZE[], unsigned short cACC[],
792 unsigned short cTC[],
793 rvec x[], rvec xprime[], rvec v[], rvec f[],
800 /* The random part of the velocity update, generated in the first
801 * half of the update, needs to be remembered for the second half.
805 int gf = 0, ga = 0, gt = 0;
806 real vn = 0, Vmh, Xmh;
814 for (n = start; n < nrend; n++)
816 ism = sqrt(invmass[n]);
830 for (d = 0; d < DIM; d++)
836 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
842 sd_X[n][d] = ism*sig[gt].X*gmx_rng_gaussian_table(gaussrand);
844 Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
845 + ism*sig[gt].Yv*gmx_rng_gaussian_table(gaussrand);
846 sd_V[n][d] = ism*sig[gt].V*gmx_rng_gaussian_table(gaussrand);
848 v[n][d] = vn*sdc[gt].em
849 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
850 + sd_V[n][d] - sdc[gt].em*Vmh;
852 xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
857 /* Correct the velocities for the constraints.
858 * This operation introduces some inaccuracy,
859 * since the velocity is determined from differences in coordinates.
862 (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
864 Xmh = sd_V[n][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
865 + ism*sig[gt].Yx*gmx_rng_gaussian_table(gaussrand);
866 sd_X[n][d] = ism*sig[gt].X*gmx_rng_gaussian_table(gaussrand);
868 xprime[n][d] += sd_X[n][d] - Xmh;
877 xprime[n][d] = x[n][d];
884 static void do_update_bd_Tconsts(double dt, real friction_coefficient,
885 int ngtc, const real ref_t[],
888 /* This is separated from the update below, because it is single threaded */
891 if (friction_coefficient != 0)
893 for (gt = 0; gt < ngtc; gt++)
895 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]/(friction_coefficient*dt));
900 for (gt = 0; gt < ngtc; gt++)
902 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]);
907 static void do_update_bd(int start, int nrend, double dt,
909 real invmass[], unsigned short ptype[],
910 unsigned short cFREEZE[], unsigned short cTC[],
911 rvec x[], rvec xprime[], rvec v[],
912 rvec f[], real friction_coefficient,
913 real *rf, gmx_rng_t gaussrand)
915 /* note -- these appear to be full step velocities . . . */
921 if (friction_coefficient != 0)
923 invfr = 1.0/friction_coefficient;
926 for (n = start; (n < nrend); n++)
936 for (d = 0; (d < DIM); d++)
938 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
940 if (friction_coefficient != 0)
942 vn = invfr*f[n][d] + rf[gt]*gmx_rng_gaussian_table(gaussrand);
946 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
947 vn = 0.5*invmass[n]*f[n][d]*dt
948 + sqrt(0.5*invmass[n])*rf[gt]*gmx_rng_gaussian_table(gaussrand);
952 xprime[n][d] = x[n][d]+vn*dt;
957 xprime[n][d] = x[n][d];
963 static void dump_it_all(FILE gmx_unused *fp, const char gmx_unused *title,
964 int gmx_unused natoms, rvec gmx_unused x[], rvec gmx_unused xp[],
965 rvec gmx_unused v[], rvec gmx_unused f[])
970 fprintf(fp, "%s\n", title);
971 pr_rvecs(fp, 0, "x", x, natoms);
972 pr_rvecs(fp, 0, "xp", xp, natoms);
973 pr_rvecs(fp, 0, "v", v, natoms);
974 pr_rvecs(fp, 0, "f", f, natoms);
979 static void calc_ke_part_normal(rvec v[], t_grpopts *opts, t_mdatoms *md,
980 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel,
981 gmx_bool bSaveEkinOld)
984 t_grp_tcstat *tcstat = ekind->tcstat;
985 t_grp_acc *grpstat = ekind->grpstat;
988 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
989 an option, but not supported now. Additionally, if we are doing iterations.
990 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
991 bSavEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't copy over the ekinh_old.
992 If FALSE, we overrwrite it.
995 /* group velocities are calculated in update_ekindata and
996 * accumulated in acumulate_groups.
997 * Now the partial global and groups ekin.
999 for (g = 0; (g < opts->ngtc); g++)
1004 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
1008 clear_mat(tcstat[g].ekinf);
1012 clear_mat(tcstat[g].ekinh);
1016 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
1019 ekind->dekindl_old = ekind->dekindl;
1021 nthread = gmx_omp_nthreads_get(emntUpdate);
1023 #pragma omp parallel for num_threads(nthread) schedule(static)
1024 for (thread = 0; thread < nthread; thread++)
1026 int start_t, end_t, n;
1034 start_t = md->start + ((thread+0)*md->homenr)/nthread;
1035 end_t = md->start + ((thread+1)*md->homenr)/nthread;
1037 ekin_sum = ekind->ekin_work[thread];
1038 dekindl_sum = ekind->dekindl_work[thread];
1040 for (gt = 0; gt < opts->ngtc; gt++)
1042 clear_mat(ekin_sum[gt]);
1048 for (n = start_t; n < end_t; n++)
1058 hm = 0.5*md->massT[n];
1060 for (d = 0; (d < DIM); d++)
1062 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
1064 for (d = 0; (d < DIM); d++)
1066 for (m = 0; (m < DIM); m++)
1068 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1069 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
1072 if (md->nMassPerturbed && md->bPerturbed[n])
1075 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1081 for (thread = 0; thread < nthread; thread++)
1083 for (g = 0; g < opts->ngtc; g++)
1087 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
1092 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1097 ekind->dekindl += *ekind->dekindl_work[thread];
1100 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1103 static void calc_ke_part_visc(matrix box, rvec x[], rvec v[],
1104 t_grpopts *opts, t_mdatoms *md,
1105 gmx_ekindata_t *ekind,
1106 t_nrnb *nrnb, gmx_bool bEkinAveVel)
1108 int start = md->start, homenr = md->homenr;
1109 int g, d, n, m, gt = 0;
1112 t_grp_tcstat *tcstat = ekind->tcstat;
1113 t_cos_acc *cosacc = &(ekind->cosacc);
1118 for (g = 0; g < opts->ngtc; g++)
1120 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1121 clear_mat(ekind->tcstat[g].ekinh);
1123 ekind->dekindl_old = ekind->dekindl;
1125 fac = 2*M_PI/box[ZZ][ZZ];
1128 for (n = start; n < start+homenr; n++)
1134 hm = 0.5*md->massT[n];
1136 /* Note that the times of x and v differ by half a step */
1137 /* MRS -- would have to be changed for VV */
1138 cosz = cos(fac*x[n][ZZ]);
1139 /* Calculate the amplitude of the new velocity profile */
1140 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1142 copy_rvec(v[n], v_corrt);
1143 /* Subtract the profile for the kinetic energy */
1144 v_corrt[XX] -= cosz*cosacc->vcos;
1145 for (d = 0; (d < DIM); d++)
1147 for (m = 0; (m < DIM); m++)
1149 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1152 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1156 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1160 if (md->nPerturbed && md->bPerturbed[n])
1162 dekindl += 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1165 ekind->dekindl = dekindl;
1166 cosacc->mvcos = mvcos;
1168 inc_nrnb(nrnb, eNR_EKIN, homenr);
1171 void calc_ke_part(t_state *state, t_grpopts *opts, t_mdatoms *md,
1172 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
1174 if (ekind->cosacc.cos_accel == 0)
1176 calc_ke_part_normal(state->v, opts, md, ekind, nrnb, bEkinAveVel, bSaveEkinOld);
1180 calc_ke_part_visc(state->box, state->x, state->v, opts, md, ekind, nrnb, bEkinAveVel);
1184 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1186 ekinstate->ekin_n = ir->opts.ngtc;
1187 snew(ekinstate->ekinh, ekinstate->ekin_n);
1188 snew(ekinstate->ekinf, ekinstate->ekin_n);
1189 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1190 snew(ekinstate->ekinscalef_nhc, ekinstate->ekin_n);
1191 snew(ekinstate->ekinscaleh_nhc, ekinstate->ekin_n);
1192 snew(ekinstate->vscale_nhc, ekinstate->ekin_n);
1193 ekinstate->dekindl = 0;
1194 ekinstate->mvcos = 0;
1197 void update_ekinstate(ekinstate_t *ekinstate, gmx_ekindata_t *ekind)
1201 for (i = 0; i < ekinstate->ekin_n; i++)
1203 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1204 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1205 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1206 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1207 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1208 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1211 copy_mat(ekind->ekin, ekinstate->ekin_total);
1212 ekinstate->dekindl = ekind->dekindl;
1213 ekinstate->mvcos = ekind->cosacc.mvcos;
1217 void restore_ekinstate_from_state(t_commrec *cr,
1218 gmx_ekindata_t *ekind, ekinstate_t *ekinstate)
1224 for (i = 0; i < ekinstate->ekin_n; i++)
1226 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1227 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1228 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1229 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1230 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1231 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1234 copy_mat(ekinstate->ekin_total, ekind->ekin);
1236 ekind->dekindl = ekinstate->dekindl;
1237 ekind->cosacc.mvcos = ekinstate->mvcos;
1238 n = ekinstate->ekin_n;
1243 gmx_bcast(sizeof(n), &n, cr);
1244 for (i = 0; i < n; i++)
1246 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1247 ekind->tcstat[i].ekinh[0], cr);
1248 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1249 ekind->tcstat[i].ekinf[0], cr);
1250 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1251 ekind->tcstat[i].ekinh_old[0], cr);
1253 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1254 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1255 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1256 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1257 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1258 &(ekind->tcstat[i].vscale_nhc), cr);
1260 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1261 ekind->ekin[0], cr);
1263 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1264 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1268 void set_deform_reference_box(gmx_update_t upd, gmx_int64_t step, matrix box)
1270 upd->deformref_step = step;
1271 copy_mat(box, upd->deformref_box);
1274 static void deform(gmx_update_t upd,
1275 int start, int homenr, rvec x[], matrix box, matrix *scale_tot,
1276 const t_inputrec *ir, gmx_int64_t step)
1278 matrix bnew, invbox, mu;
1282 elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1283 copy_mat(box, bnew);
1284 for (i = 0; i < DIM; i++)
1286 for (j = 0; j < DIM; j++)
1288 if (ir->deform[i][j] != 0)
1291 upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1295 /* We correct the off-diagonal elements,
1296 * which can grow indefinitely during shearing,
1297 * so the shifts do not get messed up.
1299 for (i = 1; i < DIM; i++)
1301 for (j = i-1; j >= 0; j--)
1303 while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1305 rvec_dec(bnew[i], bnew[j]);
1307 while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1309 rvec_inc(bnew[i], bnew[j]);
1313 m_inv_ur0(box, invbox);
1314 copy_mat(bnew, box);
1315 mmul_ur0(box, invbox, mu);
1317 for (i = start; i < start+homenr; i++)
1319 x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1320 x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1321 x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1325 /* The transposes of the scaling matrices are stored,
1326 * so we need to do matrix multiplication in the inverse order.
1328 mmul_ur0(*scale_tot, mu, *scale_tot);
1332 static void combine_forces(int nstcalclr,
1333 gmx_constr_t constr,
1334 t_inputrec *ir, t_mdatoms *md, t_idef *idef,
1337 t_state *state, gmx_bool bMolPBC,
1338 int start, int nrend,
1339 rvec f[], rvec f_lr[],
1344 /* f contains the short-range forces + the long range forces
1345 * which are stored separately in f_lr.
1348 if (constr != NULL && !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1350 /* We need to constrain the LR forces separately,
1351 * because due to the different pre-factor for the SR and LR
1352 * forces in the update algorithm, we can not determine
1353 * the constraint force for the coordinate constraining.
1354 * Constrain only the additional LR part of the force.
1356 /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1357 constrain(NULL, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
1358 state->x, f_lr, f_lr, bMolPBC, state->box, state->lambda[efptBONDED], NULL,
1359 NULL, NULL, nrnb, econqForce, ir->epc == epcMTTK, state->veta, state->veta);
1362 /* Add nstcalclr-1 times the LR force to the sum of both forces
1363 * and store the result in forces_lr.
1365 nm1 = nstcalclr - 1;
1366 for (i = start; i < nrend; i++)
1368 for (d = 0; d < DIM; d++)
1370 f_lr[i][d] = f[i][d] + nm1*f_lr[i][d];
1375 void update_tcouple(gmx_int64_t step,
1376 t_inputrec *inputrec,
1378 gmx_ekindata_t *ekind,
1384 gmx_bool bTCouple = FALSE;
1386 int i, start, end, homenr, offset;
1388 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1389 if (inputrec->etc != etcNO &&
1390 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec)))
1392 /* We should only couple after a step where energies were determined (for leapfrog versions)
1393 or the step energies are determined, for velocity verlet versions */
1395 if (EI_VV(inputrec->eI))
1403 bTCouple = (inputrec->nsttcouple == 1 ||
1404 do_per_step(step+inputrec->nsttcouple-offset,
1405 inputrec->nsttcouple));
1410 dttc = inputrec->nsttcouple*inputrec->delta_t;
1412 switch (inputrec->etc)
1417 berendsen_tcoupl(inputrec, ekind, dttc);
1420 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1421 state->nosehoover_xi, state->nosehoover_vxi, MassQ);
1424 vrescale_tcoupl(inputrec, ekind, dttc,
1425 state->therm_integral, upd->sd->gaussrand[0]);
1428 /* rescale in place here */
1429 if (EI_VV(inputrec->eI))
1431 rescale_velocities(ekind, md, md->start, md->start+md->homenr, state->v);
1436 /* Set the T scaling lambda to 1 to have no scaling */
1437 for (i = 0; (i < inputrec->opts.ngtc); i++)
1439 ekind->tcstat[i].lambda = 1.0;
1444 void update_pcouple(FILE *fplog,
1446 t_inputrec *inputrec,
1452 gmx_bool bPCouple = FALSE;
1456 /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1457 if (inputrec->epc != epcNO && (!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))))
1459 /* We should only couple after a step where energies were determined */
1460 bPCouple = (inputrec->nstpcouple == 1 ||
1461 do_per_step(step+inputrec->nstpcouple-1,
1462 inputrec->nstpcouple));
1465 clear_mat(pcoupl_mu);
1466 for (i = 0; i < DIM; i++)
1468 pcoupl_mu[i][i] = 1.0;
1475 dtpc = inputrec->nstpcouple*inputrec->delta_t;
1477 switch (inputrec->epc)
1479 /* We can always pcoupl, even if we did not sum the energies
1480 * the previous step, since state->pres_prev is only updated
1481 * when the energies have been summed.
1485 case (epcBERENDSEN):
1488 berendsen_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev, state->box,
1492 case (epcPARRINELLORAHMAN):
1493 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1494 state->box, state->box_rel, state->boxv,
1495 M, pcoupl_mu, bInitStep);
1503 static rvec *get_xprime(const t_state *state, gmx_update_t upd)
1505 if (state->nalloc > upd->xp_nalloc)
1507 upd->xp_nalloc = state->nalloc;
1508 srenew(upd->xp, upd->xp_nalloc);
1514 void update_constraints(FILE *fplog,
1516 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1517 t_inputrec *inputrec, /* input record and box stuff */
1518 gmx_ekindata_t *ekind,
1523 rvec force[], /* forces on home particles */
1528 gmx_wallcycle_t wcycle,
1530 gmx_constr_t constr,
1531 gmx_bool bFirstHalf,
1535 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1538 int start, homenr, nrend, i, n, m, g, d;
1540 rvec *vbuf, *xprime = NULL;
1547 if (bFirstHalf && !EI_VV(inputrec->eI))
1552 /* for now, SD update is here -- though it really seems like it
1553 should be reformulated as a velocity verlet method, since it has two parts */
1556 homenr = md->homenr;
1557 nrend = start+homenr;
1559 dt = inputrec->delta_t;
1564 * APPLY CONSTRAINTS:
1567 * When doing PR pressure coupling we have to constrain the
1568 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1569 * it is enough to do this once though, since the relative velocities
1570 * after this will be normal to the bond vector
1575 /* clear out constraints before applying */
1576 clear_mat(vir_part);
1578 xprime = get_xprime(state, upd);
1580 bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1581 bLog = (do_per_step(step, inputrec->nstlog) || bLastStep || (step < 0));
1582 bEner = (do_per_step(step, inputrec->nstenergy) || bLastStep);
1583 /* Constrain the coordinates xprime */
1584 wallcycle_start(wcycle, ewcCONSTR);
1585 if (EI_VV(inputrec->eI) && bFirstHalf)
1587 constrain(NULL, bLog, bEner, constr, idef,
1588 inputrec, ekind, cr, step, 1, md,
1589 state->x, state->v, state->v,
1590 bMolPBC, state->box,
1591 state->lambda[efptBONDED], dvdlambda,
1592 NULL, bCalcVir ? &vir_con : NULL, nrnb, econqVeloc,
1593 inputrec->epc == epcMTTK, state->veta, vetanew);
1597 constrain(NULL, bLog, bEner, constr, idef,
1598 inputrec, ekind, cr, step, 1, md,
1599 state->x, xprime, NULL,
1600 bMolPBC, state->box,
1601 state->lambda[efptBONDED], dvdlambda,
1602 state->v, bCalcVir ? &vir_con : NULL, nrnb, econqCoord,
1603 inputrec->epc == epcMTTK, state->veta, state->veta);
1605 wallcycle_stop(wcycle, ewcCONSTR);
1609 dump_it_all(fplog, "After Shake",
1610 state->natoms, state->x, xprime, state->v, force);
1614 if (inputrec->eI == eiSD2)
1616 /* A correction factor eph is needed for the SD constraint force */
1617 /* Here we can, unfortunately, not have proper corrections
1618 * for different friction constants, so we use the first one.
1620 for (i = 0; i < DIM; i++)
1622 for (m = 0; m < DIM; m++)
1624 vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1630 m_add(vir_part, vir_con, vir_part);
1634 pr_rvecs(debug, 0, "constraint virial", vir_part, DIM);
1640 if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1642 xprime = get_xprime(state, upd);
1644 nth = gmx_omp_nthreads_get(emntUpdate);
1646 #pragma omp parallel for num_threads(nth) schedule(static)
1647 for (th = 0; th < nth; th++)
1649 int start_th, end_th;
1651 start_th = start + ((nrend-start)* th )/nth;
1652 end_th = start + ((nrend-start)*(th+1))/nth;
1654 /* The second part of the SD integration */
1655 do_update_sd2(upd->sd, upd->sd->gaussrand[th],
1656 FALSE, start_th, end_th,
1657 inputrec->opts.acc, inputrec->opts.nFreeze,
1658 md->invmass, md->ptype,
1659 md->cFREEZE, md->cACC, md->cTC,
1660 state->x, xprime, state->v, force, state->sd_X,
1661 inputrec->opts.tau_t,
1664 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1668 /* Constrain the coordinates xprime */
1669 wallcycle_start(wcycle, ewcCONSTR);
1670 constrain(NULL, bLog, bEner, constr, idef,
1671 inputrec, NULL, cr, step, 1, md,
1672 state->x, xprime, NULL,
1673 bMolPBC, state->box,
1674 state->lambda[efptBONDED], dvdlambda,
1675 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
1676 wallcycle_stop(wcycle, ewcCONSTR);
1680 /* We must always unshift after updating coordinates; if we did not shake
1681 x was shifted in do_force */
1683 if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1685 if (graph && (graph->nnodes > 0))
1687 unshift_x(graph, state->box, state->x, upd->xp);
1688 if (TRICLINIC(state->box))
1690 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1694 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1699 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntUpdate)) schedule(static)
1700 for (i = start; i < nrend; i++)
1702 copy_rvec(upd->xp[i], state->x[i]);
1706 dump_it_all(fplog, "After unshift",
1707 state->natoms, state->x, upd->xp, state->v, force);
1709 /* ############# END the update of velocities and positions ######### */
1712 void update_box(FILE *fplog,
1714 t_inputrec *inputrec, /* input record and box stuff */
1717 rvec force[], /* forces on home particles */
1723 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE;
1726 int start, homenr, nrend, i, n, m, g;
1730 homenr = md->homenr;
1731 nrend = start+homenr;
1734 (inputrec->etc == etcNOSEHOOVER) ||
1735 (inputrec->epc == epcPARRINELLORAHMAN) ||
1736 (inputrec->epc == epcMTTK);
1738 dt = inputrec->delta_t;
1742 /* now update boxes */
1743 switch (inputrec->epc)
1747 case (epcBERENDSEN):
1748 berendsen_pscale(inputrec, pcoupl_mu, state->box, state->box_rel,
1749 start, homenr, state->x, md->cFREEZE, nrnb);
1751 case (epcPARRINELLORAHMAN):
1752 /* The box velocities were updated in do_pr_pcoupl in the update
1753 * iteration, but we dont change the box vectors until we get here
1754 * since we need to be able to shift/unshift above.
1756 for (i = 0; i < DIM; i++)
1758 for (m = 0; m <= i; m++)
1760 state->box[i][m] += dt*state->boxv[i][m];
1763 preserve_box_shape(inputrec, state->box_rel, state->box);
1765 /* Scale the coordinates */
1766 for (n = start; (n < start+homenr); n++)
1768 tmvmul_ur0(pcoupl_mu, state->x[n], state->x[n]);
1772 switch (inputrec->epct)
1774 case (epctISOTROPIC):
1775 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1776 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1777 Side length scales as exp(veta*dt) */
1779 msmul(state->box, exp(state->veta*dt), state->box);
1781 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1782 o If we assume isotropic scaling, and box length scaling
1783 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1784 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1785 determinant of B is L^DIM det(M), and the determinant
1786 of dB/dt is (dL/dT)^DIM det (M). veta will be
1787 (det(dB/dT)/det(B))^(1/3). Then since M =
1788 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1790 msmul(state->box, state->veta, state->boxv);
1800 if ((!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))) && scale_tot)
1802 /* The transposes of the scaling matrices are stored,
1803 * therefore we need to reverse the order in the multiplication.
1805 mmul_ur0(*scale_tot, pcoupl_mu, *scale_tot);
1808 if (DEFORM(*inputrec))
1810 deform(upd, start, homenr, state->x, state->box, scale_tot, inputrec, step);
1813 dump_it_all(fplog, "After update",
1814 state->natoms, state->x, upd->xp, state->v, force);
1817 void update_coords(FILE *fplog,
1819 t_inputrec *inputrec, /* input record and box stuff */
1823 rvec *f, /* forces on home particles */
1827 gmx_ekindata_t *ekind,
1832 t_commrec *cr, /* these shouldn't be here -- need to think about it */
1834 gmx_constr_t constr,
1837 gmx_bool bNH, bPR, bLastStep, bLog = FALSE, bEner = FALSE;
1839 real *imass, *imassin;
1842 int start, homenr, nrend, i, j, d, n, m, g;
1843 int blen0, blen1, iatom, jatom, nshake, nsettle, nconstr, nexpand;
1846 rvec *vcom, *xcom, *vall, *xall, *xin, *vin, *forcein, *fall, *xpall, *xprimein, *xprime;
1849 /* Running the velocity half does nothing except for velocity verlet */
1850 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1851 !EI_VV(inputrec->eI))
1853 gmx_incons("update_coords called for velocity without VV integrator");
1857 homenr = md->homenr;
1858 nrend = start+homenr;
1860 xprime = get_xprime(state, upd);
1862 dt = inputrec->delta_t;
1865 /* We need to update the NMR restraint history when time averaging is used */
1866 if (state->flags & (1<<estDISRE_RM3TAV))
1868 update_disres_history(fcd, &state->hist);
1870 if (state->flags & (1<<estORIRE_DTAV))
1872 update_orires_history(fcd, &state->hist);
1876 bNH = inputrec->etc == etcNOSEHOOVER;
1877 bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1879 if (bDoLR && inputrec->nstcalclr > 1 && !EI_VV(inputrec->eI)) /* get this working with VV? */
1881 /* Store the total force + nstcalclr-1 times the LR force
1882 * in forces_lr, so it can be used in a normal update algorithm
1883 * to produce twin time stepping.
1885 /* is this correct in the new construction? MRS */
1886 combine_forces(inputrec->nstcalclr, constr, inputrec, md, idef, cr,
1887 step, state, bMolPBC,
1888 start, nrend, f, f_lr, nrnb);
1896 /* ############# START The update of velocities and positions ######### */
1898 dump_it_all(fplog, "Before update",
1899 state->natoms, state->x, xprime, state->v, force);
1901 if (inputrec->eI == eiSD2)
1903 check_sd2_work_data_allocation(upd->sd, nrend);
1905 do_update_sd2_Tconsts(upd->sd,
1906 inputrec->opts.ngtc,
1907 inputrec->opts.tau_t,
1908 inputrec->opts.ref_t);
1910 if (inputrec->eI == eiBD)
1912 do_update_bd_Tconsts(dt, inputrec->bd_fric,
1913 inputrec->opts.ngtc, inputrec->opts.ref_t,
1917 nth = gmx_omp_nthreads_get(emntUpdate);
1919 #pragma omp parallel for num_threads(nth) schedule(static) private(alpha)
1920 for (th = 0; th < nth; th++)
1922 int start_th, end_th;
1924 start_th = start + ((nrend-start)* th )/nth;
1925 end_th = start + ((nrend-start)*(th+1))/nth;
1927 switch (inputrec->eI)
1930 if (ekind->cosacc.cos_accel == 0)
1932 do_update_md(start_th, end_th, dt,
1933 ekind->tcstat, state->nosehoover_vxi,
1934 ekind->bNEMD, ekind->grpstat, inputrec->opts.acc,
1935 inputrec->opts.nFreeze,
1936 md->invmass, md->ptype,
1937 md->cFREEZE, md->cACC, md->cTC,
1938 state->x, xprime, state->v, force, M,
1943 do_update_visc(start_th, end_th, dt,
1944 ekind->tcstat, state->nosehoover_vxi,
1945 md->invmass, md->ptype,
1946 md->cTC, state->x, xprime, state->v, force, M,
1948 ekind->cosacc.cos_accel,
1954 do_update_sd1(upd->sd, upd->sd->gaussrand[th],
1955 start_th, end_th, dt,
1956 inputrec->opts.acc, inputrec->opts.nFreeze,
1957 md->invmass, md->ptype,
1958 md->cFREEZE, md->cACC, md->cTC,
1959 state->x, xprime, state->v, force,
1960 inputrec->opts.ngtc, inputrec->opts.tau_t, inputrec->opts.ref_t);
1963 /* The SD update is done in 2 parts, because an extra constraint step
1966 do_update_sd2(upd->sd, upd->sd->gaussrand[th],
1967 bInitStep, start_th, end_th,
1968 inputrec->opts.acc, inputrec->opts.nFreeze,
1969 md->invmass, md->ptype,
1970 md->cFREEZE, md->cACC, md->cTC,
1971 state->x, xprime, state->v, force, state->sd_X,
1972 inputrec->opts.tau_t,
1976 do_update_bd(start_th, end_th, dt,
1977 inputrec->opts.nFreeze, md->invmass, md->ptype,
1978 md->cFREEZE, md->cTC,
1979 state->x, xprime, state->v, force,
1981 upd->sd->bd_rf, upd->sd->gaussrand[th]);
1985 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
1990 do_update_vv_vel(start_th, end_th, dt,
1991 inputrec->opts.acc, inputrec->opts.nFreeze,
1992 md->invmass, md->ptype,
1993 md->cFREEZE, md->cACC,
1995 (bNH || bPR), state->veta, alpha);
1998 do_update_vv_pos(start_th, end_th, dt,
1999 inputrec->opts.nFreeze,
2000 md->ptype, md->cFREEZE,
2001 state->x, xprime, state->v,
2002 (bNH || bPR), state->veta);
2007 gmx_fatal(FARGS, "Don't know how to update coordinates");
2015 void correct_ekin(FILE *log, int start, int end, rvec v[], rvec vcm, real mass[],
2016 real tmass, tensor ekin)
2019 * This is a debugging routine. It should not be called for production code
2021 * The kinetic energy should calculated according to:
2022 * Ekin = 1/2 m (v-vcm)^2
2023 * However the correction is not always applied, since vcm may not be
2024 * known in time and we compute
2025 * Ekin' = 1/2 m v^2 instead
2026 * This can be corrected afterwards by computing
2027 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
2029 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
2036 /* Local particles */
2039 /* Processor dependent part. */
2041 for (i = start; (i < end); i++)
2045 for (j = 0; (j < DIM); j++)
2051 svmul(1/tmass, vcm, vcm);
2052 svmul(0.5, vcm, hvcm);
2054 for (j = 0; (j < DIM); j++)
2056 for (k = 0; (k < DIM); k++)
2058 dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
2061 pr_rvecs(log, 0, "dekin", dekin, DIM);
2062 pr_rvecs(log, 0, " ekin", ekin, DIM);
2063 fprintf(log, "dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
2064 trace(dekin), trace(ekin), vcm[XX], vcm[YY], vcm[ZZ]);
2065 fprintf(log, "mv = (%8.4f %8.4f %8.4f)\n",
2066 mv[XX], mv[YY], mv[ZZ]);
2069 extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_int64_t step, t_mdatoms *md, t_state *state, gmx_update_t upd, t_idef *idef, gmx_constr_t constr)
2073 real rate = (ir->delta_t)/ir->opts.tau_t[0];
2074 /* proceed with andersen if 1) it's fixed probability per
2075 particle andersen or 2) it's massive andersen and it's tau_t/dt */
2076 if ((ir->etc == etcANDERSEN) || do_per_step(step, (int)(1.0/rate)))
2078 srenew(upd->randatom, state->nalloc);
2079 srenew(upd->randatom_list, state->nalloc);
2080 if (upd->randatom_list_init == FALSE)
2082 for (i = 0; i < state->nalloc; i++)
2084 upd->randatom[i] = FALSE;
2085 upd->randatom_list[i] = 0;
2087 upd->randatom_list_init = TRUE;
2089 andersen_tcoupl(ir, md, state, upd->sd->gaussrand[0], rate,
2090 (ir->etc == etcANDERSEN) ? idef : NULL,
2091 constr ? get_nblocks(constr) : 0,
2092 constr ? get_sblock(constr) : NULL,
2093 upd->randatom, upd->randatom_list,
2094 upd->sd->randomize_group, upd->sd->boltzfac);