1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GROwing Monsters And Cloning Shrimps
35 /* This file is completely threadsafe - keep it that way! */
52 #include "mtop_util.h"
53 #include "gmx_omp_nthreads.h"
55 static void init_grptcstat(int ngtc, t_grp_tcstat tcstat[])
59 for (i = 0; (i < ngtc); i++)
62 clear_mat(tcstat[i].ekinh);
63 clear_mat(tcstat[i].ekinh_old);
64 clear_mat(tcstat[i].ekinf);
68 static void init_grpstat(FILE *log,
69 gmx_mtop_t *mtop, int ngacc, t_grp_acc gstat[])
72 gmx_mtop_atomloop_all_t aloop;
78 groups = &mtop->groups;
79 aloop = gmx_mtop_atomloop_all_init(mtop);
80 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
82 grp = ggrpnr(groups, egcACC, i);
83 if ((grp < 0) && (grp >= ngacc))
85 gmx_incons("Input for acceleration groups wrong");
88 /* This will not work for integrator BD */
89 gstat[grp].mA += atom->m;
90 gstat[grp].mB += atom->mB;
95 void init_ekindata(FILE *log, gmx_mtop_t *mtop, t_grpopts *opts,
96 gmx_ekindata_t *ekind)
101 fprintf(log, "ngtc: %d, ngacc: %d, ngener: %d\n", opts->ngtc, opts->ngacc,
105 /* bNEMD tells if we should remove remove the COM velocity
106 * from the velocities during velocity scaling in T-coupling.
107 * Turn this on when we have multiple acceleration groups
108 * or one accelerated group.
110 ekind->bNEMD = (opts->ngacc > 1 || norm(opts->acc[0]) > 0);
112 ekind->ngtc = opts->ngtc;
113 snew(ekind->tcstat, opts->ngtc);
114 init_grptcstat(opts->ngtc, ekind->tcstat);
115 /* Set Berendsen tcoupl lambda's to 1,
116 * so runs without Berendsen coupling are not affected.
118 for (i = 0; i < opts->ngtc; i++)
120 ekind->tcstat[i].lambda = 1.0;
121 ekind->tcstat[i].vscale_nhc = 1.0;
122 ekind->tcstat[i].ekinscaleh_nhc = 1.0;
123 ekind->tcstat[i].ekinscalef_nhc = 1.0;
126 nthread = gmx_omp_nthreads_get(emntUpdate);
128 snew(ekind->ekin_work_alloc, nthread);
129 snew(ekind->ekin_work, nthread);
130 #pragma omp parallel for num_threads(nthread) schedule(static)
131 for (thread = 0; thread < nthread; thread++)
133 /* Allocate 2 elements extra on both sides,
134 * so in single precision we have 2*3*3*4=72 bytes buffer
135 * on both sides to avoid cache pollution.
137 snew(ekind->ekin_work_alloc[thread], ekind->ngtc+4);
138 ekind->ekin_work[thread] = ekind->ekin_work_alloc[thread] + 2;
141 ekind->ngacc = opts->ngacc;
142 snew(ekind->grpstat, opts->ngacc);
143 init_grpstat(log, mtop, opts->ngacc, ekind->grpstat);
146 void accumulate_u(t_commrec *cr, t_grpopts *opts, gmx_ekindata_t *ekind)
148 /* This routine will only be called when it's necessary */
154 for (g = 0; (g < opts->ngacc); g++)
156 add_binr(rb, DIM, ekind->grpstat[g].u);
160 for (g = 0; (g < opts->ngacc); g++)
162 extract_binr(rb, DIM*g, DIM, ekind->grpstat[g].u);
167 /* I don't think accumulate_ekin is used anymore? */
170 static void accumulate_ekin(t_commrec *cr, t_grpopts *opts,
171 gmx_ekindata_t *ekind)
177 for (g = 0; (g < opts->ngtc); g++)
179 gmx_sum(DIM*DIM, ekind->tcstat[g].ekinf[0], cr);
185 void update_ekindata(int start, int homenr, gmx_ekindata_t *ekind,
186 t_grpopts *opts, rvec v[], t_mdatoms *md, real lambda)
191 /* calculate mean velocities at whole timestep */
192 for (g = 0; (g < opts->ngtc); g++)
194 ekind->tcstat[g].T = 0;
199 for (g = 0; (g < opts->ngacc); g++)
201 clear_rvec(ekind->grpstat[g].u);
205 for (n = start; (n < start+homenr); n++)
211 for (d = 0; (d < DIM); d++)
213 mv = md->massT[n]*v[n][d];
214 ekind->grpstat[g].u[d] += mv;
218 for (g = 0; (g < opts->ngacc); g++)
220 for (d = 0; (d < DIM); d++)
222 ekind->grpstat[g].u[d] /=
223 (1-lambda)*ekind->grpstat[g].mA + lambda*ekind->grpstat[g].mB;
229 real sum_ekin(t_grpopts *opts, gmx_ekindata_t *ekind, real *dekindlambda,
230 gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld, gmx_bool bScaleEkin)
234 t_grp_tcstat *tcstat;
243 clear_mat(ekind->ekin);
245 for (i = 0; (i < ngtc); i++)
249 tcstat = &ekind->tcstat[i];
250 /* Sometimes a group does not have degrees of freedom, e.g.
251 * when it consists of shells and virtual sites, then we just
252 * set the temperatue to 0 and also neglect the kinetic
253 * energy, which should be zero anyway.
262 /* in this case, kinetic energy is from the current velocities already */
263 msmul(tcstat->ekinf, tcstat->ekinscalef_nhc, tcstat->ekinf);
269 /* Calculate the full step Ekin as the average of the half steps */
270 for (j = 0; (j < DIM); j++)
272 for (m = 0; (m < DIM); m++)
274 tcstat->ekinf[j][m] =
275 0.5*(tcstat->ekinh[j][m]*tcstat->ekinscaleh_nhc + tcstat->ekinh_old[j][m]);
279 m_add(tcstat->ekinf, ekind->ekin, ekind->ekin);
281 tcstat->Th = calc_temp(trace(tcstat->ekinh), nd);
282 tcstat->T = calc_temp(trace(tcstat->ekinf), nd);
284 /* after the scaling factors have been multiplied in, we can remove them */
287 tcstat->ekinscalef_nhc = 1.0;
291 tcstat->ekinscaleh_nhc = 1.0;
308 *dekindlambda = 0.5*(ekind->dekindl + ekind->dekindl_old);