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