1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
44 #include "types/commrec.h"
55 #include "gmx_random.h"
68 #include "gmx_wallcycle.h"
69 #include "gmx_omp_nthreads.h"
72 /*For debugging, start at v(-dt/2) for velolcity verlet -- uncomment next line */
73 /*#define STARTFROMDT2*/
93 /* The random state for ngaussrand threads.
94 * Normal thermostats need just 1 random number generator,
95 * but SD and BD with OpenMP parallelization need 1 for each thread.
103 gmx_sd_sigma_t *sdsig;
106 /* andersen temperature control stuff */
107 gmx_bool *randomize_group;
111 typedef struct gmx_update
114 /* xprime for constraint algorithms */
118 /* variable size arrays for andersen */
121 gmx_bool randatom_list_init;
123 /* Variables for the deform algorithm */
124 gmx_large_int_t deformref_step;
125 matrix deformref_box;
129 static void do_update_md(int start, int nrend, double dt,
130 t_grp_tcstat *tcstat,
132 gmx_bool bNEMD, t_grp_acc *gstat, rvec accel[],
135 unsigned short ptype[], unsigned short cFREEZE[],
136 unsigned short cACC[], unsigned short cTC[],
137 rvec x[], rvec xprime[], rvec v[],
139 gmx_bool bNH, gmx_bool bPR)
142 int gf = 0, ga = 0, gt = 0;
144 real vn, vv, va, vb, vnrel;
150 /* Update with coupling to extended ensembles, used for
151 * Nose-Hoover and Parrinello-Rahman coupling
152 * Nose-Hoover uses the reversible leap-frog integrator from
153 * Holian et al. Phys Rev E 52(3) : 2338, 1995
155 for (n = start; n < nrend; n++)
170 lg = tcstat[gt].lambda;
175 rvec_sub(v[n], gstat[ga].u, vrel);
177 for (d = 0; d < DIM; d++)
179 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
181 vnrel = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
182 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
183 /* do not scale the mean velocities u */
184 vn = gstat[ga].u[d] + accel[ga][d]*dt + vnrel;
186 xprime[n][d] = x[n][d]+vn*dt;
191 xprime[n][d] = x[n][d];
196 else if (cFREEZE != NULL ||
197 nFreeze[0][XX] || nFreeze[0][YY] || nFreeze[0][ZZ] ||
200 /* Update with Berendsen/v-rescale coupling and freeze or NEMD */
201 for (n = start; n < nrend; n++)
203 w_dt = invmass[n]*dt;
216 lg = tcstat[gt].lambda;
218 for (d = 0; d < DIM; d++)
221 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
223 vv = lg*vn + f[n][d]*w_dt;
225 /* do not scale the mean velocities u */
227 va = vv + accel[ga][d]*dt;
228 vb = va + (1.0-lg)*u;
230 xprime[n][d] = x[n][d]+vb*dt;
235 xprime[n][d] = x[n][d];
242 /* Plain update with Berendsen/v-rescale coupling */
243 for (n = start; n < nrend; n++)
245 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
247 w_dt = invmass[n]*dt;
252 lg = tcstat[gt].lambda;
254 for (d = 0; d < DIM; d++)
256 vn = lg*v[n][d] + f[n][d]*w_dt;
258 xprime[n][d] = x[n][d] + vn*dt;
263 for (d = 0; d < DIM; d++)
266 xprime[n][d] = x[n][d];
273 static void do_update_vv_vel(int start, int nrend, double dt,
274 t_grp_tcstat *tcstat, t_grp_acc *gstat,
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,
325 t_grp_tcstat *tcstat, t_grp_acc *gstat,
326 rvec accel[], ivec nFreeze[], real invmass[],
327 unsigned short ptype[], unsigned short cFREEZE[],
328 rvec x[], rvec xprime[], rvec v[],
329 rvec f[], gmx_bool bExtended, real veta, real alpha)
336 /* Would it make more sense if Parrinello-Rahman was put here? */
341 mr2 = series_sinhx(g);
349 for (n = start; n < nrend; n++)
357 for (d = 0; d < DIM; d++)
359 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
361 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
365 xprime[n][d] = x[n][d];
369 } /* do_update_vv_pos */
371 static void do_update_visc(int start, int nrend, double dt,
372 t_grp_tcstat *tcstat,
375 unsigned short ptype[], unsigned short cTC[],
376 rvec x[], rvec xprime[], rvec v[],
377 rvec f[], matrix M, matrix box, real
378 cos_accel, real vcos,
379 gmx_bool bNH, gmx_bool bPR)
384 real lg, vxi = 0, vv;
389 fac = 2*M_PI/(box[ZZ][ZZ]);
393 /* Update with coupling to extended ensembles, used for
394 * Nose-Hoover and Parrinello-Rahman coupling
396 for (n = start; n < nrend; n++)
403 lg = tcstat[gt].lambda;
404 cosz = cos(fac*x[n][ZZ]);
406 copy_rvec(v[n], vrel);
414 for (d = 0; d < DIM; d++)
418 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
420 vn = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
421 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
424 vn += vc + dt*cosz*cos_accel;
427 xprime[n][d] = x[n][d]+vn*dt;
431 xprime[n][d] = x[n][d];
438 /* Classic version of update, used with berendsen coupling */
439 for (n = start; n < nrend; n++)
441 w_dt = invmass[n]*dt;
446 lg = tcstat[gt].lambda;
447 cosz = cos(fac*x[n][ZZ]);
449 for (d = 0; d < DIM; d++)
453 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
458 /* Do not scale the cosine velocity profile */
459 vv = vc + lg*(vn - vc + f[n][d]*w_dt);
460 /* Add the cosine accelaration profile */
461 vv += dt*cosz*cos_accel;
465 vv = lg*(vn + f[n][d]*w_dt);
468 xprime[n][d] = x[n][d]+vv*dt;
473 xprime[n][d] = x[n][d];
480 /* Allocates and initializes sd->gaussrand[i] for i=1, i<sd->ngaussrand,
481 * Using seeds generated from sd->gaussrand[0].
483 static void init_multiple_gaussrand(gmx_stochd_t *sd)
488 ngr = sd->ngaussrand;
491 for (i = 1; i < ngr; i++)
493 seed[i] = gmx_rng_uniform_uint32(sd->gaussrand[0]);
496 if (ngr != gmx_omp_nthreads_get(emntUpdate))
498 gmx_incons("The number of Gaussian number generators should be equal to gmx_omp_nthreads_get(emntUpdate)");
501 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
505 th = gmx_omp_get_thread_num();
508 /* Initialize on each thread to get memory allocated thread-local */
509 sd->gaussrand[th] = gmx_rng_init(seed[th]);
516 static gmx_stochd_t *init_stochd(FILE *fplog, t_inputrec *ir, int nthreads)
525 /* Initiate random number generator for langevin type dynamics,
526 * for BD, SD or velocity rescaling temperature coupling.
528 if (ir->eI == eiBD || EI_SD(ir->eI))
530 sd->ngaussrand = nthreads;
536 snew(sd->gaussrand, sd->ngaussrand);
538 /* Initialize the first random generator */
539 sd->gaussrand[0] = gmx_rng_init(ir->ld_seed);
541 if (sd->ngaussrand > 1)
543 /* Initialize the rest of the random number generators,
544 * using the first one to generate seeds.
546 init_multiple_gaussrand(sd);
549 ngtc = ir->opts.ngtc;
553 snew(sd->bd_rf, ngtc);
555 else if (EI_SD(ir->eI))
558 snew(sd->sdsig, ngtc);
561 for (n = 0; n < ngtc; n++)
563 if (ir->opts.tau_t[n] > 0)
565 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
566 sdc[n].eph = exp(sdc[n].gdt/2);
567 sdc[n].emh = exp(-sdc[n].gdt/2);
568 sdc[n].em = exp(-sdc[n].gdt);
572 /* No friction and noise on this group */
578 if (sdc[n].gdt >= 0.05)
580 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
581 - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
582 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
583 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
588 /* Seventh order expansions for small y */
589 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
590 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))));
591 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
595 fprintf(debug, "SD const tc-grp %d: b %g c %g d %g\n",
596 n, sdc[n].b, sdc[n].c, sdc[n].d);
600 else if (ETC_ANDERSEN(ir->etc))
609 snew(sd->randomize_group, ngtc);
610 snew(sd->boltzfac, ngtc);
612 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
613 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
615 for (n = 0; n < ngtc; n++)
617 reft = max(0.0, opts->ref_t[n]);
618 if ((opts->tau_t[n] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
620 sd->randomize_group[n] = TRUE;
621 sd->boltzfac[n] = BOLTZ*opts->ref_t[n];
625 sd->randomize_group[n] = FALSE;
632 void get_stochd_state(gmx_update_t upd, t_state *state)
634 /* Note that we only get the state of the first random generator,
635 * even if there are multiple. This avoids repetition.
637 gmx_rng_get_state(upd->sd->gaussrand[0], state->ld_rng, state->ld_rngi);
640 void set_stochd_state(gmx_update_t upd, t_state *state)
647 gmx_rng_set_state(sd->gaussrand[0], state->ld_rng, state->ld_rngi[0]);
649 if (sd->ngaussrand > 1)
651 /* We only end up here with SD or BD with OpenMP.
652 * Destroy and reinitialize the rest of the random number generators,
653 * using seeds generated from the first one.
654 * Although this doesn't recover the previous state,
655 * it at least avoids repetition, which is most important.
656 * Exaclty restoring states with all MPI+OpenMP setups is difficult
657 * and as the integrator is random to start with, doesn't gain us much.
659 for (i = 1; i < sd->ngaussrand; i++)
661 gmx_rng_destroy(sd->gaussrand[i]);
664 init_multiple_gaussrand(sd);
668 gmx_update_t init_update(FILE *fplog, t_inputrec *ir)
674 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
676 upd->sd = init_stochd(fplog, ir, gmx_omp_nthreads_get(emntUpdate));
681 upd->randatom = NULL;
682 upd->randatom_list = NULL;
683 upd->randatom_list_init = FALSE; /* we have not yet cleared the data structure at this point */
688 static void do_update_sd1(gmx_stochd_t *sd,
690 int start, int nrend, double dt,
691 rvec accel[], ivec nFreeze[],
692 real invmass[], unsigned short ptype[],
693 unsigned short cFREEZE[], unsigned short cACC[],
694 unsigned short cTC[],
695 rvec x[], rvec xprime[], rvec v[], rvec f[],
697 int ngtc, real tau_t[], real ref_t[])
702 int gf = 0, ga = 0, gt = 0;
709 for (n = 0; n < ngtc; n++)
712 /* The mass is encounted for later, since this differs per atom */
713 sig[n].V = sqrt(kT*(1 - sdc[n].em*sdc[n].em));
716 for (n = start; n < nrend; n++)
718 ism = sqrt(invmass[n]);
732 for (d = 0; d < DIM; d++)
734 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
736 sd_V = ism*sig[gt].V*gmx_rng_gaussian_table(gaussrand);
738 v[n][d] = v[n][d]*sdc[gt].em
739 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
742 xprime[n][d] = x[n][d] + v[n][d]*dt;
747 xprime[n][d] = x[n][d];
753 static void check_sd2_work_data_allocation(gmx_stochd_t *sd, int nrend)
755 if (nrend > sd->sd_V_nalloc)
757 sd->sd_V_nalloc = over_alloc_dd(nrend);
758 srenew(sd->sd_V, sd->sd_V_nalloc);
762 static void do_update_sd2_Tconsts(gmx_stochd_t *sd,
767 /* This is separated from the update below, because it is single threaded */
776 for (gt = 0; gt < ngtc; gt++)
778 kT = BOLTZ*ref_t[gt];
779 /* The mass is encounted for later, since this differs per atom */
780 sig[gt].V = sqrt(kT*(1-sdc[gt].em));
781 sig[gt].X = sqrt(kT*sqr(tau_t[gt])*sdc[gt].c);
782 sig[gt].Yv = sqrt(kT*sdc[gt].b/sdc[gt].c);
783 sig[gt].Yx = sqrt(kT*sqr(tau_t[gt])*sdc[gt].b/(1-sdc[gt].em));
787 static void do_update_sd2(gmx_stochd_t *sd,
790 int start, int nrend,
791 rvec accel[], ivec nFreeze[],
792 real invmass[], unsigned short ptype[],
793 unsigned short cFREEZE[], unsigned short cACC[],
794 unsigned short cTC[],
795 rvec x[], rvec xprime[], rvec v[], rvec f[],
802 /* The random part of the velocity update, generated in the first
803 * half of the update, needs to be remembered for the second half.
807 int gf = 0, ga = 0, gt = 0;
808 real vn = 0, Vmh, Xmh;
816 for (n = start; n < nrend; n++)
818 ism = sqrt(invmass[n]);
832 for (d = 0; d < DIM; d++)
838 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
844 sd_X[n][d] = ism*sig[gt].X*gmx_rng_gaussian_table(gaussrand);
846 Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
847 + ism*sig[gt].Yv*gmx_rng_gaussian_table(gaussrand);
848 sd_V[n][d] = ism*sig[gt].V*gmx_rng_gaussian_table(gaussrand);
850 v[n][d] = vn*sdc[gt].em
851 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
852 + sd_V[n][d] - sdc[gt].em*Vmh;
854 xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
859 /* Correct the velocities for the constraints.
860 * This operation introduces some inaccuracy,
861 * since the velocity is determined from differences in coordinates.
864 (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
866 Xmh = sd_V[n][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
867 + ism*sig[gt].Yx*gmx_rng_gaussian_table(gaussrand);
868 sd_X[n][d] = ism*sig[gt].X*gmx_rng_gaussian_table(gaussrand);
870 xprime[n][d] += sd_X[n][d] - Xmh;
879 xprime[n][d] = x[n][d];
886 static void do_update_bd_Tconsts(double dt, real friction_coefficient,
887 int ngtc, const real ref_t[],
890 /* This is separated from the update below, because it is single threaded */
893 if (friction_coefficient != 0)
895 for (gt = 0; gt < ngtc; gt++)
897 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]/(friction_coefficient*dt));
902 for (gt = 0; gt < ngtc; gt++)
904 rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]);
909 static void do_update_bd(int start, int nrend, double dt,
911 real invmass[], unsigned short ptype[],
912 unsigned short cFREEZE[], unsigned short cTC[],
913 rvec x[], rvec xprime[], rvec v[],
914 rvec f[], real friction_coefficient,
915 real *rf, gmx_rng_t gaussrand)
917 /* note -- these appear to be full step velocities . . . */
923 if (friction_coefficient != 0)
925 invfr = 1.0/friction_coefficient;
928 for (n = start; (n < nrend); n++)
938 for (d = 0; (d < DIM); d++)
940 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
942 if (friction_coefficient != 0)
944 vn = invfr*f[n][d] + rf[gt]*gmx_rng_gaussian_table(gaussrand);
948 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
949 vn = 0.5*invmass[n]*f[n][d]*dt
950 + sqrt(0.5*invmass[n])*rf[gt]*gmx_rng_gaussian_table(gaussrand);
954 xprime[n][d] = x[n][d]+vn*dt;
959 xprime[n][d] = x[n][d];
965 static void dump_it_all(FILE *fp, const char *title,
966 int natoms, rvec x[], rvec xp[], rvec v[], rvec f[])
971 fprintf(fp, "%s\n", title);
972 pr_rvecs(fp, 0, "x", x, natoms);
973 pr_rvecs(fp, 0, "xp", xp, natoms);
974 pr_rvecs(fp, 0, "v", v, natoms);
975 pr_rvecs(fp, 0, "f", f, natoms);
980 static void calc_ke_part_normal(rvec v[], t_grpopts *opts, t_mdatoms *md,
981 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel,
982 gmx_bool bSaveEkinOld)
985 t_grp_tcstat *tcstat = ekind->tcstat;
986 t_grp_acc *grpstat = ekind->grpstat;
989 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
990 an option, but not supported now. Additionally, if we are doing iterations.
991 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
992 bSavEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't copy over the ekinh_old.
993 If FALSE, we overrwrite it.
996 /* group velocities are calculated in update_ekindata and
997 * accumulated in acumulate_groups.
998 * Now the partial global and groups ekin.
1000 for (g = 0; (g < opts->ngtc); g++)
1005 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
1009 clear_mat(tcstat[g].ekinf);
1013 clear_mat(tcstat[g].ekinh);
1017 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
1020 ekind->dekindl_old = ekind->dekindl;
1022 nthread = gmx_omp_nthreads_get(emntUpdate);
1024 #pragma omp parallel for num_threads(nthread) schedule(static)
1025 for (thread = 0; thread < nthread; thread++)
1027 int start_t, end_t, n;
1035 start_t = md->start + ((thread+0)*md->homenr)/nthread;
1036 end_t = md->start + ((thread+1)*md->homenr)/nthread;
1038 ekin_sum = ekind->ekin_work[thread];
1039 dekindl_sum = ekind->dekindl_work[thread];
1041 for (gt = 0; gt < opts->ngtc; gt++)
1043 clear_mat(ekin_sum[gt]);
1049 for (n = start_t; n < end_t; n++)
1059 hm = 0.5*md->massT[n];
1061 for (d = 0; (d < DIM); d++)
1063 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
1065 for (d = 0; (d < DIM); d++)
1067 for (m = 0; (m < DIM); m++)
1069 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1070 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
1073 if (md->nMassPerturbed && md->bPerturbed[n])
1076 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1082 for (thread = 0; thread < nthread; thread++)
1084 for (g = 0; g < opts->ngtc; g++)
1088 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
1093 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1098 ekind->dekindl += *ekind->dekindl_work[thread];
1101 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1104 static void calc_ke_part_visc(matrix box, rvec x[], rvec v[],
1105 t_grpopts *opts, t_mdatoms *md,
1106 gmx_ekindata_t *ekind,
1107 t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
1109 int start = md->start, homenr = md->homenr;
1110 int g, d, n, m, gt = 0;
1113 t_grp_tcstat *tcstat = ekind->tcstat;
1114 t_cos_acc *cosacc = &(ekind->cosacc);
1119 for (g = 0; g < opts->ngtc; g++)
1121 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1122 clear_mat(ekind->tcstat[g].ekinh);
1124 ekind->dekindl_old = ekind->dekindl;
1126 fac = 2*M_PI/box[ZZ][ZZ];
1129 for (n = start; n < start+homenr; n++)
1135 hm = 0.5*md->massT[n];
1137 /* Note that the times of x and v differ by half a step */
1138 /* MRS -- would have to be changed for VV */
1139 cosz = cos(fac*x[n][ZZ]);
1140 /* Calculate the amplitude of the new velocity profile */
1141 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1143 copy_rvec(v[n], v_corrt);
1144 /* Subtract the profile for the kinetic energy */
1145 v_corrt[XX] -= cosz*cosacc->vcos;
1146 for (d = 0; (d < DIM); d++)
1148 for (m = 0; (m < DIM); m++)
1150 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1153 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1157 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1161 if (md->nPerturbed && md->bPerturbed[n])
1163 dekindl += 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1166 ekind->dekindl = dekindl;
1167 cosacc->mvcos = mvcos;
1169 inc_nrnb(nrnb, eNR_EKIN, homenr);
1172 void calc_ke_part(t_state *state, t_grpopts *opts, t_mdatoms *md,
1173 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
1175 if (ekind->cosacc.cos_accel == 0)
1177 calc_ke_part_normal(state->v, opts, md, ekind, nrnb, bEkinAveVel, bSaveEkinOld);
1181 calc_ke_part_visc(state->box, state->x, state->v, opts, md, ekind, nrnb, bEkinAveVel, bSaveEkinOld);
1185 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1187 ekinstate->ekin_n = ir->opts.ngtc;
1188 snew(ekinstate->ekinh, ekinstate->ekin_n);
1189 snew(ekinstate->ekinf, ekinstate->ekin_n);
1190 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1191 snew(ekinstate->ekinscalef_nhc, ekinstate->ekin_n);
1192 snew(ekinstate->ekinscaleh_nhc, ekinstate->ekin_n);
1193 snew(ekinstate->vscale_nhc, ekinstate->ekin_n);
1194 ekinstate->dekindl = 0;
1195 ekinstate->mvcos = 0;
1198 void update_ekinstate(ekinstate_t *ekinstate, gmx_ekindata_t *ekind)
1202 for (i = 0; i < ekinstate->ekin_n; i++)
1204 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1205 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1206 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1207 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1208 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1209 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1212 copy_mat(ekind->ekin, ekinstate->ekin_total);
1213 ekinstate->dekindl = ekind->dekindl;
1214 ekinstate->mvcos = ekind->cosacc.mvcos;
1218 void restore_ekinstate_from_state(t_commrec *cr,
1219 gmx_ekindata_t *ekind, ekinstate_t *ekinstate)
1225 for (i = 0; i < ekinstate->ekin_n; i++)
1227 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1228 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1229 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1230 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1231 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1232 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1235 copy_mat(ekinstate->ekin_total, ekind->ekin);
1237 ekind->dekindl = ekinstate->dekindl;
1238 ekind->cosacc.mvcos = ekinstate->mvcos;
1239 n = ekinstate->ekin_n;
1244 gmx_bcast(sizeof(n), &n, cr);
1245 for (i = 0; i < n; i++)
1247 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1248 ekind->tcstat[i].ekinh[0], cr);
1249 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1250 ekind->tcstat[i].ekinf[0], cr);
1251 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1252 ekind->tcstat[i].ekinh_old[0], cr);
1254 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1255 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1256 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1257 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1258 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1259 &(ekind->tcstat[i].vscale_nhc), cr);
1261 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1262 ekind->ekin[0], cr);
1264 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1265 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1269 void set_deform_reference_box(gmx_update_t upd, gmx_large_int_t step, matrix box)
1271 upd->deformref_step = step;
1272 copy_mat(box, upd->deformref_box);
1275 static void deform(gmx_update_t upd,
1276 int start, int homenr, rvec x[], matrix box, matrix *scale_tot,
1277 const t_inputrec *ir, gmx_large_int_t step)
1279 matrix bnew, invbox, mu;
1283 elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1284 copy_mat(box, bnew);
1285 for (i = 0; i < DIM; i++)
1287 for (j = 0; j < DIM; j++)
1289 if (ir->deform[i][j] != 0)
1292 upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1296 /* We correct the off-diagonal elements,
1297 * which can grow indefinitely during shearing,
1298 * so the shifts do not get messed up.
1300 for (i = 1; i < DIM; i++)
1302 for (j = i-1; j >= 0; j--)
1304 while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1306 rvec_dec(bnew[i], bnew[j]);
1308 while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1310 rvec_inc(bnew[i], bnew[j]);
1314 m_inv_ur0(box, invbox);
1315 copy_mat(bnew, box);
1316 mmul_ur0(box, invbox, mu);
1318 for (i = start; i < start+homenr; i++)
1320 x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1321 x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1322 x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1326 /* The transposes of the scaling matrices are stored,
1327 * so we need to do matrix multiplication in the inverse order.
1329 mmul_ur0(*scale_tot, mu, *scale_tot);
1333 static void combine_forces(int nstcalclr,
1334 gmx_constr_t constr,
1335 t_inputrec *ir, t_mdatoms *md, t_idef *idef,
1337 gmx_large_int_t step,
1338 t_state *state, gmx_bool bMolPBC,
1339 int start, int nrend,
1340 rvec f[], rvec f_lr[],
1345 /* f contains the short-range forces + the long range forces
1346 * which are stored separately in f_lr.
1349 if (constr != NULL && !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1351 /* We need to constrain the LR forces separately,
1352 * because due to the different pre-factor for the SR and LR
1353 * forces in the update algorithm, we can not determine
1354 * the constraint force for the coordinate constraining.
1355 * Constrain only the additional LR part of the force.
1357 /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1358 constrain(NULL, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
1359 state->x, f_lr, f_lr, bMolPBC, state->box, state->lambda[efptBONDED], NULL,
1360 NULL, NULL, nrnb, econqForce, ir->epc == epcMTTK, state->veta, state->veta);
1363 /* Add nstcalclr-1 times the LR force to the sum of both forces
1364 * and store the result in forces_lr.
1366 nm1 = nstcalclr - 1;
1367 for (i = start; i < nrend; i++)
1369 for (d = 0; d < DIM; d++)
1371 f_lr[i][d] = f[i][d] + nm1*f_lr[i][d];
1376 void update_tcouple(FILE *fplog,
1377 gmx_large_int_t step,
1378 t_inputrec *inputrec,
1380 gmx_ekindata_t *ekind,
1381 gmx_wallcycle_t wcycle,
1387 gmx_bool bTCouple = FALSE;
1389 int i, start, end, homenr, offset;
1391 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1392 if (inputrec->etc != etcNO &&
1393 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec)))
1395 /* We should only couple after a step where energies were determined (for leapfrog versions)
1396 or the step energies are determined, for velocity verlet versions */
1398 if (EI_VV(inputrec->eI))
1406 bTCouple = (inputrec->nsttcouple == 1 ||
1407 do_per_step(step+inputrec->nsttcouple-offset,
1408 inputrec->nsttcouple));
1413 dttc = inputrec->nsttcouple*inputrec->delta_t;
1415 switch (inputrec->etc)
1420 berendsen_tcoupl(inputrec, ekind, dttc);
1423 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1424 state->nosehoover_xi, state->nosehoover_vxi, MassQ);
1427 vrescale_tcoupl(inputrec, ekind, dttc,
1428 state->therm_integral, upd->sd->gaussrand[0]);
1431 /* rescale in place here */
1432 if (EI_VV(inputrec->eI))
1434 rescale_velocities(ekind, md, md->start, md->start+md->homenr, state->v);
1439 /* Set the T scaling lambda to 1 to have no scaling */
1440 for (i = 0; (i < inputrec->opts.ngtc); i++)
1442 ekind->tcstat[i].lambda = 1.0;
1447 void update_pcouple(FILE *fplog,
1448 gmx_large_int_t step,
1449 t_inputrec *inputrec,
1453 gmx_wallcycle_t wcycle,
1457 gmx_bool bPCouple = FALSE;
1461 /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1462 if (inputrec->epc != epcNO && (!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))))
1464 /* We should only couple after a step where energies were determined */
1465 bPCouple = (inputrec->nstpcouple == 1 ||
1466 do_per_step(step+inputrec->nstpcouple-1,
1467 inputrec->nstpcouple));
1470 clear_mat(pcoupl_mu);
1471 for (i = 0; i < DIM; i++)
1473 pcoupl_mu[i][i] = 1.0;
1480 dtpc = inputrec->nstpcouple*inputrec->delta_t;
1482 switch (inputrec->epc)
1484 /* We can always pcoupl, even if we did not sum the energies
1485 * the previous step, since state->pres_prev is only updated
1486 * when the energies have been summed.
1490 case (epcBERENDSEN):
1493 berendsen_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev, state->box,
1497 case (epcPARRINELLORAHMAN):
1498 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1499 state->box, state->box_rel, state->boxv,
1500 M, pcoupl_mu, bInitStep);
1508 static rvec *get_xprime(const t_state *state, gmx_update_t upd)
1510 if (state->nalloc > upd->xp_nalloc)
1512 upd->xp_nalloc = state->nalloc;
1513 srenew(upd->xp, upd->xp_nalloc);
1519 void update_constraints(FILE *fplog,
1520 gmx_large_int_t step,
1521 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1522 t_inputrec *inputrec, /* input record and box stuff */
1523 gmx_ekindata_t *ekind,
1528 rvec force[], /* forces on home particles */
1531 tensor vir, /* tensors for virial and ekin, needed for computing */
1534 gmx_wallcycle_t wcycle,
1536 gmx_constr_t constr,
1538 gmx_bool bFirstHalf,
1542 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1545 int start, homenr, nrend, i, n, m, g, d;
1547 rvec *vbuf, *xprime = NULL;
1554 if (bFirstHalf && !EI_VV(inputrec->eI))
1559 /* for now, SD update is here -- though it really seems like it
1560 should be reformulated as a velocity verlet method, since it has two parts */
1563 homenr = md->homenr;
1564 nrend = start+homenr;
1566 dt = inputrec->delta_t;
1571 * APPLY CONSTRAINTS:
1574 * When doing PR pressure coupling we have to constrain the
1575 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1576 * it is enough to do this once though, since the relative velocities
1577 * after this will be normal to the bond vector
1582 /* clear out constraints before applying */
1583 clear_mat(vir_part);
1585 xprime = get_xprime(state, upd);
1587 bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1588 bLog = (do_per_step(step, inputrec->nstlog) || bLastStep || (step < 0));
1589 bEner = (do_per_step(step, inputrec->nstenergy) || bLastStep);
1590 /* Constrain the coordinates xprime */
1591 wallcycle_start(wcycle, ewcCONSTR);
1592 if (EI_VV(inputrec->eI) && bFirstHalf)
1594 constrain(NULL, bLog, bEner, constr, idef,
1595 inputrec, ekind, cr, step, 1, md,
1596 state->x, state->v, state->v,
1597 bMolPBC, state->box,
1598 state->lambda[efptBONDED], dvdlambda,
1599 NULL, bCalcVir ? &vir_con : NULL, nrnb, econqVeloc,
1600 inputrec->epc == epcMTTK, state->veta, vetanew);
1604 constrain(NULL, bLog, bEner, constr, idef,
1605 inputrec, ekind, cr, step, 1, md,
1606 state->x, xprime, NULL,
1607 bMolPBC, state->box,
1608 state->lambda[efptBONDED], dvdlambda,
1609 state->v, bCalcVir ? &vir_con : NULL, nrnb, econqCoord,
1610 inputrec->epc == epcMTTK, state->veta, state->veta);
1612 wallcycle_stop(wcycle, ewcCONSTR);
1616 dump_it_all(fplog, "After Shake",
1617 state->natoms, state->x, xprime, state->v, force);
1621 if (inputrec->eI == eiSD2)
1623 /* A correction factor eph is needed for the SD constraint force */
1624 /* Here we can, unfortunately, not have proper corrections
1625 * for different friction constants, so we use the first one.
1627 for (i = 0; i < DIM; i++)
1629 for (m = 0; m < DIM; m++)
1631 vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1637 m_add(vir_part, vir_con, vir_part);
1641 pr_rvecs(debug, 0, "constraint virial", vir_part, DIM);
1647 if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1649 xprime = get_xprime(state, upd);
1651 nth = gmx_omp_nthreads_get(emntUpdate);
1653 #pragma omp parallel for num_threads(nth) schedule(static)
1654 for (th = 0; th < nth; th++)
1656 int start_th, end_th;
1658 start_th = start + ((nrend-start)* th )/nth;
1659 end_th = start + ((nrend-start)*(th+1))/nth;
1661 /* The second part of the SD integration */
1662 do_update_sd2(upd->sd, upd->sd->gaussrand[th],
1663 FALSE, start_th, end_th,
1664 inputrec->opts.acc, inputrec->opts.nFreeze,
1665 md->invmass, md->ptype,
1666 md->cFREEZE, md->cACC, md->cTC,
1667 state->x, xprime, state->v, force, state->sd_X,
1668 inputrec->opts.tau_t,
1671 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1675 /* Constrain the coordinates xprime */
1676 wallcycle_start(wcycle, ewcCONSTR);
1677 constrain(NULL, bLog, bEner, constr, idef,
1678 inputrec, NULL, cr, step, 1, md,
1679 state->x, xprime, NULL,
1680 bMolPBC, state->box,
1681 state->lambda[efptBONDED], dvdlambda,
1682 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
1683 wallcycle_stop(wcycle, ewcCONSTR);
1687 /* We must always unshift after updating coordinates; if we did not shake
1688 x was shifted in do_force */
1690 if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1692 if (graph && (graph->nnodes > 0))
1694 unshift_x(graph, state->box, state->x, upd->xp);
1695 if (TRICLINIC(state->box))
1697 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1701 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1706 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntUpdate)) schedule(static)
1707 for (i = start; i < nrend; i++)
1709 copy_rvec(upd->xp[i], state->x[i]);
1713 dump_it_all(fplog, "After unshift",
1714 state->natoms, state->x, upd->xp, state->v, force);
1716 /* ############# END the update of velocities and positions ######### */
1719 void update_box(FILE *fplog,
1720 gmx_large_int_t step,
1721 t_inputrec *inputrec, /* input record and box stuff */
1725 rvec force[], /* forces on home particles */
1729 gmx_wallcycle_t wcycle,
1732 gmx_bool bFirstHalf)
1734 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE;
1737 int start, homenr, nrend, i, n, m, g;
1741 homenr = md->homenr;
1742 nrend = start+homenr;
1745 (inputrec->etc == etcNOSEHOOVER) ||
1746 (inputrec->epc == epcPARRINELLORAHMAN) ||
1747 (inputrec->epc == epcMTTK);
1749 dt = inputrec->delta_t;
1753 /* now update boxes */
1754 switch (inputrec->epc)
1758 case (epcBERENDSEN):
1759 berendsen_pscale(inputrec, pcoupl_mu, state->box, state->box_rel,
1760 start, homenr, state->x, md->cFREEZE, nrnb);
1762 case (epcPARRINELLORAHMAN):
1763 /* The box velocities were updated in do_pr_pcoupl in the update
1764 * iteration, but we dont change the box vectors until we get here
1765 * since we need to be able to shift/unshift above.
1767 for (i = 0; i < DIM; i++)
1769 for (m = 0; m <= i; m++)
1771 state->box[i][m] += dt*state->boxv[i][m];
1774 preserve_box_shape(inputrec, state->box_rel, state->box);
1776 /* Scale the coordinates */
1777 for (n = start; (n < start+homenr); n++)
1779 tmvmul_ur0(pcoupl_mu, state->x[n], state->x[n]);
1783 switch (inputrec->epct)
1785 case (epctISOTROPIC):
1786 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1787 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1788 Side length scales as exp(veta*dt) */
1790 msmul(state->box, exp(state->veta*dt), state->box);
1792 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1793 o If we assume isotropic scaling, and box length scaling
1794 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1795 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1796 determinant of B is L^DIM det(M), and the determinant
1797 of dB/dt is (dL/dT)^DIM det (M). veta will be
1798 (det(dB/dT)/det(B))^(1/3). Then since M =
1799 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1801 msmul(state->box, state->veta, state->boxv);
1811 if ((!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))) && scale_tot)
1813 /* The transposes of the scaling matrices are stored,
1814 * therefore we need to reverse the order in the multiplication.
1816 mmul_ur0(*scale_tot, pcoupl_mu, *scale_tot);
1819 if (DEFORM(*inputrec))
1821 deform(upd, start, homenr, state->x, state->box, scale_tot, inputrec, step);
1824 dump_it_all(fplog, "After update",
1825 state->natoms, state->x, upd->xp, state->v, force);
1828 void update_coords(FILE *fplog,
1829 gmx_large_int_t step,
1830 t_inputrec *inputrec, /* input record and box stuff */
1834 rvec *f, /* forces on home particles */
1838 gmx_ekindata_t *ekind,
1840 gmx_wallcycle_t wcycle,
1844 t_commrec *cr, /* these shouldn't be here -- need to think about it */
1846 gmx_constr_t constr,
1849 gmx_bool bNH, bPR, bLastStep, bLog = FALSE, bEner = FALSE;
1851 real *imass, *imassin;
1854 int start, homenr, nrend, i, j, d, n, m, g;
1855 int blen0, blen1, iatom, jatom, nshake, nsettle, nconstr, nexpand;
1858 rvec *vcom, *xcom, *vall, *xall, *xin, *vin, *forcein, *fall, *xpall, *xprimein, *xprime;
1861 /* Running the velocity half does nothing except for velocity verlet */
1862 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1863 !EI_VV(inputrec->eI))
1865 gmx_incons("update_coords called for velocity without VV integrator");
1869 homenr = md->homenr;
1870 nrend = start+homenr;
1872 xprime = get_xprime(state, upd);
1874 dt = inputrec->delta_t;
1877 /* We need to update the NMR restraint history when time averaging is used */
1878 if (state->flags & (1<<estDISRE_RM3TAV))
1880 update_disres_history(fcd, &state->hist);
1882 if (state->flags & (1<<estORIRE_DTAV))
1884 update_orires_history(fcd, &state->hist);
1888 bNH = inputrec->etc == etcNOSEHOOVER;
1889 bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1891 if (bDoLR && inputrec->nstcalclr > 1 && !EI_VV(inputrec->eI)) /* get this working with VV? */
1893 /* Store the total force + nstcalclr-1 times the LR force
1894 * in forces_lr, so it can be used in a normal update algorithm
1895 * to produce twin time stepping.
1897 /* is this correct in the new construction? MRS */
1898 combine_forces(inputrec->nstcalclr, constr, inputrec, md, idef, cr,
1899 step, state, bMolPBC,
1900 start, nrend, f, f_lr, nrnb);
1908 /* ############# START The update of velocities and positions ######### */
1910 dump_it_all(fplog, "Before update",
1911 state->natoms, state->x, xprime, state->v, force);
1913 if (inputrec->eI == eiSD2)
1915 check_sd2_work_data_allocation(upd->sd, nrend);
1917 do_update_sd2_Tconsts(upd->sd,
1918 inputrec->opts.ngtc,
1919 inputrec->opts.tau_t,
1920 inputrec->opts.ref_t);
1922 if (inputrec->eI == eiBD)
1924 do_update_bd_Tconsts(dt, inputrec->bd_fric,
1925 inputrec->opts.ngtc, inputrec->opts.ref_t,
1929 nth = gmx_omp_nthreads_get(emntUpdate);
1931 #pragma omp parallel for num_threads(nth) schedule(static) private(alpha)
1932 for (th = 0; th < nth; th++)
1934 int start_th, end_th;
1936 start_th = start + ((nrend-start)* th )/nth;
1937 end_th = start + ((nrend-start)*(th+1))/nth;
1939 switch (inputrec->eI)
1942 if (ekind->cosacc.cos_accel == 0)
1944 do_update_md(start_th, end_th, dt,
1945 ekind->tcstat, state->nosehoover_vxi,
1946 ekind->bNEMD, ekind->grpstat, inputrec->opts.acc,
1947 inputrec->opts.nFreeze,
1948 md->invmass, md->ptype,
1949 md->cFREEZE, md->cACC, md->cTC,
1950 state->x, xprime, state->v, force, M,
1955 do_update_visc(start_th, end_th, dt,
1956 ekind->tcstat, state->nosehoover_vxi,
1957 md->invmass, md->ptype,
1958 md->cTC, state->x, xprime, state->v, force, M,
1960 ekind->cosacc.cos_accel,
1966 do_update_sd1(upd->sd, upd->sd->gaussrand[th],
1967 start_th, end_th, dt,
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.ngtc, inputrec->opts.tau_t, inputrec->opts.ref_t);
1975 /* The SD update is done in 2 parts, because an extra constraint step
1978 do_update_sd2(upd->sd, upd->sd->gaussrand[th],
1979 bInitStep, start_th, end_th,
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, state->sd_X,
1984 inputrec->opts.tau_t,
1988 do_update_bd(start_th, end_th, dt,
1989 inputrec->opts.nFreeze, md->invmass, md->ptype,
1990 md->cFREEZE, md->cTC,
1991 state->x, xprime, state->v, force,
1993 upd->sd->bd_rf, upd->sd->gaussrand[th]);
1997 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
2002 do_update_vv_vel(start_th, end_th, dt,
2003 ekind->tcstat, ekind->grpstat,
2004 inputrec->opts.acc, inputrec->opts.nFreeze,
2005 md->invmass, md->ptype,
2006 md->cFREEZE, md->cACC,
2008 (bNH || bPR), state->veta, alpha);
2011 do_update_vv_pos(start_th, end_th, dt,
2012 ekind->tcstat, ekind->grpstat,
2013 inputrec->opts.acc, inputrec->opts.nFreeze,
2014 md->invmass, md->ptype, md->cFREEZE,
2015 state->x, xprime, state->v, force,
2016 (bNH || bPR), state->veta, alpha);
2021 gmx_fatal(FARGS, "Don't know how to update coordinates");
2029 void correct_ekin(FILE *log, int start, int end, rvec v[], rvec vcm, real mass[],
2030 real tmass, tensor ekin)
2033 * This is a debugging routine. It should not be called for production code
2035 * The kinetic energy should calculated according to:
2036 * Ekin = 1/2 m (v-vcm)^2
2037 * However the correction is not always applied, since vcm may not be
2038 * known in time and we compute
2039 * Ekin' = 1/2 m v^2 instead
2040 * This can be corrected afterwards by computing
2041 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
2043 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
2050 /* Local particles */
2053 /* Processor dependent part. */
2055 for (i = start; (i < end); i++)
2059 for (j = 0; (j < DIM); j++)
2065 svmul(1/tmass, vcm, vcm);
2066 svmul(0.5, vcm, hvcm);
2068 for (j = 0; (j < DIM); j++)
2070 for (k = 0; (k < DIM); k++)
2072 dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
2075 pr_rvecs(log, 0, "dekin", dekin, DIM);
2076 pr_rvecs(log, 0, " ekin", ekin, DIM);
2077 fprintf(log, "dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
2078 trace(dekin), trace(ekin), vcm[XX], vcm[YY], vcm[ZZ]);
2079 fprintf(log, "mv = (%8.4f %8.4f %8.4f)\n",
2080 mv[XX], mv[YY], mv[ZZ]);
2083 extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_large_int_t step, t_mdatoms *md, t_state *state, gmx_update_t upd, t_idef *idef, gmx_constr_t constr)
2087 real rate = (ir->delta_t)/ir->opts.tau_t[0];
2088 /* proceed with andersen if 1) it's fixed probability per
2089 particle andersen or 2) it's massive andersen and it's tau_t/dt */
2090 if ((ir->etc == etcANDERSEN) || do_per_step(step, (int)(1.0/rate)))
2092 srenew(upd->randatom, state->nalloc);
2093 srenew(upd->randatom_list, state->nalloc);
2094 if (upd->randatom_list_init == FALSE)
2096 for (i = 0; i < state->nalloc; i++)
2098 upd->randatom[i] = FALSE;
2099 upd->randatom_list[i] = 0;
2101 upd->randatom_list_init = TRUE;
2103 andersen_tcoupl(ir, md, state, upd->sd->gaussrand[0], rate,
2104 (ir->etc == etcANDERSEN) ? idef : NULL,
2105 constr ? get_nblocks(constr) : 0,
2106 constr ? get_sblock(constr) : NULL,
2107 upd->randatom, upd->randatom_list,
2108 upd->sd->randomize_group, upd->sd->boltzfac);