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