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,2015,2016,2017,2018,2019, 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! */
44 #include "gromacs/gmxlib/network.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/mdlib/gmx_omp_nthreads.h"
47 #include "gromacs/mdlib/rbin.h"
48 #include "gromacs/mdlib/update.h"
49 #include "gromacs/mdtypes/group.h"
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/mdtypes/mdatom.h"
52 #include "gromacs/topology/mtop_util.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/utility/exceptions.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/futil.h"
57 #include "gromacs/utility/smalloc.h"
59 static void init_grpstat(const gmx_mtop_t *mtop, int ngacc, t_grp_acc gstat[])
63 const gmx_groups_t &groups = mtop->groups;
64 for (const AtomProxy &atomP : AtomRange(*mtop))
66 const t_atom &local = atomP.atom();
67 int i = atomP.globalAtomNumber();
68 int grp = getGroupType(groups, egcACC, i);
69 if ((grp < 0) && (grp >= ngacc))
71 gmx_incons("Input for acceleration groups wrong");
74 /* This will not work for integrator BD */
75 gstat[grp].mA += local.m;
76 gstat[grp].mB += local.mB;
81 void init_ekindata(FILE gmx_unused *log, const gmx_mtop_t *mtop, const t_grpopts *opts,
82 gmx_ekindata_t *ekind)
87 /* bNEMD tells if we should remove remove the COM velocity
88 * from the velocities during velocity scaling in T-coupling.
89 * Turn this on when we have multiple acceleration groups
90 * or one accelerated group.
92 ekind->bNEMD = (opts->ngacc > 1 || norm2(opts->acc[0]) > 0);
94 ekind->ngtc = opts->ngtc;
95 ekind->tcstat.resize(opts->ngtc);
96 /* Set Berendsen tcoupl lambda's to 1,
97 * so runs without Berendsen coupling are not affected.
99 for (i = 0; i < opts->ngtc; i++)
101 ekind->tcstat[i].lambda = 1.0;
102 ekind->tcstat[i].vscale_nhc = 1.0;
103 ekind->tcstat[i].ekinscaleh_nhc = 1.0;
104 ekind->tcstat[i].ekinscalef_nhc = 1.0;
107 nthread = gmx_omp_nthreads_get(emntUpdate);
108 ekind->nthreads = nthread;
109 snew(ekind->ekin_work_alloc, nthread);
110 snew(ekind->ekin_work, nthread);
111 snew(ekind->dekindl_work, nthread);
112 #pragma omp parallel for num_threads(nthread) schedule(static)
113 for (thread = 0; thread < nthread; thread++)
117 #define EKIN_WORK_BUFFER_SIZE 2
118 /* Allocate 2 extra elements on both sides, so in single
120 * EKIN_WORK_BUFFER_SIZE*DIM*DIM*sizeof(real) = 72/144 bytes
121 * buffer on both sides to avoid cache pollution.
123 snew(ekind->ekin_work_alloc[thread], ekind->ngtc+2*EKIN_WORK_BUFFER_SIZE);
124 ekind->ekin_work[thread] = ekind->ekin_work_alloc[thread] + EKIN_WORK_BUFFER_SIZE;
125 /* Nasty hack so we can have the per-thread accumulation
126 * variable for dekindl in the same thread-local cache lines
127 * as the per-thread accumulation tensors for ekin[fh],
128 * because they are accumulated in the same loop. */
129 ekind->dekindl_work[thread] = &(ekind->ekin_work[thread][ekind->ngtc][0][0]);
130 #undef EKIN_WORK_BUFFER_SIZE
132 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
135 ekind->ngacc = opts->ngacc;
136 ekind->grpstat.resize(opts->ngacc);
137 init_grpstat(mtop, opts->ngacc, ekind->grpstat.data());
140 void accumulate_u(const t_commrec *cr, const t_grpopts *opts, gmx_ekindata_t *ekind)
142 /* This routine will only be called when it's necessary */
148 for (g = 0; (g < opts->ngacc); g++)
150 add_binr(rb, DIM, ekind->grpstat[g].u);
154 for (g = 0; (g < opts->ngacc); g++)
156 extract_binr(rb, DIM*g, DIM, ekind->grpstat[g].u);
161 void update_ekindata(int start, int homenr, gmx_ekindata_t *ekind,
162 const t_grpopts *opts, const rvec v[], const t_mdatoms *md, real lambda)
167 /* calculate mean velocities at whole timestep */
168 for (g = 0; (g < opts->ngtc); g++)
170 ekind->tcstat[g].T = 0;
175 for (g = 0; (g < opts->ngacc); g++)
177 clear_rvec(ekind->grpstat[g].u);
181 for (n = start; (n < start+homenr); n++)
187 for (d = 0; (d < DIM); d++)
189 mv = md->massT[n]*v[n][d];
190 ekind->grpstat[g].u[d] += mv;
194 for (g = 0; (g < opts->ngacc); g++)
196 for (d = 0; (d < DIM); d++)
198 ekind->grpstat[g].u[d] /=
199 (1-lambda)*ekind->grpstat[g].mA + lambda*ekind->grpstat[g].mB;
205 real sum_ekin(const t_grpopts *opts, gmx_ekindata_t *ekind, real *dekindlambda,
206 gmx_bool bEkinAveVel, gmx_bool bScaleEkin)
210 t_grp_tcstat *tcstat;
219 clear_mat(ekind->ekin);
221 for (i = 0; (i < ngtc); i++)
225 tcstat = &ekind->tcstat[i];
226 /* Sometimes a group does not have degrees of freedom, e.g.
227 * when it consists of shells and virtual sites, then we just
228 * set the temperatue to 0 and also neglect the kinetic
229 * energy, which should be zero anyway.
238 /* in this case, kinetic energy is from the current velocities already */
239 msmul(tcstat->ekinf, tcstat->ekinscalef_nhc, tcstat->ekinf);
244 /* Calculate the full step Ekin as the average of the half steps */
245 for (j = 0; (j < DIM); j++)
247 for (m = 0; (m < DIM); m++)
249 tcstat->ekinf[j][m] =
250 0.5*(tcstat->ekinh[j][m]*tcstat->ekinscaleh_nhc + tcstat->ekinh_old[j][m]);
254 m_add(tcstat->ekinf, ekind->ekin, ekind->ekin);
256 tcstat->Th = calc_temp(trace(tcstat->ekinh), nd);
257 tcstat->T = calc_temp(trace(tcstat->ekinf), nd);
259 /* after the scaling factors have been multiplied in, we can remove them */
262 tcstat->ekinscalef_nhc = 1.0;
266 tcstat->ekinscaleh_nhc = 1.0;
285 *dekindlambda = ekind->dekindl;
289 *dekindlambda = 0.5*(ekind->dekindl + ekind->dekindl_old);