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