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