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