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"
69 #include "gmx_wallcycle.h"
70 #include "gmx_omp_nthreads.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_large_int_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 t_grp_tcstat *tcstat, t_grp_acc *gstat,
276 rvec accel[], ivec nFreeze[], real invmass[],
277 unsigned short ptype[], unsigned short cFREEZE[],
278 unsigned short cACC[], rvec v[], rvec f[],
279 gmx_bool bExtended, real veta, real alpha)
284 real u, vn, vv, va, vb, vnrel;
290 g = 0.25*dt*veta*alpha;
292 mv2 = series_sinhx(g);
299 for (n = start; n < nrend; n++)
301 w_dt = invmass[n]*dt;
311 for (d = 0; d < DIM; d++)
313 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
315 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
323 } /* do_update_vv_vel */
325 static void do_update_vv_pos(int start, int nrend, double dt,
326 t_grp_tcstat *tcstat, t_grp_acc *gstat,
327 rvec accel[], ivec nFreeze[], real invmass[],
328 unsigned short ptype[], unsigned short cFREEZE[],
329 rvec x[], rvec xprime[], rvec v[],
330 rvec f[], gmx_bool bExtended, real veta, real alpha)
337 /* Would it make more sense if Parrinello-Rahman was put here? */
342 mr2 = series_sinhx(g);
350 for (n = start; n < nrend; n++)
358 for (d = 0; d < DIM; d++)
360 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
362 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
366 xprime[n][d] = x[n][d];
370 } /* do_update_vv_pos */
372 static void do_update_visc(int start, int nrend, double dt,
373 t_grp_tcstat *tcstat,
376 unsigned short ptype[], unsigned short cTC[],
377 rvec x[], rvec xprime[], rvec v[],
378 rvec f[], matrix M, matrix box, real
379 cos_accel, real vcos,
380 gmx_bool bNH, gmx_bool bPR)
385 real lg, vxi = 0, vv;
390 fac = 2*M_PI/(box[ZZ][ZZ]);
394 /* Update with coupling to extended ensembles, used for
395 * Nose-Hoover and Parrinello-Rahman coupling
397 for (n = start; n < nrend; n++)
404 lg = tcstat[gt].lambda;
405 cosz = cos(fac*x[n][ZZ]);
407 copy_rvec(v[n], vrel);
415 for (d = 0; d < DIM; d++)
419 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
421 vn = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
422 - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
425 vn += vc + dt*cosz*cos_accel;
428 xprime[n][d] = x[n][d]+vn*dt;
432 xprime[n][d] = x[n][d];
439 /* Classic version of update, used with berendsen coupling */
440 for (n = start; n < nrend; n++)
442 w_dt = invmass[n]*dt;
447 lg = tcstat[gt].lambda;
448 cosz = cos(fac*x[n][ZZ]);
450 for (d = 0; d < DIM; d++)
454 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
459 /* Do not scale the cosine velocity profile */
460 vv = vc + lg*(vn - vc + f[n][d]*w_dt);
461 /* Add the cosine accelaration profile */
462 vv += dt*cosz*cos_accel;
466 vv = lg*(vn + f[n][d]*w_dt);
469 xprime[n][d] = x[n][d]+vv*dt;
474 xprime[n][d] = x[n][d];
481 /* Allocates and initializes sd->gaussrand[i] for i=1, i<sd->ngaussrand,
482 * Using seeds generated from sd->gaussrand[0].
484 static void init_multiple_gaussrand(gmx_stochd_t *sd)
489 ngr = sd->ngaussrand;
492 for (i = 1; i < ngr; i++)
494 seed[i] = gmx_rng_uniform_uint32(sd->gaussrand[0]);
497 #pragma omp parallel num_threads(ngr)
501 th = gmx_omp_get_thread_num();
504 /* Initialize on each thread to have thread-local memory alloced */
505 sd->gaussrand[th] = gmx_rng_init(seed[th]);
512 static gmx_stochd_t *init_stochd(FILE *fplog, t_inputrec *ir, int nthreads)
521 /* Initiate random number generator for langevin type dynamics,
522 * for BD, SD or velocity rescaling temperature coupling.
524 if (ir->eI == eiBD || EI_SD(ir->eI))
526 sd->ngaussrand = nthreads;
532 snew(sd->gaussrand, sd->ngaussrand);
534 /* Initialize the first random generator */
535 sd->gaussrand[0] = gmx_rng_init(ir->ld_seed);
537 if (sd->ngaussrand > 1)
539 /* Initialize the rest of the random number generators,
540 * using the first one to generate seeds.
542 init_multiple_gaussrand(sd);
545 ngtc = ir->opts.ngtc;
549 snew(sd->bd_rf, ngtc);
551 else if (EI_SD(ir->eI))
554 snew(sd->sdsig, ngtc);
557 for (n = 0; n < ngtc; n++)
559 if (ir->opts.tau_t[n] > 0)
561 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
562 sdc[n].eph = exp(sdc[n].gdt/2);
563 sdc[n].emh = exp(-sdc[n].gdt/2);
564 sdc[n].em = exp(-sdc[n].gdt);
568 /* No friction and noise on this group */
574 if (sdc[n].gdt >= 0.05)
576 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
577 - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
578 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
579 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
584 /* Seventh order expansions for small y */
585 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
586 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))));
587 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
591 fprintf(debug, "SD const tc-grp %d: b %g c %g d %g\n",
592 n, sdc[n].b, sdc[n].c, sdc[n].d);
596 else if (ETC_ANDERSEN(ir->etc))
605 snew(sd->randomize_group, ngtc);
606 snew(sd->boltzfac, ngtc);
608 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
609 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
611 for (n = 0; n < ngtc; n++)
613 reft = max(0.0, opts->ref_t[n]);
614 if ((opts->tau_t[n] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
616 sd->randomize_group[n] = TRUE;
617 sd->boltzfac[n] = BOLTZ*opts->ref_t[n];
621 sd->randomize_group[n] = FALSE;
628 void get_stochd_state(gmx_update_t upd, t_state *state)
630 /* Note that we only get the state of the first random generator,
631 * even if there are multiple. This avoids repetition.
633 gmx_rng_get_state(upd->sd->gaussrand[0], state->ld_rng, state->ld_rngi);
636 void set_stochd_state(gmx_update_t upd, t_state *state)
643 gmx_rng_set_state(sd->gaussrand[0], state->ld_rng, state->ld_rngi[0]);
645 if (sd->ngaussrand > 1)
647 /* We only end up here with SD or BD with OpenMP.
648 * Destroy and reinitialize the rest of the random number generators,
649 * using seeds generated from the first one.
650 * Although this doesn't recover the previous state,
651 * it at least avoids repetition, which is most important.
652 * Exaclty restoring states with all MPI+OpenMP setups is difficult
653 * and as the integrator is random to start with, doesn't gain us much.
655 for (i = 1; i < sd->ngaussrand; i++)
657 gmx_rng_destroy(sd->gaussrand[i]);
660 init_multiple_gaussrand(sd);
664 gmx_update_t init_update(FILE *fplog, t_inputrec *ir)
670 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
672 upd->sd = init_stochd(fplog, ir, gmx_omp_nthreads_get(emntUpdate));
677 upd->randatom = NULL;
678 upd->randatom_list = NULL;
679 upd->randatom_list_init = FALSE; /* we have not yet cleared the data structure at this point */
684 static void do_update_sd1(gmx_stochd_t *sd,
686 int start, int nrend, double dt,
687 rvec accel[], ivec nFreeze[],
688 real invmass[], unsigned short ptype[],
689 unsigned short cFREEZE[], unsigned short cACC[],
690 unsigned short cTC[],
691 rvec x[], rvec xprime[], rvec v[], rvec f[],
693 int ngtc, real tau_t[], real ref_t[])
698 int gf = 0, ga = 0, gt = 0;
705 for (n = 0; n < ngtc; n++)
708 /* The mass is encounted for later, since this differs per atom */
709 sig[n].V = sqrt(kT*(1 - sdc[n].em*sdc[n].em));
712 for (n = start; n < nrend; n++)
714 ism = sqrt(invmass[n]);
728 for (d = 0; d < DIM; d++)
730 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
732 sd_V = ism*sig[gt].V*gmx_rng_gaussian_table(gaussrand);
734 v[n][d] = v[n][d]*sdc[gt].em
735 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
738 xprime[n][d] = x[n][d] + v[n][d]*dt;
743 xprime[n][d] = x[n][d];
749 static void check_sd2_work_data_allocation(gmx_stochd_t *sd, int nrend)
751 if (nrend > sd->sd_V_nalloc)
753 sd->sd_V_nalloc = over_alloc_dd(nrend);
754 srenew(sd->sd_V, sd->sd_V_nalloc);
758 static void do_update_sd2(gmx_stochd_t *sd,
761 int start, int nrend,
762 rvec accel[], ivec nFreeze[],
763 real invmass[], unsigned short ptype[],
764 unsigned short cFREEZE[], unsigned short cACC[],
765 unsigned short cTC[],
766 rvec x[], rvec xprime[], rvec v[], rvec f[],
768 int ngtc, real tau_t[], real ref_t[],
773 /* The random part of the velocity update, generated in the first
774 * half of the update, needs to be remembered for the second half.
778 int gf = 0, ga = 0, gt = 0;
779 real vn = 0, Vmh, Xmh;
789 for (n = 0; n < ngtc; n++)
792 /* The mass is encounted for later, since this differs per atom */
793 sig[n].V = sqrt(kT*(1-sdc[n].em));
794 sig[n].X = sqrt(kT*sqr(tau_t[n])*sdc[n].c);
795 sig[n].Yv = sqrt(kT*sdc[n].b/sdc[n].c);
796 sig[n].Yx = sqrt(kT*sqr(tau_t[n])*sdc[n].b/(1-sdc[n].em));
800 for (n = start; n < nrend; n++)
802 ism = sqrt(invmass[n]);
816 for (d = 0; d < DIM; d++)
822 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
828 sd_X[n][d] = ism*sig[gt].X*gmx_rng_gaussian_table(gaussrand);
830 Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
831 + ism*sig[gt].Yv*gmx_rng_gaussian_table(gaussrand);
832 sd_V[n][d] = ism*sig[gt].V*gmx_rng_gaussian_table(gaussrand);
834 v[n][d] = vn*sdc[gt].em
835 + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
836 + sd_V[n][d] - sdc[gt].em*Vmh;
838 xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
843 /* Correct the velocities for the constraints.
844 * This operation introduces some inaccuracy,
845 * since the velocity is determined from differences in coordinates.
848 (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
850 Xmh = sd_V[n][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
851 + ism*sig[gt].Yx*gmx_rng_gaussian_table(gaussrand);
852 sd_X[n][d] = ism*sig[gt].X*gmx_rng_gaussian_table(gaussrand);
854 xprime[n][d] += sd_X[n][d] - Xmh;
863 xprime[n][d] = x[n][d];
870 static void do_update_bd(int start, int nrend, double dt,
872 real invmass[], unsigned short ptype[],
873 unsigned short cFREEZE[], unsigned short cTC[],
874 rvec x[], rvec xprime[], rvec v[],
875 rvec f[], real friction_coefficient,
876 int ngtc, real tau_t[], real ref_t[],
877 real *rf, gmx_rng_t gaussrand)
879 /* note -- these appear to be full step velocities . . . */
885 if (friction_coefficient != 0)
887 invfr = 1.0/friction_coefficient;
888 for (n = 0; n < ngtc; n++)
890 rf[n] = sqrt(2.0*BOLTZ*ref_t[n]/(friction_coefficient*dt));
895 for (n = 0; n < ngtc; n++)
897 rf[n] = sqrt(2.0*BOLTZ*ref_t[n]);
900 for (n = start; (n < nrend); n++)
910 for (d = 0; (d < DIM); d++)
912 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
914 if (friction_coefficient != 0)
916 vn = invfr*f[n][d] + rf[gt]*gmx_rng_gaussian_table(gaussrand);
920 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
921 vn = 0.5*invmass[n]*f[n][d]*dt
922 + sqrt(0.5*invmass[n])*rf[gt]*gmx_rng_gaussian_table(gaussrand);
926 xprime[n][d] = x[n][d]+vn*dt;
931 xprime[n][d] = x[n][d];
937 static void dump_it_all(FILE *fp, const char *title,
938 int natoms, rvec x[], rvec xp[], rvec v[], rvec f[])
943 fprintf(fp, "%s\n", title);
944 pr_rvecs(fp, 0, "x", x, natoms);
945 pr_rvecs(fp, 0, "xp", xp, natoms);
946 pr_rvecs(fp, 0, "v", v, natoms);
947 pr_rvecs(fp, 0, "f", f, natoms);
952 static void calc_ke_part_normal(rvec v[], t_grpopts *opts, t_mdatoms *md,
953 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel,
954 gmx_bool bSaveEkinOld)
957 t_grp_tcstat *tcstat = ekind->tcstat;
958 t_grp_acc *grpstat = ekind->grpstat;
961 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
962 an option, but not supported now. Additionally, if we are doing iterations.
963 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
964 bSavEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't copy over the ekinh_old.
965 If FALSE, we overrwrite it.
968 /* group velocities are calculated in update_ekindata and
969 * accumulated in acumulate_groups.
970 * Now the partial global and groups ekin.
972 for (g = 0; (g < opts->ngtc); g++)
977 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
981 clear_mat(tcstat[g].ekinf);
985 clear_mat(tcstat[g].ekinh);
989 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
992 ekind->dekindl_old = ekind->dekindl;
994 nthread = gmx_omp_nthreads_get(emntUpdate);
996 #pragma omp parallel for num_threads(nthread) schedule(static)
997 for (thread = 0; thread < nthread; thread++)
999 int start_t, end_t, n;
1007 start_t = md->start + ((thread+0)*md->homenr)/nthread;
1008 end_t = md->start + ((thread+1)*md->homenr)/nthread;
1010 ekin_sum = ekind->ekin_work[thread];
1011 dekindl_sum = &ekind->ekin_work[thread][opts->ngtc][0][0];
1013 for (gt = 0; gt < opts->ngtc; gt++)
1015 clear_mat(ekin_sum[gt]);
1020 for (n = start_t; n < end_t; n++)
1030 hm = 0.5*md->massT[n];
1032 for (d = 0; (d < DIM); d++)
1034 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
1036 for (d = 0; (d < DIM); d++)
1038 for (m = 0; (m < DIM); m++)
1040 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1041 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
1044 if (md->nMassPerturbed && md->bPerturbed[n])
1047 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1053 for (thread = 0; thread < nthread; thread++)
1055 for (g = 0; g < opts->ngtc; g++)
1059 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
1064 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1069 ekind->dekindl += ekind->ekin_work[thread][opts->ngtc][0][0];
1072 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1075 static void calc_ke_part_visc(matrix box, rvec x[], rvec v[],
1076 t_grpopts *opts, t_mdatoms *md,
1077 gmx_ekindata_t *ekind,
1078 t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
1080 int start = md->start, homenr = md->homenr;
1081 int g, d, n, m, gt = 0;
1084 t_grp_tcstat *tcstat = ekind->tcstat;
1085 t_cos_acc *cosacc = &(ekind->cosacc);
1090 for (g = 0; g < opts->ngtc; g++)
1092 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1093 clear_mat(ekind->tcstat[g].ekinh);
1095 ekind->dekindl_old = ekind->dekindl;
1097 fac = 2*M_PI/box[ZZ][ZZ];
1100 for (n = start; n < start+homenr; n++)
1106 hm = 0.5*md->massT[n];
1108 /* Note that the times of x and v differ by half a step */
1109 /* MRS -- would have to be changed for VV */
1110 cosz = cos(fac*x[n][ZZ]);
1111 /* Calculate the amplitude of the new velocity profile */
1112 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1114 copy_rvec(v[n], v_corrt);
1115 /* Subtract the profile for the kinetic energy */
1116 v_corrt[XX] -= cosz*cosacc->vcos;
1117 for (d = 0; (d < DIM); d++)
1119 for (m = 0; (m < DIM); m++)
1121 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1124 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1128 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1132 if (md->nPerturbed && md->bPerturbed[n])
1134 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1137 ekind->dekindl = dekindl;
1138 cosacc->mvcos = mvcos;
1140 inc_nrnb(nrnb, eNR_EKIN, homenr);
1143 void calc_ke_part(t_state *state, t_grpopts *opts, t_mdatoms *md,
1144 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
1146 if (ekind->cosacc.cos_accel == 0)
1148 calc_ke_part_normal(state->v, opts, md, ekind, nrnb, bEkinAveVel, bSaveEkinOld);
1152 calc_ke_part_visc(state->box, state->x, state->v, opts, md, ekind, nrnb, bEkinAveVel, bSaveEkinOld);
1156 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1158 ekinstate->ekin_n = ir->opts.ngtc;
1159 snew(ekinstate->ekinh, ekinstate->ekin_n);
1160 snew(ekinstate->ekinf, ekinstate->ekin_n);
1161 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1162 snew(ekinstate->ekinscalef_nhc, ekinstate->ekin_n);
1163 snew(ekinstate->ekinscaleh_nhc, ekinstate->ekin_n);
1164 snew(ekinstate->vscale_nhc, ekinstate->ekin_n);
1165 ekinstate->dekindl = 0;
1166 ekinstate->mvcos = 0;
1169 void update_ekinstate(ekinstate_t *ekinstate, gmx_ekindata_t *ekind)
1173 for (i = 0; i < ekinstate->ekin_n; i++)
1175 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1176 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1177 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1178 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1179 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1180 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1183 copy_mat(ekind->ekin, ekinstate->ekin_total);
1184 ekinstate->dekindl = ekind->dekindl;
1185 ekinstate->mvcos = ekind->cosacc.mvcos;
1189 void restore_ekinstate_from_state(t_commrec *cr,
1190 gmx_ekindata_t *ekind, ekinstate_t *ekinstate)
1196 for (i = 0; i < ekinstate->ekin_n; i++)
1198 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1199 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1200 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1201 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1202 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1203 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1206 copy_mat(ekinstate->ekin_total, ekind->ekin);
1208 ekind->dekindl = ekinstate->dekindl;
1209 ekind->cosacc.mvcos = ekinstate->mvcos;
1210 n = ekinstate->ekin_n;
1215 gmx_bcast(sizeof(n), &n, cr);
1216 for (i = 0; i < n; i++)
1218 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1219 ekind->tcstat[i].ekinh[0], cr);
1220 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1221 ekind->tcstat[i].ekinf[0], cr);
1222 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1223 ekind->tcstat[i].ekinh_old[0], cr);
1225 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1226 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1227 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1228 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1229 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1230 &(ekind->tcstat[i].vscale_nhc), cr);
1232 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1233 ekind->ekin[0], cr);
1235 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1236 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1240 void set_deform_reference_box(gmx_update_t upd, gmx_large_int_t step, matrix box)
1242 upd->deformref_step = step;
1243 copy_mat(box, upd->deformref_box);
1246 static void deform(gmx_update_t upd,
1247 int start, int homenr, rvec x[], matrix box, matrix *scale_tot,
1248 const t_inputrec *ir, gmx_large_int_t step)
1250 matrix bnew, invbox, mu;
1254 elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1255 copy_mat(box, bnew);
1256 for (i = 0; i < DIM; i++)
1258 for (j = 0; j < DIM; j++)
1260 if (ir->deform[i][j] != 0)
1263 upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1267 /* We correct the off-diagonal elements,
1268 * which can grow indefinitely during shearing,
1269 * so the shifts do not get messed up.
1271 for (i = 1; i < DIM; i++)
1273 for (j = i-1; j >= 0; j--)
1275 while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1277 rvec_dec(bnew[i], bnew[j]);
1279 while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1281 rvec_inc(bnew[i], bnew[j]);
1285 m_inv_ur0(box, invbox);
1286 copy_mat(bnew, box);
1287 mmul_ur0(box, invbox, mu);
1289 for (i = start; i < start+homenr; i++)
1291 x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1292 x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1293 x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1297 /* The transposes of the scaling matrices are stored,
1298 * so we need to do matrix multiplication in the inverse order.
1300 mmul_ur0(*scale_tot, mu, *scale_tot);
1304 static void combine_forces(int nstcalclr,
1305 gmx_constr_t constr,
1306 t_inputrec *ir, t_mdatoms *md, t_idef *idef,
1308 gmx_large_int_t step,
1309 t_state *state, gmx_bool bMolPBC,
1310 int start, int nrend,
1311 rvec f[], rvec f_lr[],
1316 /* f contains the short-range forces + the long range forces
1317 * which are stored separately in f_lr.
1320 if (constr != NULL && !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1322 /* We need to constrain the LR forces separately,
1323 * because due to the different pre-factor for the SR and LR
1324 * forces in the update algorithm, we can not determine
1325 * the constraint force for the coordinate constraining.
1326 * Constrain only the additional LR part of the force.
1328 /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1329 constrain(NULL, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
1330 state->x, f_lr, f_lr, bMolPBC, state->box, state->lambda[efptBONDED], NULL,
1331 NULL, NULL, nrnb, econqForce, ir->epc == epcMTTK, state->veta, state->veta);
1334 /* Add nstcalclr-1 times the LR force to the sum of both forces
1335 * and store the result in forces_lr.
1337 nm1 = nstcalclr - 1;
1338 for (i = start; i < nrend; i++)
1340 for (d = 0; d < DIM; d++)
1342 f_lr[i][d] = f[i][d] + nm1*f_lr[i][d];
1347 void update_tcouple(FILE *fplog,
1348 gmx_large_int_t step,
1349 t_inputrec *inputrec,
1351 gmx_ekindata_t *ekind,
1352 gmx_wallcycle_t wcycle,
1358 gmx_bool bTCouple = FALSE;
1360 int i, start, end, homenr, offset;
1362 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1363 if (inputrec->etc != etcNO &&
1364 !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec)))
1366 /* We should only couple after a step where energies were determined (for leapfrog versions)
1367 or the step energies are determined, for velocity verlet versions */
1369 if (EI_VV(inputrec->eI))
1377 bTCouple = (inputrec->nsttcouple == 1 ||
1378 do_per_step(step+inputrec->nsttcouple-offset,
1379 inputrec->nsttcouple));
1384 dttc = inputrec->nsttcouple*inputrec->delta_t;
1386 switch (inputrec->etc)
1391 berendsen_tcoupl(inputrec, ekind, dttc);
1394 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1395 state->nosehoover_xi, state->nosehoover_vxi, MassQ);
1398 vrescale_tcoupl(inputrec, ekind, dttc,
1399 state->therm_integral, upd->sd->gaussrand[0]);
1402 /* rescale in place here */
1403 if (EI_VV(inputrec->eI))
1405 rescale_velocities(ekind, md, md->start, md->start+md->homenr, state->v);
1410 /* Set the T scaling lambda to 1 to have no scaling */
1411 for (i = 0; (i < inputrec->opts.ngtc); i++)
1413 ekind->tcstat[i].lambda = 1.0;
1418 void update_pcouple(FILE *fplog,
1419 gmx_large_int_t step,
1420 t_inputrec *inputrec,
1424 gmx_wallcycle_t wcycle,
1428 gmx_bool bPCouple = FALSE;
1432 /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1433 if (inputrec->epc != epcNO && (!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))))
1435 /* We should only couple after a step where energies were determined */
1436 bPCouple = (inputrec->nstpcouple == 1 ||
1437 do_per_step(step+inputrec->nstpcouple-1,
1438 inputrec->nstpcouple));
1441 clear_mat(pcoupl_mu);
1442 for (i = 0; i < DIM; i++)
1444 pcoupl_mu[i][i] = 1.0;
1451 dtpc = inputrec->nstpcouple*inputrec->delta_t;
1453 switch (inputrec->epc)
1455 /* We can always pcoupl, even if we did not sum the energies
1456 * the previous step, since state->pres_prev is only updated
1457 * when the energies have been summed.
1461 case (epcBERENDSEN):
1464 berendsen_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev, state->box,
1468 case (epcPARRINELLORAHMAN):
1469 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1470 state->box, state->box_rel, state->boxv,
1471 M, pcoupl_mu, bInitStep);
1479 static rvec *get_xprime(const t_state *state, gmx_update_t upd)
1481 if (state->nalloc > upd->xp_nalloc)
1483 upd->xp_nalloc = state->nalloc;
1484 srenew(upd->xp, upd->xp_nalloc);
1490 void update_constraints(FILE *fplog,
1491 gmx_large_int_t step,
1492 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1493 t_inputrec *inputrec, /* input record and box stuff */
1494 gmx_ekindata_t *ekind,
1499 rvec force[], /* forces on home particles */
1502 tensor vir, /* tensors for virial and ekin, needed for computing */
1505 gmx_wallcycle_t wcycle,
1507 gmx_constr_t constr,
1509 gmx_bool bFirstHalf,
1513 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1516 int start, homenr, nrend, i, n, m, g, d;
1518 rvec *vbuf, *xprime = NULL;
1525 if (bFirstHalf && !EI_VV(inputrec->eI))
1530 /* for now, SD update is here -- though it really seems like it
1531 should be reformulated as a velocity verlet method, since it has two parts */
1534 homenr = md->homenr;
1535 nrend = start+homenr;
1537 dt = inputrec->delta_t;
1542 * APPLY CONSTRAINTS:
1545 * When doing PR pressure coupling we have to constrain the
1546 * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1547 * it is enough to do this once though, since the relative velocities
1548 * after this will be normal to the bond vector
1553 /* clear out constraints before applying */
1554 clear_mat(vir_part);
1556 xprime = get_xprime(state, upd);
1558 bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1559 bLog = (do_per_step(step, inputrec->nstlog) || bLastStep || (step < 0));
1560 bEner = (do_per_step(step, inputrec->nstenergy) || bLastStep);
1561 /* Constrain the coordinates xprime */
1562 wallcycle_start(wcycle, ewcCONSTR);
1563 if (EI_VV(inputrec->eI) && bFirstHalf)
1565 constrain(NULL, bLog, bEner, constr, idef,
1566 inputrec, ekind, cr, step, 1, md,
1567 state->x, state->v, state->v,
1568 bMolPBC, state->box,
1569 state->lambda[efptBONDED], dvdlambda,
1570 NULL, bCalcVir ? &vir_con : NULL, nrnb, econqVeloc,
1571 inputrec->epc == epcMTTK, state->veta, vetanew);
1575 constrain(NULL, bLog, bEner, constr, idef,
1576 inputrec, ekind, cr, step, 1, md,
1577 state->x, xprime, NULL,
1578 bMolPBC, state->box,
1579 state->lambda[efptBONDED], dvdlambda,
1580 state->v, bCalcVir ? &vir_con : NULL, nrnb, econqCoord,
1581 inputrec->epc == epcMTTK, state->veta, state->veta);
1583 wallcycle_stop(wcycle, ewcCONSTR);
1587 dump_it_all(fplog, "After Shake",
1588 state->natoms, state->x, xprime, state->v, force);
1592 if (inputrec->eI == eiSD2)
1594 /* A correction factor eph is needed for the SD constraint force */
1595 /* Here we can, unfortunately, not have proper corrections
1596 * for different friction constants, so we use the first one.
1598 for (i = 0; i < DIM; i++)
1600 for (m = 0; m < DIM; m++)
1602 vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1608 m_add(vir_part, vir_con, vir_part);
1612 pr_rvecs(debug, 0, "constraint virial", vir_part, DIM);
1618 if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1620 xprime = get_xprime(state, upd);
1622 nth = gmx_omp_nthreads_get(emntUpdate);
1624 #pragma omp parallel for num_threads(nth) schedule(static)
1625 for (th = 0; th < nth; th++)
1627 int start_th, end_th;
1629 start_th = start + ((nrend-start)* th )/nth;
1630 end_th = start + ((nrend-start)*(th+1))/nth;
1632 /* The second part of the SD integration */
1633 do_update_sd2(upd->sd, upd->sd->gaussrand[th],
1634 FALSE, start_th, end_th,
1635 inputrec->opts.acc, inputrec->opts.nFreeze,
1636 md->invmass, md->ptype,
1637 md->cFREEZE, md->cACC, md->cTC,
1638 state->x, xprime, state->v, force, state->sd_X,
1639 inputrec->opts.ngtc, inputrec->opts.tau_t,
1640 inputrec->opts.ref_t, FALSE);
1642 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1646 /* Constrain the coordinates xprime */
1647 wallcycle_start(wcycle, ewcCONSTR);
1648 constrain(NULL, bLog, bEner, constr, idef,
1649 inputrec, NULL, cr, step, 1, md,
1650 state->x, xprime, NULL,
1651 bMolPBC, state->box,
1652 state->lambda[efptBONDED], dvdlambda,
1653 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
1654 wallcycle_stop(wcycle, ewcCONSTR);
1658 /* We must always unshift after updating coordinates; if we did not shake
1659 x was shifted in do_force */
1661 if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1663 if (graph && (graph->nnodes > 0))
1665 unshift_x(graph, state->box, state->x, upd->xp);
1666 if (TRICLINIC(state->box))
1668 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1672 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1677 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntUpdate)) schedule(static)
1678 for (i = start; i < nrend; i++)
1680 copy_rvec(upd->xp[i], state->x[i]);
1684 dump_it_all(fplog, "After unshift",
1685 state->natoms, state->x, upd->xp, state->v, force);
1687 /* ############# END the update of velocities and positions ######### */
1690 void update_box(FILE *fplog,
1691 gmx_large_int_t step,
1692 t_inputrec *inputrec, /* input record and box stuff */
1696 rvec force[], /* forces on home particles */
1700 gmx_wallcycle_t wcycle,
1703 gmx_bool bFirstHalf)
1705 gmx_bool bExtended, bLastStep, bLog = FALSE, bEner = FALSE;
1708 int start, homenr, nrend, i, n, m, g;
1712 homenr = md->homenr;
1713 nrend = start+homenr;
1716 (inputrec->etc == etcNOSEHOOVER) ||
1717 (inputrec->epc == epcPARRINELLORAHMAN) ||
1718 (inputrec->epc == epcMTTK);
1720 dt = inputrec->delta_t;
1724 /* now update boxes */
1725 switch (inputrec->epc)
1729 case (epcBERENDSEN):
1730 berendsen_pscale(inputrec, pcoupl_mu, state->box, state->box_rel,
1731 start, homenr, state->x, md->cFREEZE, nrnb);
1733 case (epcPARRINELLORAHMAN):
1734 /* The box velocities were updated in do_pr_pcoupl in the update
1735 * iteration, but we dont change the box vectors until we get here
1736 * since we need to be able to shift/unshift above.
1738 for (i = 0; i < DIM; i++)
1740 for (m = 0; m <= i; m++)
1742 state->box[i][m] += dt*state->boxv[i][m];
1745 preserve_box_shape(inputrec, state->box_rel, state->box);
1747 /* Scale the coordinates */
1748 for (n = start; (n < start+homenr); n++)
1750 tmvmul_ur0(pcoupl_mu, state->x[n], state->x[n]);
1754 switch (inputrec->epct)
1756 case (epctISOTROPIC):
1757 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1758 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1759 Side length scales as exp(veta*dt) */
1761 msmul(state->box, exp(state->veta*dt), state->box);
1763 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1764 o If we assume isotropic scaling, and box length scaling
1765 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1766 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1767 determinant of B is L^DIM det(M), and the determinant
1768 of dB/dt is (dL/dT)^DIM det (M). veta will be
1769 (det(dB/dT)/det(B))^(1/3). Then since M =
1770 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1772 msmul(state->box, state->veta, state->boxv);
1782 if ((!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))) && scale_tot)
1784 /* The transposes of the scaling matrices are stored,
1785 * therefore we need to reverse the order in the multiplication.
1787 mmul_ur0(*scale_tot, pcoupl_mu, *scale_tot);
1790 if (DEFORM(*inputrec))
1792 deform(upd, start, homenr, state->x, state->box, scale_tot, inputrec, step);
1795 dump_it_all(fplog, "After update",
1796 state->natoms, state->x, upd->xp, state->v, force);
1799 void update_coords(FILE *fplog,
1800 gmx_large_int_t step,
1801 t_inputrec *inputrec, /* input record and box stuff */
1805 rvec *f, /* forces on home particles */
1809 gmx_ekindata_t *ekind,
1811 gmx_wallcycle_t wcycle,
1815 t_commrec *cr, /* these shouldn't be here -- need to think about it */
1817 gmx_constr_t constr,
1820 gmx_bool bNH, bPR, bLastStep, bLog = FALSE, bEner = FALSE;
1822 real *imass, *imassin;
1825 int start, homenr, nrend, i, j, d, n, m, g;
1826 int blen0, blen1, iatom, jatom, nshake, nsettle, nconstr, nexpand;
1829 rvec *vcom, *xcom, *vall, *xall, *xin, *vin, *forcein, *fall, *xpall, *xprimein, *xprime;
1832 /* Running the velocity half does nothing except for velocity verlet */
1833 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1834 !EI_VV(inputrec->eI))
1836 gmx_incons("update_coords called for velocity without VV integrator");
1840 homenr = md->homenr;
1841 nrend = start+homenr;
1843 xprime = get_xprime(state, upd);
1845 dt = inputrec->delta_t;
1848 /* We need to update the NMR restraint history when time averaging is used */
1849 if (state->flags & (1<<estDISRE_RM3TAV))
1851 update_disres_history(fcd, &state->hist);
1853 if (state->flags & (1<<estORIRE_DTAV))
1855 update_orires_history(fcd, &state->hist);
1859 bNH = inputrec->etc == etcNOSEHOOVER;
1860 bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1862 if (bDoLR && inputrec->nstcalclr > 1 && !EI_VV(inputrec->eI)) /* get this working with VV? */
1864 /* Store the total force + nstcalclr-1 times the LR force
1865 * in forces_lr, so it can be used in a normal update algorithm
1866 * to produce twin time stepping.
1868 /* is this correct in the new construction? MRS */
1869 combine_forces(inputrec->nstcalclr, constr, inputrec, md, idef, cr,
1870 step, state, bMolPBC,
1871 start, nrend, f, f_lr, nrnb);
1879 /* ############# START The update of velocities and positions ######### */
1881 dump_it_all(fplog, "Before update",
1882 state->natoms, state->x, xprime, state->v, force);
1884 if (EI_RANDOM(inputrec->eI))
1886 /* We still need to take care of generating random seeds properly
1887 * when multi-threading.
1893 nth = gmx_omp_nthreads_get(emntUpdate);
1896 if (inputrec->eI == eiSD2)
1898 check_sd2_work_data_allocation(upd->sd, nrend);
1901 #pragma omp parallel for num_threads(nth) schedule(static) private(alpha)
1902 for (th = 0; th < nth; th++)
1904 int start_th, end_th;
1906 start_th = start + ((nrend-start)* th )/nth;
1907 end_th = start + ((nrend-start)*(th+1))/nth;
1909 switch (inputrec->eI)
1912 if (ekind->cosacc.cos_accel == 0)
1914 do_update_md(start_th, end_th, dt,
1915 ekind->tcstat, state->nosehoover_vxi,
1916 ekind->bNEMD, ekind->grpstat, inputrec->opts.acc,
1917 inputrec->opts.nFreeze,
1918 md->invmass, md->ptype,
1919 md->cFREEZE, md->cACC, md->cTC,
1920 state->x, xprime, state->v, force, M,
1925 do_update_visc(start_th, end_th, dt,
1926 ekind->tcstat, state->nosehoover_vxi,
1927 md->invmass, md->ptype,
1928 md->cTC, state->x, xprime, state->v, force, M,
1930 ekind->cosacc.cos_accel,
1936 do_update_sd1(upd->sd, upd->sd->gaussrand[th],
1937 start_th, end_th, dt,
1938 inputrec->opts.acc, inputrec->opts.nFreeze,
1939 md->invmass, md->ptype,
1940 md->cFREEZE, md->cACC, md->cTC,
1941 state->x, xprime, state->v, force, state->sd_X,
1942 inputrec->opts.ngtc, inputrec->opts.tau_t, inputrec->opts.ref_t);
1945 /* The SD update is done in 2 parts, because an extra constraint step
1948 do_update_sd2(upd->sd, upd->sd->gaussrand[th],
1949 bInitStep, start_th, end_th,
1950 inputrec->opts.acc, inputrec->opts.nFreeze,
1951 md->invmass, md->ptype,
1952 md->cFREEZE, md->cACC, md->cTC,
1953 state->x, xprime, state->v, force, state->sd_X,
1954 inputrec->opts.ngtc, inputrec->opts.tau_t, inputrec->opts.ref_t,
1958 do_update_bd(start_th, end_th, dt,
1959 inputrec->opts.nFreeze, md->invmass, md->ptype,
1960 md->cFREEZE, md->cTC,
1961 state->x, xprime, state->v, force,
1963 inputrec->opts.ngtc, inputrec->opts.tau_t, inputrec->opts.ref_t,
1964 upd->sd->bd_rf, upd->sd->gaussrand[th]);
1968 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
1973 do_update_vv_vel(start_th, end_th, dt,
1974 ekind->tcstat, ekind->grpstat,
1975 inputrec->opts.acc, inputrec->opts.nFreeze,
1976 md->invmass, md->ptype,
1977 md->cFREEZE, md->cACC,
1979 (bNH || bPR), state->veta, alpha);
1982 do_update_vv_pos(start_th, end_th, dt,
1983 ekind->tcstat, ekind->grpstat,
1984 inputrec->opts.acc, inputrec->opts.nFreeze,
1985 md->invmass, md->ptype, md->cFREEZE,
1986 state->x, xprime, state->v, force,
1987 (bNH || bPR), state->veta, alpha);
1992 gmx_fatal(FARGS, "Don't know how to update coordinates");
2000 void correct_ekin(FILE *log, int start, int end, rvec v[], rvec vcm, real mass[],
2001 real tmass, tensor ekin)
2004 * This is a debugging routine. It should not be called for production code
2006 * The kinetic energy should calculated according to:
2007 * Ekin = 1/2 m (v-vcm)^2
2008 * However the correction is not always applied, since vcm may not be
2009 * known in time and we compute
2010 * Ekin' = 1/2 m v^2 instead
2011 * This can be corrected afterwards by computing
2012 * Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
2014 * Ekin = Ekin' - m v vcm + 1/2 m vcm^2
2021 /* Local particles */
2024 /* Processor dependent part. */
2026 for (i = start; (i < end); i++)
2030 for (j = 0; (j < DIM); j++)
2036 svmul(1/tmass, vcm, vcm);
2037 svmul(0.5, vcm, hvcm);
2039 for (j = 0; (j < DIM); j++)
2041 for (k = 0; (k < DIM); k++)
2043 dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
2046 pr_rvecs(log, 0, "dekin", dekin, DIM);
2047 pr_rvecs(log, 0, " ekin", ekin, DIM);
2048 fprintf(log, "dekin = %g, ekin = %g vcm = (%8.4f %8.4f %8.4f)\n",
2049 trace(dekin), trace(ekin), vcm[XX], vcm[YY], vcm[ZZ]);
2050 fprintf(log, "mv = (%8.4f %8.4f %8.4f)\n",
2051 mv[XX], mv[YY], mv[ZZ]);
2054 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)
2058 real rate = (ir->delta_t)/ir->opts.tau_t[0];
2059 /* proceed with andersen if 1) it's fixed probability per
2060 particle andersen or 2) it's massive andersen and it's tau_t/dt */
2061 if ((ir->etc == etcANDERSEN) || do_per_step(step, (int)(1.0/rate)))
2063 srenew(upd->randatom, state->nalloc);
2064 srenew(upd->randatom_list, state->nalloc);
2065 if (upd->randatom_list_init == FALSE)
2067 for (i = 0; i < state->nalloc; i++)
2069 upd->randatom[i] = FALSE;
2070 upd->randatom_list[i] = 0;
2072 upd->randatom_list_init = TRUE;
2074 andersen_tcoupl(ir, md, state, upd->sd->gaussrand[0], rate,
2075 (ir->etc == etcANDERSEN) ? idef : NULL,
2076 constr ? get_nblocks(constr) : 0,
2077 constr ? get_sblock(constr) : NULL,
2078 upd->randatom, upd->randatom_list,
2079 upd->sd->randomize_group, upd->sd->boltzfac);