Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / update.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41
42 #include <stdio.h>
43 #include <math.h>
44
45 #include "types/commrec.h"
46 #include "typedefs.h"
47 #include "nrnb.h"
48 #include "gromacs/math/units.h"
49 #include "macros.h"
50 #include "gromacs/math/vec.h"
51 #include "update.h"
52 #include "gromacs/random/random.h"
53 #include "tgroup.h"
54 #include "force.h"
55 #include "names.h"
56 #include "txtdump.h"
57 #include "mdrun.h"
58 #include "constr.h"
59 #include "disre.h"
60 #include "orires.h"
61 #include "gmx_omp_nthreads.h"
62
63 #include "gromacs/fileio/confio.h"
64 #include "gromacs/pbcutil/mshift.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/pulling/pull.h"
67 #include "gromacs/timing/wallcycle.h"
68 #include "gromacs/utility/futil.h"
69 #include "gromacs/utility/gmxomp.h"
70 #include "gromacs/utility/smalloc.h"
71
72 /*For debugging, start at v(-dt/2) for velolcity verlet -- uncomment next line */
73 /*#define STARTFROMDT2*/
74
75 typedef struct {
76     double gdt;
77     double eph;
78     double emh;
79     double em;
80     double b;
81     double c;
82     double d;
83 } gmx_sd_const_t;
84
85 typedef struct {
86     real V;
87     real X;
88     real Yv;
89     real Yx;
90 } gmx_sd_sigma_t;
91
92 typedef struct {
93     /* BD stuff */
94     real           *bd_rf;
95     /* SD stuff */
96     gmx_sd_const_t *sdc;
97     gmx_sd_sigma_t *sdsig;
98     rvec           *sd_V;
99     int             sd_V_nalloc;
100     /* andersen temperature control stuff */
101     gmx_bool       *randomize_group;
102     real           *boltzfac;
103 } gmx_stochd_t;
104
105 typedef struct gmx_update
106 {
107     gmx_stochd_t *sd;
108     /* xprime for constraint algorithms */
109     rvec         *xp;
110     int           xp_nalloc;
111
112     /* Variables for the deform algorithm */
113     gmx_int64_t     deformref_step;
114     matrix          deformref_box;
115 } t_gmx_update;
116
117
118 static void do_update_md(int start, int nrend, double dt,
119                          t_grp_tcstat *tcstat,
120                          double nh_vxi[],
121                          gmx_bool bNEMD, t_grp_acc *gstat, rvec accel[],
122                          ivec nFreeze[],
123                          real invmass[],
124                          unsigned short ptype[], unsigned short cFREEZE[],
125                          unsigned short cACC[], unsigned short cTC[],
126                          rvec x[], rvec xprime[], rvec v[],
127                          rvec f[], matrix M,
128                          gmx_bool bNH, gmx_bool bPR)
129 {
130     double imass, w_dt;
131     int    gf = 0, ga = 0, gt = 0;
132     rvec   vrel;
133     real   vn, vv, va, vb, vnrel;
134     real   lg, vxi = 0, u;
135     int    n, d;
136
137     if (bNH || bPR)
138     {
139         /* Update with coupling to extended ensembles, used for
140          * Nose-Hoover and Parrinello-Rahman coupling
141          * Nose-Hoover uses the reversible leap-frog integrator from
142          * Holian et al. Phys Rev E 52(3) : 2338, 1995
143          */
144         for (n = start; n < nrend; n++)
145         {
146             imass = invmass[n];
147             if (cFREEZE)
148             {
149                 gf   = cFREEZE[n];
150             }
151             if (cACC)
152             {
153                 ga   = cACC[n];
154             }
155             if (cTC)
156             {
157                 gt   = cTC[n];
158             }
159             lg   = tcstat[gt].lambda;
160             if (bNH)
161             {
162                 vxi   = nh_vxi[gt];
163             }
164             rvec_sub(v[n], gstat[ga].u, vrel);
165
166             for (d = 0; d < DIM; d++)
167             {
168                 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
169                 {
170                     vnrel = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
171                                               - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
172                     /* do not scale the mean velocities u */
173                     vn             = gstat[ga].u[d] + accel[ga][d]*dt + vnrel;
174                     v[n][d]        = vn;
175                     xprime[n][d]   = x[n][d]+vn*dt;
176                 }
177                 else
178                 {
179                     v[n][d]        = 0.0;
180                     xprime[n][d]   = x[n][d];
181                 }
182             }
183         }
184     }
185     else if (cFREEZE != NULL ||
186              nFreeze[0][XX] || nFreeze[0][YY] || nFreeze[0][ZZ] ||
187              bNEMD)
188     {
189         /* Update with Berendsen/v-rescale coupling and freeze or NEMD */
190         for (n = start; n < nrend; n++)
191         {
192             w_dt = invmass[n]*dt;
193             if (cFREEZE)
194             {
195                 gf   = cFREEZE[n];
196             }
197             if (cACC)
198             {
199                 ga   = cACC[n];
200             }
201             if (cTC)
202             {
203                 gt   = cTC[n];
204             }
205             lg   = tcstat[gt].lambda;
206
207             for (d = 0; d < DIM; d++)
208             {
209                 vn             = v[n][d];
210                 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
211                 {
212                     vv             = lg*vn + f[n][d]*w_dt;
213
214                     /* do not scale the mean velocities u */
215                     u              = gstat[ga].u[d];
216                     va             = vv + accel[ga][d]*dt;
217                     vb             = va + (1.0-lg)*u;
218                     v[n][d]        = vb;
219                     xprime[n][d]   = x[n][d]+vb*dt;
220                 }
221                 else
222                 {
223                     v[n][d]        = 0.0;
224                     xprime[n][d]   = x[n][d];
225                 }
226             }
227         }
228     }
229     else
230     {
231         /* Plain update with Berendsen/v-rescale coupling */
232         for (n = start; n < nrend; n++)
233         {
234             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
235             {
236                 w_dt = invmass[n]*dt;
237                 if (cTC)
238                 {
239                     gt = cTC[n];
240                 }
241                 lg = tcstat[gt].lambda;
242
243                 for (d = 0; d < DIM; d++)
244                 {
245                     vn           = lg*v[n][d] + f[n][d]*w_dt;
246                     v[n][d]      = vn;
247                     xprime[n][d] = x[n][d] + vn*dt;
248                 }
249             }
250             else
251             {
252                 for (d = 0; d < DIM; d++)
253                 {
254                     v[n][d]        = 0.0;
255                     xprime[n][d]   = x[n][d];
256                 }
257             }
258         }
259     }
260 }
261
262 static void do_update_vv_vel(int start, int nrend, double dt,
263                              rvec accel[], ivec nFreeze[], real invmass[],
264                              unsigned short ptype[], unsigned short cFREEZE[],
265                              unsigned short cACC[], rvec v[], rvec f[],
266                              gmx_bool bExtended, real veta, real alpha)
267 {
268     double imass, w_dt;
269     int    gf = 0, ga = 0;
270     rvec   vrel;
271     real   u, vn, vv, va, vb, vnrel;
272     int    n, d;
273     double g, mv1, mv2;
274
275     if (bExtended)
276     {
277         g        = 0.25*dt*veta*alpha;
278         mv1      = exp(-g);
279         mv2      = series_sinhx(g);
280     }
281     else
282     {
283         mv1      = 1.0;
284         mv2      = 1.0;
285     }
286     for (n = start; n < nrend; n++)
287     {
288         w_dt = invmass[n]*dt;
289         if (cFREEZE)
290         {
291             gf   = cFREEZE[n];
292         }
293         if (cACC)
294         {
295             ga   = cACC[n];
296         }
297
298         for (d = 0; d < DIM; d++)
299         {
300             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
301             {
302                 v[n][d]             = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
303             }
304             else
305             {
306                 v[n][d]        = 0.0;
307             }
308         }
309     }
310 } /* do_update_vv_vel */
311
312 static void do_update_vv_pos(int start, int nrend, double dt,
313                              ivec nFreeze[],
314                              unsigned short ptype[], unsigned short cFREEZE[],
315                              rvec x[], rvec xprime[], rvec v[],
316                              gmx_bool bExtended, real veta)
317 {
318     double imass, w_dt;
319     int    gf = 0;
320     int    n, d;
321     double g, mr1, mr2;
322
323     /* Would it make more sense if Parrinello-Rahman was put here? */
324     if (bExtended)
325     {
326         g        = 0.5*dt*veta;
327         mr1      = exp(g);
328         mr2      = series_sinhx(g);
329     }
330     else
331     {
332         mr1      = 1.0;
333         mr2      = 1.0;
334     }
335
336     for (n = start; n < nrend; n++)
337     {
338
339         if (cFREEZE)
340         {
341             gf   = cFREEZE[n];
342         }
343
344         for (d = 0; d < DIM; d++)
345         {
346             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
347             {
348                 xprime[n][d]   = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
349             }
350             else
351             {
352                 xprime[n][d]   = x[n][d];
353             }
354         }
355     }
356 } /* do_update_vv_pos */
357
358 static void do_update_visc(int start, int nrend, double dt,
359                            t_grp_tcstat *tcstat,
360                            double nh_vxi[],
361                            real invmass[],
362                            unsigned short ptype[], unsigned short cTC[],
363                            rvec x[], rvec xprime[], rvec v[],
364                            rvec f[], matrix M, matrix box, real
365                            cos_accel, real vcos,
366                            gmx_bool bNH, gmx_bool bPR)
367 {
368     double imass, w_dt;
369     int    gt = 0;
370     real   vn, vc;
371     real   lg, vxi = 0, vv;
372     real   fac, cosz;
373     rvec   vrel;
374     int    n, d;
375
376     fac = 2*M_PI/(box[ZZ][ZZ]);
377
378     if (bNH || bPR)
379     {
380         /* Update with coupling to extended ensembles, used for
381          * Nose-Hoover and Parrinello-Rahman coupling
382          */
383         for (n = start; n < nrend; n++)
384         {
385             imass = invmass[n];
386             if (cTC)
387             {
388                 gt   = cTC[n];
389             }
390             lg   = tcstat[gt].lambda;
391             cosz = cos(fac*x[n][ZZ]);
392
393             copy_rvec(v[n], vrel);
394
395             vc            = cosz*vcos;
396             vrel[XX]     -= vc;
397             if (bNH)
398             {
399                 vxi        = nh_vxi[gt];
400             }
401             for (d = 0; d < DIM; d++)
402             {
403                 vn             = v[n][d];
404
405                 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
406                 {
407                     vn  = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
408                                             - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
409                     if (d == XX)
410                     {
411                         vn += vc + dt*cosz*cos_accel;
412                     }
413                     v[n][d]        = vn;
414                     xprime[n][d]   = x[n][d]+vn*dt;
415                 }
416                 else
417                 {
418                     xprime[n][d]   = x[n][d];
419                 }
420             }
421         }
422     }
423     else
424     {
425         /* Classic version of update, used with berendsen coupling */
426         for (n = start; n < nrend; n++)
427         {
428             w_dt = invmass[n]*dt;
429             if (cTC)
430             {
431                 gt   = cTC[n];
432             }
433             lg   = tcstat[gt].lambda;
434             cosz = cos(fac*x[n][ZZ]);
435
436             for (d = 0; d < DIM; d++)
437             {
438                 vn             = v[n][d];
439
440                 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
441                 {
442                     if (d == XX)
443                     {
444                         vc           = cosz*vcos;
445                         /* Do not scale the cosine velocity profile */
446                         vv           = vc + lg*(vn - vc + f[n][d]*w_dt);
447                         /* Add the cosine accelaration profile */
448                         vv          += dt*cosz*cos_accel;
449                     }
450                     else
451                     {
452                         vv           = lg*(vn + f[n][d]*w_dt);
453                     }
454                     v[n][d]        = vv;
455                     xprime[n][d]   = x[n][d]+vv*dt;
456                 }
457                 else
458                 {
459                     v[n][d]        = 0.0;
460                     xprime[n][d]   = x[n][d];
461                 }
462             }
463         }
464     }
465 }
466
467 static gmx_stochd_t *init_stochd(t_inputrec *ir)
468 {
469     gmx_stochd_t   *sd;
470     gmx_sd_const_t *sdc;
471     int             ngtc, n, th;
472     real            y;
473
474     snew(sd, 1);
475
476     ngtc = ir->opts.ngtc;
477
478     if (ir->eI == eiBD)
479     {
480         snew(sd->bd_rf, ngtc);
481     }
482     else if (EI_SD(ir->eI))
483     {
484         snew(sd->sdc, ngtc);
485         snew(sd->sdsig, ngtc);
486
487         sdc = sd->sdc;
488         for (n = 0; n < ngtc; n++)
489         {
490             if (ir->opts.tau_t[n] > 0)
491             {
492                 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
493                 sdc[n].eph = exp(sdc[n].gdt/2);
494                 sdc[n].emh = exp(-sdc[n].gdt/2);
495                 sdc[n].em  = exp(-sdc[n].gdt);
496             }
497             else
498             {
499                 /* No friction and noise on this group */
500                 sdc[n].gdt = 0;
501                 sdc[n].eph = 1;
502                 sdc[n].emh = 1;
503                 sdc[n].em  = 1;
504             }
505             if (sdc[n].gdt >= 0.05)
506             {
507                 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
508                     - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
509                 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
510                 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
511             }
512             else
513             {
514                 y = sdc[n].gdt/2;
515                 /* Seventh order expansions for small y */
516                 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
517                 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))));
518                 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
519             }
520             if (debug)
521             {
522                 fprintf(debug, "SD const tc-grp %d: b %g  c %g  d %g\n",
523                         n, sdc[n].b, sdc[n].c, sdc[n].d);
524             }
525         }
526     }
527     else if (ETC_ANDERSEN(ir->etc))
528     {
529         int        ngtc;
530         t_grpopts *opts;
531         real       reft;
532
533         opts = &ir->opts;
534         ngtc = opts->ngtc;
535
536         snew(sd->randomize_group, ngtc);
537         snew(sd->boltzfac, ngtc);
538
539         /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
540         /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
541
542         for (n = 0; n < ngtc; n++)
543         {
544             reft = max(0.0, opts->ref_t[n]);
545             if ((opts->tau_t[n] > 0) && (reft > 0))  /* tau_t or ref_t = 0 means that no randomization is done */
546             {
547                 sd->randomize_group[n] = TRUE;
548                 sd->boltzfac[n]        = BOLTZ*opts->ref_t[n];
549             }
550             else
551             {
552                 sd->randomize_group[n] = FALSE;
553             }
554         }
555     }
556     return sd;
557 }
558
559 gmx_update_t init_update(t_inputrec *ir)
560 {
561     t_gmx_update *upd;
562
563     snew(upd, 1);
564
565     if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
566     {
567         upd->sd    = init_stochd(ir);
568     }
569
570     upd->xp        = NULL;
571     upd->xp_nalloc = 0;
572
573     return upd;
574 }
575
576 static void do_update_sd1(gmx_stochd_t *sd,
577                           int start, int nrend, double dt,
578                           rvec accel[], ivec nFreeze[],
579                           real invmass[], unsigned short ptype[],
580                           unsigned short cFREEZE[], unsigned short cACC[],
581                           unsigned short cTC[],
582                           rvec x[], rvec xprime[], rvec v[], rvec f[],
583                           int ngtc, real tau_t[], real ref_t[],
584                           gmx_int64_t step, int seed, int* gatindex)
585 {
586     gmx_sd_const_t *sdc;
587     gmx_sd_sigma_t *sig;
588     real            kT;
589     int             gf = 0, ga = 0, gt = 0;
590     real            ism, sd_V;
591     int             n, d;
592
593     sdc = sd->sdc;
594     sig = sd->sdsig;
595
596     for (n = 0; n < ngtc; n++)
597     {
598         kT = BOLTZ*ref_t[n];
599         /* The mass is encounted for later, since this differs per atom */
600         sig[n].V  = sqrt(kT*(1 - sdc[n].em*sdc[n].em));
601     }
602
603     for (n = start; n < nrend; n++)
604     {
605         real rnd[3];
606         int  ng  = gatindex ? gatindex[n] : n;
607         ism = sqrt(invmass[n]);
608         if (cFREEZE)
609         {
610             gf  = cFREEZE[n];
611         }
612         if (cACC)
613         {
614             ga  = cACC[n];
615         }
616         if (cTC)
617         {
618             gt  = cTC[n];
619         }
620
621         gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
622         for (d = 0; d < DIM; d++)
623         {
624             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
625             {
626                 sd_V = ism*sig[gt].V*rnd[d];
627
628                 v[n][d] = v[n][d]*sdc[gt].em
629                     + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
630                     + sd_V;
631
632                 xprime[n][d] = x[n][d] + v[n][d]*dt;
633             }
634             else
635             {
636                 v[n][d]      = 0.0;
637                 xprime[n][d] = x[n][d];
638             }
639         }
640     }
641 }
642
643 static void check_sd2_work_data_allocation(gmx_stochd_t *sd, int nrend)
644 {
645     if (nrend > sd->sd_V_nalloc)
646     {
647         sd->sd_V_nalloc = over_alloc_dd(nrend);
648         srenew(sd->sd_V, sd->sd_V_nalloc);
649     }
650 }
651
652 static void do_update_sd2_Tconsts(gmx_stochd_t *sd,
653                                   int           ngtc,
654                                   const real    tau_t[],
655                                   const real    ref_t[])
656 {
657     /* This is separated from the update below, because it is single threaded */
658     gmx_sd_const_t *sdc;
659     gmx_sd_sigma_t *sig;
660     int             gt;
661     real            kT;
662
663     sdc = sd->sdc;
664     sig = sd->sdsig;
665
666     for (gt = 0; gt < ngtc; gt++)
667     {
668         kT = BOLTZ*ref_t[gt];
669         /* The mass is encounted for later, since this differs per atom */
670         sig[gt].V  = sqrt(kT*(1-sdc[gt].em));
671         sig[gt].X  = sqrt(kT*sqr(tau_t[gt])*sdc[gt].c);
672         sig[gt].Yv = sqrt(kT*sdc[gt].b/sdc[gt].c);
673         sig[gt].Yx = sqrt(kT*sqr(tau_t[gt])*sdc[gt].b/(1-sdc[gt].em));
674     }
675 }
676
677 static void do_update_sd2(gmx_stochd_t *sd,
678                           gmx_bool bInitStep,
679                           int start, int nrend,
680                           rvec accel[], ivec nFreeze[],
681                           real invmass[], unsigned short ptype[],
682                           unsigned short cFREEZE[], unsigned short cACC[],
683                           unsigned short cTC[],
684                           rvec x[], rvec xprime[], rvec v[], rvec f[],
685                           rvec sd_X[],
686                           const real tau_t[],
687                           gmx_bool bFirstHalf, gmx_int64_t step, int seed,
688                           int* gatindex)
689 {
690     gmx_sd_const_t *sdc;
691     gmx_sd_sigma_t *sig;
692     /* The random part of the velocity update, generated in the first
693      * half of the update, needs to be remembered for the second half.
694      */
695     rvec  *sd_V;
696     real   kT;
697     int    gf = 0, ga = 0, gt = 0;
698     real   vn = 0, Vmh, Xmh;
699     real   ism;
700     int    n, d, ng;
701
702     sdc  = sd->sdc;
703     sig  = sd->sdsig;
704     sd_V = sd->sd_V;
705
706     for (n = start; n < nrend; n++)
707     {
708         real rnd[6], rndi[3];
709         ng  = gatindex ? gatindex[n] : n;
710         ism = sqrt(invmass[n]);
711         if (cFREEZE)
712         {
713             gf  = cFREEZE[n];
714         }
715         if (cACC)
716         {
717             ga  = cACC[n];
718         }
719         if (cTC)
720         {
721             gt  = cTC[n];
722         }
723
724         gmx_rng_cycle_6gaussian_table(step*2+(bFirstHalf ? 1 : 2), ng, seed, RND_SEED_UPDATE, rnd);
725         if (bInitStep)
726         {
727             gmx_rng_cycle_3gaussian_table(step*2, ng, seed, RND_SEED_UPDATE, rndi);
728         }
729         for (d = 0; d < DIM; d++)
730         {
731             if (bFirstHalf)
732             {
733                 vn             = v[n][d];
734             }
735             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
736             {
737                 if (bFirstHalf)
738                 {
739                     if (bInitStep)
740                     {
741                         sd_X[n][d] = ism*sig[gt].X*rndi[d];
742                     }
743                     Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
744                         + ism*sig[gt].Yv*rnd[d*2];
745                     sd_V[n][d] = ism*sig[gt].V*rnd[d*2+1];
746
747                     v[n][d] = vn*sdc[gt].em
748                         + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
749                         + sd_V[n][d] - sdc[gt].em*Vmh;
750
751                     xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
752                 }
753                 else
754                 {
755                     /* Correct the velocities for the constraints.
756                      * This operation introduces some inaccuracy,
757                      * since the velocity is determined from differences in coordinates.
758                      */
759                     v[n][d] =
760                         (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
761
762                     Xmh = sd_V[n][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
763                         + ism*sig[gt].Yx*rnd[d*2];
764                     sd_X[n][d] = ism*sig[gt].X*rnd[d*2+1];
765
766                     xprime[n][d] += sd_X[n][d] - Xmh;
767
768                 }
769             }
770             else
771             {
772                 if (bFirstHalf)
773                 {
774                     v[n][d]        = 0.0;
775                     xprime[n][d]   = x[n][d];
776                 }
777             }
778         }
779     }
780 }
781
782 static void do_update_bd_Tconsts(double dt, real friction_coefficient,
783                                  int ngtc, const real ref_t[],
784                                  real *rf)
785 {
786     /* This is separated from the update below, because it is single threaded */
787     int gt;
788
789     if (friction_coefficient != 0)
790     {
791         for (gt = 0; gt < ngtc; gt++)
792         {
793             rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]/(friction_coefficient*dt));
794         }
795     }
796     else
797     {
798         for (gt = 0; gt < ngtc; gt++)
799         {
800             rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]);
801         }
802     }
803 }
804
805 static void do_update_bd(int start, int nrend, double dt,
806                          ivec nFreeze[],
807                          real invmass[], unsigned short ptype[],
808                          unsigned short cFREEZE[], unsigned short cTC[],
809                          rvec x[], rvec xprime[], rvec v[],
810                          rvec f[], real friction_coefficient,
811                          real *rf, gmx_int64_t step, int seed,
812                          int* gatindex)
813 {
814     /* note -- these appear to be full step velocities . . .  */
815     int    gf = 0, gt = 0;
816     real   vn;
817     real   invfr = 0;
818     int    n, d;
819
820     if (friction_coefficient != 0)
821     {
822         invfr = 1.0/friction_coefficient;
823     }
824
825     for (n = start; (n < nrend); n++)
826     {
827         real rnd[3];
828         int  ng  = gatindex ? gatindex[n] : n;
829
830         if (cFREEZE)
831         {
832             gf = cFREEZE[n];
833         }
834         if (cTC)
835         {
836             gt = cTC[n];
837         }
838         gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
839         for (d = 0; (d < DIM); d++)
840         {
841             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
842             {
843                 if (friction_coefficient != 0)
844                 {
845                     vn = invfr*f[n][d] + rf[gt]*rnd[d];
846                 }
847                 else
848                 {
849                     /* NOTE: invmass = 2/(mass*friction_constant*dt) */
850                     vn = 0.5*invmass[n]*f[n][d]*dt
851                         + sqrt(0.5*invmass[n])*rf[gt]*rnd[d];
852                 }
853
854                 v[n][d]      = vn;
855                 xprime[n][d] = x[n][d]+vn*dt;
856             }
857             else
858             {
859                 v[n][d]      = 0.0;
860                 xprime[n][d] = x[n][d];
861             }
862         }
863     }
864 }
865
866 static void dump_it_all(FILE gmx_unused *fp, const char gmx_unused *title,
867                         int gmx_unused natoms, rvec gmx_unused x[], rvec gmx_unused xp[],
868                         rvec gmx_unused v[], rvec gmx_unused f[])
869 {
870 #ifdef DEBUG
871     if (fp)
872     {
873         fprintf(fp, "%s\n", title);
874         pr_rvecs(fp, 0, "x", x, natoms);
875         pr_rvecs(fp, 0, "xp", xp, natoms);
876         pr_rvecs(fp, 0, "v", v, natoms);
877         pr_rvecs(fp, 0, "f", f, natoms);
878     }
879 #endif
880 }
881
882 static void calc_ke_part_normal(rvec v[], t_grpopts *opts, t_mdatoms *md,
883                                 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel,
884                                 gmx_bool bSaveEkinOld)
885 {
886     int           g;
887     t_grp_tcstat *tcstat  = ekind->tcstat;
888     t_grp_acc    *grpstat = ekind->grpstat;
889     int           nthread, thread;
890
891     /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin.  Leap with AveVel is also
892        an option, but not supported now.  Additionally, if we are doing iterations.
893        bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
894        bSavEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't copy over the ekinh_old.
895        If FALSE, we overrwrite it.
896      */
897
898     /* group velocities are calculated in update_ekindata and
899      * accumulated in acumulate_groups.
900      * Now the partial global and groups ekin.
901      */
902     for (g = 0; (g < opts->ngtc); g++)
903     {
904
905         if (!bSaveEkinOld)
906         {
907             copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
908         }
909         if (bEkinAveVel)
910         {
911             clear_mat(tcstat[g].ekinf);
912         }
913         else
914         {
915             clear_mat(tcstat[g].ekinh);
916         }
917         if (bEkinAveVel)
918         {
919             tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
920         }
921     }
922     ekind->dekindl_old = ekind->dekindl;
923
924     nthread = gmx_omp_nthreads_get(emntUpdate);
925
926 #pragma omp parallel for num_threads(nthread) schedule(static)
927     for (thread = 0; thread < nthread; thread++)
928     {
929         int     start_t, end_t, n;
930         int     ga, gt;
931         rvec    v_corrt;
932         real    hm;
933         int     d, m;
934         matrix *ekin_sum;
935         real   *dekindl_sum;
936
937         start_t = ((thread+0)*md->homenr)/nthread;
938         end_t   = ((thread+1)*md->homenr)/nthread;
939
940         ekin_sum    = ekind->ekin_work[thread];
941         dekindl_sum = ekind->dekindl_work[thread];
942
943         for (gt = 0; gt < opts->ngtc; gt++)
944         {
945             clear_mat(ekin_sum[gt]);
946         }
947         *dekindl_sum = 0.0;
948
949         ga = 0;
950         gt = 0;
951         for (n = start_t; n < end_t; n++)
952         {
953             if (md->cACC)
954             {
955                 ga = md->cACC[n];
956             }
957             if (md->cTC)
958             {
959                 gt = md->cTC[n];
960             }
961             hm   = 0.5*md->massT[n];
962
963             for (d = 0; (d < DIM); d++)
964             {
965                 v_corrt[d]  = v[n][d]  - grpstat[ga].u[d];
966             }
967             for (d = 0; (d < DIM); d++)
968             {
969                 for (m = 0; (m < DIM); m++)
970                 {
971                     /* if we're computing a full step velocity, v_corrt[d] has v(t).  Otherwise, v(t+dt/2) */
972                     ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
973                 }
974             }
975             if (md->nMassPerturbed && md->bPerturbed[n])
976             {
977                 *dekindl_sum +=
978                     0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
979             }
980         }
981     }
982
983     ekind->dekindl = 0;
984     for (thread = 0; thread < nthread; thread++)
985     {
986         for (g = 0; g < opts->ngtc; g++)
987         {
988             if (bEkinAveVel)
989             {
990                 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
991                       tcstat[g].ekinf);
992             }
993             else
994             {
995                 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
996                       tcstat[g].ekinh);
997             }
998         }
999
1000         ekind->dekindl += *ekind->dekindl_work[thread];
1001     }
1002
1003     inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1004 }
1005
1006 static void calc_ke_part_visc(matrix box, rvec x[], rvec v[],
1007                               t_grpopts *opts, t_mdatoms *md,
1008                               gmx_ekindata_t *ekind,
1009                               t_nrnb *nrnb, gmx_bool bEkinAveVel)
1010 {
1011     int           start = 0, homenr = md->homenr;
1012     int           g, d, n, m, gt = 0;
1013     rvec          v_corrt;
1014     real          hm;
1015     t_grp_tcstat *tcstat = ekind->tcstat;
1016     t_cos_acc    *cosacc = &(ekind->cosacc);
1017     real          dekindl;
1018     real          fac, cosz;
1019     double        mvcos;
1020
1021     for (g = 0; g < opts->ngtc; g++)
1022     {
1023         copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1024         clear_mat(ekind->tcstat[g].ekinh);
1025     }
1026     ekind->dekindl_old = ekind->dekindl;
1027
1028     fac     = 2*M_PI/box[ZZ][ZZ];
1029     mvcos   = 0;
1030     dekindl = 0;
1031     for (n = start; n < start+homenr; n++)
1032     {
1033         if (md->cTC)
1034         {
1035             gt = md->cTC[n];
1036         }
1037         hm   = 0.5*md->massT[n];
1038
1039         /* Note that the times of x and v differ by half a step */
1040         /* MRS -- would have to be changed for VV */
1041         cosz         = cos(fac*x[n][ZZ]);
1042         /* Calculate the amplitude of the new velocity profile */
1043         mvcos       += 2*cosz*md->massT[n]*v[n][XX];
1044
1045         copy_rvec(v[n], v_corrt);
1046         /* Subtract the profile for the kinetic energy */
1047         v_corrt[XX] -= cosz*cosacc->vcos;
1048         for (d = 0; (d < DIM); d++)
1049         {
1050             for (m = 0; (m < DIM); m++)
1051             {
1052                 /* if we're computing a full step velocity, v_corrt[d] has v(t).  Otherwise, v(t+dt/2) */
1053                 if (bEkinAveVel)
1054                 {
1055                     tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1056                 }
1057                 else
1058                 {
1059                     tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1060                 }
1061             }
1062         }
1063         if (md->nPerturbed && md->bPerturbed[n])
1064         {
1065             /* The minus sign here might be confusing.
1066              * The kinetic contribution from dH/dl doesn't come from
1067              * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
1068              * where p are the momenta. The difference is only a minus sign.
1069              */
1070             dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1071         }
1072     }
1073     ekind->dekindl = dekindl;
1074     cosacc->mvcos  = mvcos;
1075
1076     inc_nrnb(nrnb, eNR_EKIN, homenr);
1077 }
1078
1079 void calc_ke_part(t_state *state, t_grpopts *opts, t_mdatoms *md,
1080                   gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
1081 {
1082     if (ekind->cosacc.cos_accel == 0)
1083     {
1084         calc_ke_part_normal(state->v, opts, md, ekind, nrnb, bEkinAveVel, bSaveEkinOld);
1085     }
1086     else
1087     {
1088         calc_ke_part_visc(state->box, state->x, state->v, opts, md, ekind, nrnb, bEkinAveVel);
1089     }
1090 }
1091
1092 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1093 {
1094     ekinstate->ekin_n = ir->opts.ngtc;
1095     snew(ekinstate->ekinh, ekinstate->ekin_n);
1096     snew(ekinstate->ekinf, ekinstate->ekin_n);
1097     snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1098     snew(ekinstate->ekinscalef_nhc, ekinstate->ekin_n);
1099     snew(ekinstate->ekinscaleh_nhc, ekinstate->ekin_n);
1100     snew(ekinstate->vscale_nhc, ekinstate->ekin_n);
1101     ekinstate->dekindl = 0;
1102     ekinstate->mvcos   = 0;
1103 }
1104
1105 void update_ekinstate(ekinstate_t *ekinstate, gmx_ekindata_t *ekind)
1106 {
1107     int i;
1108
1109     for (i = 0; i < ekinstate->ekin_n; i++)
1110     {
1111         copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1112         copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1113         copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1114         ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1115         ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1116         ekinstate->vscale_nhc[i]     = ekind->tcstat[i].vscale_nhc;
1117     }
1118
1119     copy_mat(ekind->ekin, ekinstate->ekin_total);
1120     ekinstate->dekindl = ekind->dekindl;
1121     ekinstate->mvcos   = ekind->cosacc.mvcos;
1122
1123 }
1124
1125 void restore_ekinstate_from_state(t_commrec *cr,
1126                                   gmx_ekindata_t *ekind, ekinstate_t *ekinstate)
1127 {
1128     int i, n;
1129
1130     if (MASTER(cr))
1131     {
1132         for (i = 0; i < ekinstate->ekin_n; i++)
1133         {
1134             copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1135             copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1136             copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1137             ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1138             ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1139             ekind->tcstat[i].vscale_nhc     = ekinstate->vscale_nhc[i];
1140         }
1141
1142         copy_mat(ekinstate->ekin_total, ekind->ekin);
1143
1144         ekind->dekindl      = ekinstate->dekindl;
1145         ekind->cosacc.mvcos = ekinstate->mvcos;
1146         n                   = ekinstate->ekin_n;
1147     }
1148
1149     if (PAR(cr))
1150     {
1151         gmx_bcast(sizeof(n), &n, cr);
1152         for (i = 0; i < n; i++)
1153         {
1154             gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1155                       ekind->tcstat[i].ekinh[0], cr);
1156             gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1157                       ekind->tcstat[i].ekinf[0], cr);
1158             gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1159                       ekind->tcstat[i].ekinh_old[0], cr);
1160
1161             gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1162                       &(ekind->tcstat[i].ekinscalef_nhc), cr);
1163             gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1164                       &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1165             gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1166                       &(ekind->tcstat[i].vscale_nhc), cr);
1167         }
1168         gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1169                   ekind->ekin[0], cr);
1170
1171         gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1172         gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1173     }
1174 }
1175
1176 void set_deform_reference_box(gmx_update_t upd, gmx_int64_t step, matrix box)
1177 {
1178     upd->deformref_step = step;
1179     copy_mat(box, upd->deformref_box);
1180 }
1181
1182 static void deform(gmx_update_t upd,
1183                    int start, int homenr, rvec x[], matrix box, matrix *scale_tot,
1184                    const t_inputrec *ir, gmx_int64_t step)
1185 {
1186     matrix bnew, invbox, mu;
1187     real   elapsed_time;
1188     int    i, j;
1189
1190     elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1191     copy_mat(box, bnew);
1192     for (i = 0; i < DIM; i++)
1193     {
1194         for (j = 0; j < DIM; j++)
1195         {
1196             if (ir->deform[i][j] != 0)
1197             {
1198                 bnew[i][j] =
1199                     upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1200             }
1201         }
1202     }
1203     /* We correct the off-diagonal elements,
1204      * which can grow indefinitely during shearing,
1205      * so the shifts do not get messed up.
1206      */
1207     for (i = 1; i < DIM; i++)
1208     {
1209         for (j = i-1; j >= 0; j--)
1210         {
1211             while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1212             {
1213                 rvec_dec(bnew[i], bnew[j]);
1214             }
1215             while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1216             {
1217                 rvec_inc(bnew[i], bnew[j]);
1218             }
1219         }
1220     }
1221     m_inv_ur0(box, invbox);
1222     copy_mat(bnew, box);
1223     mmul_ur0(box, invbox, mu);
1224
1225     for (i = start; i < start+homenr; i++)
1226     {
1227         x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1228         x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1229         x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1230     }
1231     if (*scale_tot)
1232     {
1233         /* The transposes of the scaling matrices are stored,
1234          * so we need to do matrix multiplication in the inverse order.
1235          */
1236         mmul_ur0(*scale_tot, mu, *scale_tot);
1237     }
1238 }
1239
1240 void update_tcouple(gmx_int64_t       step,
1241                     t_inputrec       *inputrec,
1242                     t_state          *state,
1243                     gmx_ekindata_t   *ekind,
1244                     t_extmass        *MassQ,
1245                     t_mdatoms        *md)
1246
1247 {
1248     gmx_bool   bTCouple = FALSE;
1249     real       dttc;
1250     int        i, start, end, homenr, offset;
1251
1252     /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1253     if (inputrec->etc != etcNO &&
1254         !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec)))
1255     {
1256         /* We should only couple after a step where energies were determined (for leapfrog versions)
1257            or the step energies are determined, for velocity verlet versions */
1258
1259         if (EI_VV(inputrec->eI))
1260         {
1261             offset = 0;
1262         }
1263         else
1264         {
1265             offset = 1;
1266         }
1267         bTCouple = (inputrec->nsttcouple == 1 ||
1268                     do_per_step(step+inputrec->nsttcouple-offset,
1269                                 inputrec->nsttcouple));
1270     }
1271
1272     if (bTCouple)
1273     {
1274         dttc = inputrec->nsttcouple*inputrec->delta_t;
1275
1276         switch (inputrec->etc)
1277         {
1278             case etcNO:
1279                 break;
1280             case etcBERENDSEN:
1281                 berendsen_tcoupl(inputrec, ekind, dttc);
1282                 break;
1283             case etcNOSEHOOVER:
1284                 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1285                                   state->nosehoover_xi, state->nosehoover_vxi, MassQ);
1286                 break;
1287             case etcVRESCALE:
1288                 vrescale_tcoupl(inputrec, step, ekind, dttc,
1289                                 state->therm_integral);
1290                 break;
1291         }
1292         /* rescale in place here */
1293         if (EI_VV(inputrec->eI))
1294         {
1295             rescale_velocities(ekind, md, 0, md->homenr, state->v);
1296         }
1297     }
1298     else
1299     {
1300         /* Set the T scaling lambda to 1 to have no scaling */
1301         for (i = 0; (i < inputrec->opts.ngtc); i++)
1302         {
1303             ekind->tcstat[i].lambda = 1.0;
1304         }
1305     }
1306 }
1307
1308 void update_pcouple(FILE             *fplog,
1309                     gmx_int64_t       step,
1310                     t_inputrec       *inputrec,
1311                     t_state          *state,
1312                     matrix            pcoupl_mu,
1313                     matrix            M,
1314                     gmx_bool          bInitStep)
1315 {
1316     gmx_bool   bPCouple = FALSE;
1317     real       dtpc     = 0;
1318     int        i;
1319
1320     /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1321     if (inputrec->epc != epcNO && (!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))))
1322     {
1323         /* We should only couple after a step where energies were determined */
1324         bPCouple = (inputrec->nstpcouple == 1 ||
1325                     do_per_step(step+inputrec->nstpcouple-1,
1326                                 inputrec->nstpcouple));
1327     }
1328
1329     clear_mat(pcoupl_mu);
1330     for (i = 0; i < DIM; i++)
1331     {
1332         pcoupl_mu[i][i] = 1.0;
1333     }
1334
1335     clear_mat(M);
1336
1337     if (bPCouple)
1338     {
1339         dtpc = inputrec->nstpcouple*inputrec->delta_t;
1340
1341         switch (inputrec->epc)
1342         {
1343             /* We can always pcoupl, even if we did not sum the energies
1344              * the previous step, since state->pres_prev is only updated
1345              * when the energies have been summed.
1346              */
1347             case (epcNO):
1348                 break;
1349             case (epcBERENDSEN):
1350                 if (!bInitStep)
1351                 {
1352                     berendsen_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev, state->box,
1353                                      pcoupl_mu);
1354                 }
1355                 break;
1356             case (epcPARRINELLORAHMAN):
1357                 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1358                                         state->box, state->box_rel, state->boxv,
1359                                         M, pcoupl_mu, bInitStep);
1360                 break;
1361             default:
1362                 break;
1363         }
1364     }
1365 }
1366
1367 static rvec *get_xprime(const t_state *state, gmx_update_t upd)
1368 {
1369     if (state->nalloc > upd->xp_nalloc)
1370     {
1371         upd->xp_nalloc = state->nalloc;
1372         srenew(upd->xp, upd->xp_nalloc);
1373     }
1374
1375     return upd->xp;
1376 }
1377
1378 static void combine_forces(gmx_update_t upd,
1379                            int nstcalclr,
1380                            gmx_constr_t constr,
1381                            t_inputrec *ir, t_mdatoms *md, t_idef *idef,
1382                            t_commrec *cr,
1383                            gmx_int64_t step,
1384                            t_state *state, gmx_bool bMolPBC,
1385                            int start, int nrend,
1386                            rvec f[], rvec f_lr[],
1387                            tensor *vir_lr_constr,
1388                            t_nrnb *nrnb)
1389 {
1390     int  i, d;
1391
1392     /* f contains the short-range forces + the long range forces
1393      * which are stored separately in f_lr.
1394      */
1395
1396     if (constr != NULL && vir_lr_constr != NULL &&
1397         !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1398     {
1399         /* We need to constrain the LR forces separately,
1400          * because due to the different pre-factor for the SR and LR
1401          * forces in the update algorithm, we have to correct
1402          * the constraint virial for the nstcalclr-1 extra f_lr.
1403          * Constrain only the additional LR part of the force.
1404          */
1405         /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1406         rvec *xp;
1407         real  fac;
1408         int   gf = 0;
1409
1410         xp  = get_xprime(state, upd);
1411
1412         fac = (nstcalclr - 1)*ir->delta_t*ir->delta_t;
1413
1414         for (i = 0; i < md->homenr; i++)
1415         {
1416             if (md->cFREEZE != NULL)
1417             {
1418                 gf = md->cFREEZE[i];
1419             }
1420             for (d = 0; d < DIM; d++)
1421             {
1422                 if ((md->ptype[i] != eptVSite) &&
1423                     (md->ptype[i] != eptShell) &&
1424                     !ir->opts.nFreeze[gf][d])
1425                 {
1426                     xp[i][d] = state->x[i][d] + fac*f_lr[i][d]*md->invmass[i];
1427                 }
1428             }
1429         }
1430         constrain(NULL, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
1431                   state->x, xp, xp, bMolPBC, state->box, state->lambda[efptBONDED], NULL,
1432                   NULL, vir_lr_constr, nrnb, econqCoord, ir->epc == epcMTTK, state->veta, state->veta);
1433     }
1434
1435     /* Add nstcalclr-1 times the LR force to the sum of both forces
1436      * and store the result in forces_lr.
1437      */
1438     for (i = start; i < nrend; i++)
1439     {
1440         for (d = 0; d < DIM; d++)
1441         {
1442             f_lr[i][d] = f[i][d] + (nstcalclr - 1)*f_lr[i][d];
1443         }
1444     }
1445 }
1446
1447 void update_constraints(FILE             *fplog,
1448                         gmx_int64_t       step,
1449                         real             *dvdlambda, /* the contribution to be added to the bonded interactions */
1450                         t_inputrec       *inputrec,  /* input record and box stuff      */
1451                         gmx_ekindata_t   *ekind,
1452                         t_mdatoms        *md,
1453                         t_state          *state,
1454                         gmx_bool          bMolPBC,
1455                         t_graph          *graph,
1456                         rvec              force[],   /* forces on home particles */
1457                         t_idef           *idef,
1458                         tensor            vir_part,
1459                         t_commrec        *cr,
1460                         t_nrnb           *nrnb,
1461                         gmx_wallcycle_t   wcycle,
1462                         gmx_update_t      upd,
1463                         gmx_constr_t      constr,
1464                         gmx_bool          bFirstHalf,
1465                         gmx_bool          bCalcVir,
1466                         real              vetanew)
1467 {
1468     gmx_bool             bExtended, bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1469     double               dt;
1470     real                 dt_1;
1471     int                  start, homenr, nrend, i, n, m, g, d;
1472     tensor               vir_con;
1473     rvec                *vbuf, *xprime = NULL;
1474     int                  nth, th;
1475
1476     if (constr)
1477     {
1478         bDoConstr = TRUE;
1479     }
1480     if (bFirstHalf && !EI_VV(inputrec->eI))
1481     {
1482         bDoConstr = FALSE;
1483     }
1484
1485     /* for now, SD update is here -- though it really seems like it
1486        should be reformulated as a velocity verlet method, since it has two parts */
1487
1488     start  = 0;
1489     homenr = md->homenr;
1490     nrend  = start+homenr;
1491
1492     dt   = inputrec->delta_t;
1493     dt_1 = 1.0/dt;
1494
1495     /*
1496      *  Steps (7C, 8C)
1497      *  APPLY CONSTRAINTS:
1498      *  BLOCK SHAKE
1499
1500      * When doing PR pressure coupling we have to constrain the
1501      * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1502      * it is enough to do this once though, since the relative velocities
1503      * after this will be normal to the bond vector
1504      */
1505
1506     if (bDoConstr)
1507     {
1508         /* clear out constraints before applying */
1509         clear_mat(vir_part);
1510
1511         xprime = get_xprime(state, upd);
1512
1513         bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1514         bLog      = (do_per_step(step, inputrec->nstlog) || bLastStep || (step < 0));
1515         bEner     = (do_per_step(step, inputrec->nstenergy) || bLastStep);
1516         /* Constrain the coordinates xprime */
1517         wallcycle_start(wcycle, ewcCONSTR);
1518         if (EI_VV(inputrec->eI) && bFirstHalf)
1519         {
1520             constrain(NULL, bLog, bEner, constr, idef,
1521                       inputrec, ekind, cr, step, 1, md,
1522                       state->x, state->v, state->v,
1523                       bMolPBC, state->box,
1524                       state->lambda[efptBONDED], dvdlambda,
1525                       NULL, bCalcVir ? &vir_con : NULL, nrnb, econqVeloc,
1526                       inputrec->epc == epcMTTK, state->veta, vetanew);
1527         }
1528         else
1529         {
1530             constrain(NULL, bLog, bEner, constr, idef,
1531                       inputrec, ekind, cr, step, 1, md,
1532                       state->x, xprime, NULL,
1533                       bMolPBC, state->box,
1534                       state->lambda[efptBONDED], dvdlambda,
1535                       state->v, bCalcVir ? &vir_con : NULL, nrnb, econqCoord,
1536                       inputrec->epc == epcMTTK, state->veta, state->veta);
1537         }
1538         wallcycle_stop(wcycle, ewcCONSTR);
1539
1540         where();
1541
1542         dump_it_all(fplog, "After Shake",
1543                     state->natoms, state->x, xprime, state->v, force);
1544
1545         if (bCalcVir)
1546         {
1547             if (inputrec->eI == eiSD2)
1548             {
1549                 /* A correction factor eph is needed for the SD constraint force */
1550                 /* Here we can, unfortunately, not have proper corrections
1551                  * for different friction constants, so we use the first one.
1552                  */
1553                 for (i = 0; i < DIM; i++)
1554                 {
1555                     for (m = 0; m < DIM; m++)
1556                     {
1557                         vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1558                     }
1559                 }
1560             }
1561             else
1562             {
1563                 m_add(vir_part, vir_con, vir_part);
1564             }
1565             if (debug)
1566             {
1567                 pr_rvecs(debug, 0, "constraint virial", vir_part, DIM);
1568             }
1569         }
1570     }
1571
1572     where();
1573     if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1574     {
1575         xprime = get_xprime(state, upd);
1576
1577         nth = gmx_omp_nthreads_get(emntUpdate);
1578
1579 #pragma omp parallel for num_threads(nth) schedule(static)
1580         for (th = 0; th < nth; th++)
1581         {
1582             int start_th, end_th;
1583
1584             start_th = start + ((nrend-start)* th   )/nth;
1585             end_th   = start + ((nrend-start)*(th+1))/nth;
1586
1587             /* The second part of the SD integration */
1588             do_update_sd2(upd->sd,
1589                           FALSE, start_th, end_th,
1590                           inputrec->opts.acc, inputrec->opts.nFreeze,
1591                           md->invmass, md->ptype,
1592                           md->cFREEZE, md->cACC, md->cTC,
1593                           state->x, xprime, state->v, force, state->sd_X,
1594                           inputrec->opts.tau_t,
1595                           FALSE, step, inputrec->ld_seed,
1596                           DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1597         }
1598         inc_nrnb(nrnb, eNR_UPDATE, homenr);
1599
1600         if (bDoConstr)
1601         {
1602             /* Constrain the coordinates xprime */
1603             wallcycle_start(wcycle, ewcCONSTR);
1604             constrain(NULL, bLog, bEner, constr, idef,
1605                       inputrec, NULL, cr, step, 1, md,
1606                       state->x, xprime, NULL,
1607                       bMolPBC, state->box,
1608                       state->lambda[efptBONDED], dvdlambda,
1609                       NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
1610             wallcycle_stop(wcycle, ewcCONSTR);
1611         }
1612     }
1613
1614     /* We must always unshift after updating coordinates; if we did not shake
1615        x was shifted in do_force */
1616
1617     if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1618     {
1619         if (graph && (graph->nnodes > 0))
1620         {
1621             unshift_x(graph, state->box, state->x, upd->xp);
1622             if (TRICLINIC(state->box))
1623             {
1624                 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1625             }
1626             else
1627             {
1628                 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1629             }
1630         }
1631         else
1632         {
1633 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntUpdate)) schedule(static)
1634             for (i = start; i < nrend; i++)
1635             {
1636                 copy_rvec(upd->xp[i], state->x[i]);
1637             }
1638         }
1639
1640         dump_it_all(fplog, "After unshift",
1641                     state->natoms, state->x, upd->xp, state->v, force);
1642     }
1643 /* ############# END the update of velocities and positions ######### */
1644 }
1645
1646 void update_box(FILE             *fplog,
1647                 gmx_int64_t       step,
1648                 t_inputrec       *inputrec,  /* input record and box stuff      */
1649                 t_mdatoms        *md,
1650                 t_state          *state,
1651                 rvec              force[],   /* forces on home particles */
1652                 matrix           *scale_tot,
1653                 matrix            pcoupl_mu,
1654                 t_nrnb           *nrnb,
1655                 gmx_update_t      upd)
1656 {
1657     gmx_bool             bExtended, bLastStep, bLog = FALSE, bEner = FALSE;
1658     double               dt;
1659     real                 dt_1;
1660     int                  start, homenr, nrend, i, n, m, g;
1661     tensor               vir_con;
1662
1663     start  = 0;
1664     homenr = md->homenr;
1665     nrend  = start+homenr;
1666
1667     bExtended =
1668         (inputrec->etc == etcNOSEHOOVER) ||
1669         (inputrec->epc == epcPARRINELLORAHMAN) ||
1670         (inputrec->epc == epcMTTK);
1671
1672     dt = inputrec->delta_t;
1673
1674     where();
1675
1676     /* now update boxes */
1677     switch (inputrec->epc)
1678     {
1679         case (epcNO):
1680             break;
1681         case (epcBERENDSEN):
1682             berendsen_pscale(inputrec, pcoupl_mu, state->box, state->box_rel,
1683                              start, homenr, state->x, md->cFREEZE, nrnb);
1684             break;
1685         case (epcPARRINELLORAHMAN):
1686             /* The box velocities were updated in do_pr_pcoupl in the update
1687              * iteration, but we dont change the box vectors until we get here
1688              * since we need to be able to shift/unshift above.
1689              */
1690             for (i = 0; i < DIM; i++)
1691             {
1692                 for (m = 0; m <= i; m++)
1693                 {
1694                     state->box[i][m] += dt*state->boxv[i][m];
1695                 }
1696             }
1697             preserve_box_shape(inputrec, state->box_rel, state->box);
1698
1699             /* Scale the coordinates */
1700             for (n = start; (n < start+homenr); n++)
1701             {
1702                 tmvmul_ur0(pcoupl_mu, state->x[n], state->x[n]);
1703             }
1704             break;
1705         case (epcMTTK):
1706             switch (inputrec->epct)
1707             {
1708                 case (epctISOTROPIC):
1709                     /* DIM * eta = ln V.  so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1710                        ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1711                        Side length scales as exp(veta*dt) */
1712
1713                     msmul(state->box, exp(state->veta*dt), state->box);
1714
1715                     /* Relate veta to boxv.  veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1716                        o               If we assume isotropic scaling, and box length scaling
1717                        factor L, then V = L^DIM (det(M)).  So dV/dt = DIM
1718                        L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt.  The
1719                        determinant of B is L^DIM det(M), and the determinant
1720                        of dB/dt is (dL/dT)^DIM det (M).  veta will be
1721                        (det(dB/dT)/det(B))^(1/3).  Then since M =
1722                        B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1723
1724                     msmul(state->box, state->veta, state->boxv);
1725                     break;
1726                 default:
1727                     break;
1728             }
1729             break;
1730         default:
1731             break;
1732     }
1733
1734     if ((!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))) && scale_tot)
1735     {
1736         /* The transposes of the scaling matrices are stored,
1737          * therefore we need to reverse the order in the multiplication.
1738          */
1739         mmul_ur0(*scale_tot, pcoupl_mu, *scale_tot);
1740     }
1741
1742     if (DEFORM(*inputrec))
1743     {
1744         deform(upd, start, homenr, state->x, state->box, scale_tot, inputrec, step);
1745     }
1746     where();
1747     dump_it_all(fplog, "After update",
1748                 state->natoms, state->x, upd->xp, state->v, force);
1749 }
1750
1751 void update_coords(FILE             *fplog,
1752                    gmx_int64_t       step,
1753                    t_inputrec       *inputrec,  /* input record and box stuff   */
1754                    t_mdatoms        *md,
1755                    t_state          *state,
1756                    gmx_bool          bMolPBC,
1757                    rvec             *f,    /* forces on home particles */
1758                    gmx_bool          bDoLR,
1759                    rvec             *f_lr,
1760                    tensor           *vir_lr_constr,
1761                    t_fcdata         *fcd,
1762                    gmx_ekindata_t   *ekind,
1763                    matrix            M,
1764                    gmx_update_t      upd,
1765                    gmx_bool          bInitStep,
1766                    int               UpdatePart,
1767                    t_commrec        *cr, /* these shouldn't be here -- need to think about it */
1768                    t_nrnb           *nrnb,
1769                    gmx_constr_t      constr,
1770                    t_idef           *idef)
1771 {
1772     gmx_bool          bNH, bPR, bLastStep, bLog = FALSE, bEner = FALSE;
1773     double            dt, alpha;
1774     real             *imass, *imassin;
1775     rvec             *force;
1776     real              dt_1;
1777     int               start, homenr, nrend, i, j, d, n, m, g;
1778     int               blen0, blen1, iatom, jatom, nshake, nsettle, nconstr, nexpand;
1779     int              *icom = NULL;
1780     tensor            vir_con;
1781     rvec             *vcom, *xcom, *vall, *xall, *xin, *vin, *forcein, *fall, *xpall, *xprimein, *xprime;
1782     int               nth, th;
1783
1784     /* Running the velocity half does nothing except for velocity verlet */
1785     if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1786         !EI_VV(inputrec->eI))
1787     {
1788         gmx_incons("update_coords called for velocity without VV integrator");
1789     }
1790
1791     start  = 0;
1792     homenr = md->homenr;
1793     nrend  = start+homenr;
1794
1795     xprime = get_xprime(state, upd);
1796
1797     dt   = inputrec->delta_t;
1798     dt_1 = 1.0/dt;
1799
1800     /* We need to update the NMR restraint history when time averaging is used */
1801     if (state->flags & (1<<estDISRE_RM3TAV))
1802     {
1803         update_disres_history(fcd, &state->hist);
1804     }
1805     if (state->flags & (1<<estORIRE_DTAV))
1806     {
1807         update_orires_history(fcd, &state->hist);
1808     }
1809
1810
1811     bNH = inputrec->etc == etcNOSEHOOVER;
1812     bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1813
1814     if (bDoLR && inputrec->nstcalclr > 1 && !EI_VV(inputrec->eI))  /* get this working with VV? */
1815     {
1816         /* Store the total force + nstcalclr-1 times the LR force
1817          * in forces_lr, so it can be used in a normal update algorithm
1818          * to produce twin time stepping.
1819          */
1820         /* is this correct in the new construction? MRS */
1821         combine_forces(upd,
1822                        inputrec->nstcalclr, constr, inputrec, md, idef, cr,
1823                        step, state, bMolPBC,
1824                        start, nrend, f, f_lr, vir_lr_constr, nrnb);
1825         force = f_lr;
1826     }
1827     else
1828     {
1829         force = f;
1830     }
1831
1832     /* ############# START The update of velocities and positions ######### */
1833     where();
1834     dump_it_all(fplog, "Before update",
1835                 state->natoms, state->x, xprime, state->v, force);
1836
1837     if (inputrec->eI == eiSD2)
1838     {
1839         check_sd2_work_data_allocation(upd->sd, nrend);
1840
1841         do_update_sd2_Tconsts(upd->sd,
1842                               inputrec->opts.ngtc,
1843                               inputrec->opts.tau_t,
1844                               inputrec->opts.ref_t);
1845     }
1846     if (inputrec->eI == eiBD)
1847     {
1848         do_update_bd_Tconsts(dt, inputrec->bd_fric,
1849                              inputrec->opts.ngtc, inputrec->opts.ref_t,
1850                              upd->sd->bd_rf);
1851     }
1852
1853     nth = gmx_omp_nthreads_get(emntUpdate);
1854
1855 #pragma omp parallel for num_threads(nth) schedule(static) private(alpha)
1856     for (th = 0; th < nth; th++)
1857     {
1858         int start_th, end_th;
1859
1860         start_th = start + ((nrend-start)* th   )/nth;
1861         end_th   = start + ((nrend-start)*(th+1))/nth;
1862
1863         switch (inputrec->eI)
1864         {
1865             case (eiMD):
1866                 if (ekind->cosacc.cos_accel == 0)
1867                 {
1868                     do_update_md(start_th, end_th, dt,
1869                                  ekind->tcstat, state->nosehoover_vxi,
1870                                  ekind->bNEMD, ekind->grpstat, inputrec->opts.acc,
1871                                  inputrec->opts.nFreeze,
1872                                  md->invmass, md->ptype,
1873                                  md->cFREEZE, md->cACC, md->cTC,
1874                                  state->x, xprime, state->v, force, M,
1875                                  bNH, bPR);
1876                 }
1877                 else
1878                 {
1879                     do_update_visc(start_th, end_th, dt,
1880                                    ekind->tcstat, state->nosehoover_vxi,
1881                                    md->invmass, md->ptype,
1882                                    md->cTC, state->x, xprime, state->v, force, M,
1883                                    state->box,
1884                                    ekind->cosacc.cos_accel,
1885                                    ekind->cosacc.vcos,
1886                                    bNH, bPR);
1887                 }
1888                 break;
1889             case (eiSD1):
1890                 do_update_sd1(upd->sd,
1891                               start_th, end_th, dt,
1892                               inputrec->opts.acc, inputrec->opts.nFreeze,
1893                               md->invmass, md->ptype,
1894                               md->cFREEZE, md->cACC, md->cTC,
1895                               state->x, xprime, state->v, force,
1896                               inputrec->opts.ngtc, inputrec->opts.tau_t, inputrec->opts.ref_t,
1897                               step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1898                 break;
1899             case (eiSD2):
1900                 /* The SD update is done in 2 parts, because an extra constraint step
1901                  * is needed
1902                  */
1903                 do_update_sd2(upd->sd,
1904                               bInitStep, start_th, end_th,
1905                               inputrec->opts.acc, inputrec->opts.nFreeze,
1906                               md->invmass, md->ptype,
1907                               md->cFREEZE, md->cACC, md->cTC,
1908                               state->x, xprime, state->v, force, state->sd_X,
1909                               inputrec->opts.tau_t,
1910                               TRUE, step, inputrec->ld_seed,
1911                               DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1912                 break;
1913             case (eiBD):
1914                 do_update_bd(start_th, end_th, dt,
1915                              inputrec->opts.nFreeze, md->invmass, md->ptype,
1916                              md->cFREEZE, md->cTC,
1917                              state->x, xprime, state->v, force,
1918                              inputrec->bd_fric,
1919                              upd->sd->bd_rf,
1920                              step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1921                 break;
1922             case (eiVV):
1923             case (eiVVAK):
1924                 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
1925                 switch (UpdatePart)
1926                 {
1927                     case etrtVELOCITY1:
1928                     case etrtVELOCITY2:
1929                         do_update_vv_vel(start_th, end_th, dt,
1930                                          inputrec->opts.acc, inputrec->opts.nFreeze,
1931                                          md->invmass, md->ptype,
1932                                          md->cFREEZE, md->cACC,
1933                                          state->v, force,
1934                                          (bNH || bPR), state->veta, alpha);
1935                         break;
1936                     case etrtPOSITION:
1937                         do_update_vv_pos(start_th, end_th, dt,
1938                                          inputrec->opts.nFreeze,
1939                                          md->ptype, md->cFREEZE,
1940                                          state->x, xprime, state->v,
1941                                          (bNH || bPR), state->veta);
1942                         break;
1943                 }
1944                 break;
1945             default:
1946                 gmx_fatal(FARGS, "Don't know how to update coordinates");
1947                 break;
1948         }
1949     }
1950
1951 }
1952
1953
1954 void correct_ekin(FILE *log, int start, int end, rvec v[], rvec vcm, real mass[],
1955                   real tmass, tensor ekin)
1956 {
1957     /*
1958      * This is a debugging routine. It should not be called for production code
1959      *
1960      * The kinetic energy should calculated according to:
1961      *   Ekin = 1/2 m (v-vcm)^2
1962      * However the correction is not always applied, since vcm may not be
1963      * known in time and we compute
1964      *   Ekin' = 1/2 m v^2 instead
1965      * This can be corrected afterwards by computing
1966      *   Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
1967      * or in hsorthand:
1968      *   Ekin = Ekin' - m v vcm + 1/2 m vcm^2
1969      */
1970     int    i, j, k;
1971     real   m, tm;
1972     rvec   hvcm, mv;
1973     tensor dekin;
1974
1975     /* Local particles */
1976     clear_rvec(mv);
1977
1978     /* Processor dependent part. */
1979     tm = 0;
1980     for (i = start; (i < end); i++)
1981     {
1982         m      = mass[i];
1983         tm    += m;
1984         for (j = 0; (j < DIM); j++)
1985         {
1986             mv[j] += m*v[i][j];
1987         }
1988     }
1989     /* Shortcut */
1990     svmul(1/tmass, vcm, vcm);
1991     svmul(0.5, vcm, hvcm);
1992     clear_mat(dekin);
1993     for (j = 0; (j < DIM); j++)
1994     {
1995         for (k = 0; (k < DIM); k++)
1996         {
1997             dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
1998         }
1999     }
2000     pr_rvecs(log, 0, "dekin", dekin, DIM);
2001     pr_rvecs(log, 0, " ekin", ekin, DIM);
2002     fprintf(log, "dekin = %g, ekin = %g  vcm = (%8.4f %8.4f %8.4f)\n",
2003             trace(dekin), trace(ekin), vcm[XX], vcm[YY], vcm[ZZ]);
2004     fprintf(log, "mv = (%8.4f %8.4f %8.4f)\n",
2005             mv[XX], mv[YY], mv[ZZ]);
2006 }
2007
2008 extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_int64_t step, const t_commrec *cr,
2009                                             t_mdatoms *md, t_state *state, gmx_update_t upd, gmx_constr_t constr)
2010 {
2011
2012     int  i;
2013     real rate = (ir->delta_t)/ir->opts.tau_t[0];
2014
2015     if (ir->etc == etcANDERSEN && constr != NULL)
2016     {
2017         gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
2018     }
2019
2020     /* proceed with andersen if 1) it's fixed probability per
2021        particle andersen or 2) it's massive andersen and it's tau_t/dt */
2022     if ((ir->etc == etcANDERSEN) || do_per_step(step, (int)(1.0/rate)))
2023     {
2024         andersen_tcoupl(ir, step, cr, md, state, rate,
2025                         upd->sd->randomize_group, upd->sd->boltzfac);
2026         return TRUE;
2027     }
2028     return FALSE;
2029 }