2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
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.
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.
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.
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.
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.
37 /* This file is completely threadsafe - keep it that way! */
45 #include "gromacs/utility/smalloc.h"
46 #include "gromacs/utility/futil.h"
48 #include "gromacs/math/vec.h"
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(gmx_mtop_t *mtop, int ngacc, t_grp_acc gstat[])
71 gmx_mtop_atomloop_all_t aloop;
77 groups = &mtop->groups;
78 aloop = gmx_mtop_atomloop_all_init(mtop);
79 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
81 grp = ggrpnr(groups, egcACC, i);
82 if ((grp < 0) && (grp >= ngacc))
84 gmx_incons("Input for acceleration groups wrong");
87 /* This will not work for integrator BD */
88 gstat[grp].mA += atom->m;
89 gstat[grp].mB += atom->mB;
94 void init_ekindata(FILE gmx_unused *log, gmx_mtop_t *mtop, t_grpopts *opts,
95 gmx_ekindata_t *ekind)
100 fprintf(log, "ngtc: %d, ngacc: %d, ngener: %d\n", opts->ngtc, opts->ngacc,
104 /* bNEMD tells if we should remove remove the COM velocity
105 * from the velocities during velocity scaling in T-coupling.
106 * Turn this on when we have multiple acceleration groups
107 * or one accelerated group.
109 ekind->bNEMD = (opts->ngacc > 1 || norm(opts->acc[0]) > 0);
111 ekind->ngtc = opts->ngtc;
112 snew(ekind->tcstat, opts->ngtc);
113 init_grptcstat(opts->ngtc, ekind->tcstat);
114 /* Set Berendsen tcoupl lambda's to 1,
115 * so runs without Berendsen coupling are not affected.
117 for (i = 0; i < opts->ngtc; i++)
119 ekind->tcstat[i].lambda = 1.0;
120 ekind->tcstat[i].vscale_nhc = 1.0;
121 ekind->tcstat[i].ekinscaleh_nhc = 1.0;
122 ekind->tcstat[i].ekinscalef_nhc = 1.0;
125 nthread = gmx_omp_nthreads_get(emntUpdate);
127 snew(ekind->ekin_work_alloc, nthread);
128 snew(ekind->ekin_work, nthread);
129 snew(ekind->dekindl_work, nthread);
130 #pragma omp parallel for num_threads(nthread) schedule(static)
131 for (thread = 0; thread < nthread; thread++)
133 #define EKIN_WORK_BUFFER_SIZE 2
134 /* Allocate 2 extra elements on both sides, so in single
136 * EKIN_WORK_BUFFER_SIZE*DIM*DIM*sizeof(real) = 72/144 bytes
137 * buffer on both sides to avoid cache pollution.
139 snew(ekind->ekin_work_alloc[thread], ekind->ngtc+2*EKIN_WORK_BUFFER_SIZE);
140 ekind->ekin_work[thread] = ekind->ekin_work_alloc[thread] + EKIN_WORK_BUFFER_SIZE;
141 /* Nasty hack so we can have the per-thread accumulation
142 * variable for dekindl in the same thread-local cache lines
143 * as the per-thread accumulation tensors for ekin[fh],
144 * because they are accumulated in the same loop. */
145 ekind->dekindl_work[thread] = &(ekind->ekin_work[thread][ekind->ngtc][0][0]);
146 #undef EKIN_WORK_BUFFER_SIZE
149 ekind->ngacc = opts->ngacc;
150 snew(ekind->grpstat, opts->ngacc);
151 init_grpstat(mtop, opts->ngacc, ekind->grpstat);
154 void accumulate_u(t_commrec *cr, t_grpopts *opts, gmx_ekindata_t *ekind)
156 /* This routine will only be called when it's necessary */
162 for (g = 0; (g < opts->ngacc); g++)
164 add_binr(rb, DIM, ekind->grpstat[g].u);
168 for (g = 0; (g < opts->ngacc); g++)
170 extract_binr(rb, DIM*g, DIM, ekind->grpstat[g].u);
175 /* I don't think accumulate_ekin is used anymore? */
178 static void accumulate_ekin(t_commrec *cr, t_grpopts *opts,
179 gmx_ekindata_t *ekind)
185 for (g = 0; (g < opts->ngtc); g++)
187 gmx_sum(DIM*DIM, ekind->tcstat[g].ekinf[0], cr);
193 void update_ekindata(int start, int homenr, gmx_ekindata_t *ekind,
194 t_grpopts *opts, rvec v[], t_mdatoms *md, real lambda)
199 /* calculate mean velocities at whole timestep */
200 for (g = 0; (g < opts->ngtc); g++)
202 ekind->tcstat[g].T = 0;
207 for (g = 0; (g < opts->ngacc); g++)
209 clear_rvec(ekind->grpstat[g].u);
213 for (n = start; (n < start+homenr); n++)
219 for (d = 0; (d < DIM); d++)
221 mv = md->massT[n]*v[n][d];
222 ekind->grpstat[g].u[d] += mv;
226 for (g = 0; (g < opts->ngacc); g++)
228 for (d = 0; (d < DIM); d++)
230 ekind->grpstat[g].u[d] /=
231 (1-lambda)*ekind->grpstat[g].mA + lambda*ekind->grpstat[g].mB;
237 real sum_ekin(t_grpopts *opts, gmx_ekindata_t *ekind, real *dekindlambda,
238 gmx_bool bEkinAveVel, gmx_bool bScaleEkin)
242 t_grp_tcstat *tcstat;
251 clear_mat(ekind->ekin);
253 for (i = 0; (i < ngtc); i++)
257 tcstat = &ekind->tcstat[i];
258 /* Sometimes a group does not have degrees of freedom, e.g.
259 * when it consists of shells and virtual sites, then we just
260 * set the temperatue to 0 and also neglect the kinetic
261 * energy, which should be zero anyway.
270 /* in this case, kinetic energy is from the current velocities already */
271 msmul(tcstat->ekinf, tcstat->ekinscalef_nhc, tcstat->ekinf);
276 /* Calculate the full step Ekin as the average of the half steps */
277 for (j = 0; (j < DIM); j++)
279 for (m = 0; (m < DIM); m++)
281 tcstat->ekinf[j][m] =
282 0.5*(tcstat->ekinh[j][m]*tcstat->ekinscaleh_nhc + tcstat->ekinh_old[j][m]);
286 m_add(tcstat->ekinf, ekind->ekin, ekind->ekin);
288 tcstat->Th = calc_temp(trace(tcstat->ekinh), nd);
289 tcstat->T = calc_temp(trace(tcstat->ekinf), nd);
291 /* after the scaling factors have been multiplied in, we can remove them */
294 tcstat->ekinscalef_nhc = 1.0;
298 tcstat->ekinscaleh_nhc = 1.0;
317 *dekindlambda = ekind->dekindl;
321 *dekindlambda = 0.5*(ekind->dekindl + ekind->dekindl_old);