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, 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_grptcstat(int ngtc
, t_grp_tcstat tcstat
[])
63 for (i
= 0; (i
< ngtc
); i
++)
66 clear_mat(tcstat
[i
].ekinh
);
67 clear_mat(tcstat
[i
].ekinh_old
);
68 clear_mat(tcstat
[i
].ekinf
);
72 static void init_grpstat(const gmx_mtop_t
*mtop
, int ngacc
, t_grp_acc gstat
[])
74 const gmx_groups_t
*groups
;
75 gmx_mtop_atomloop_all_t aloop
;
81 groups
= &mtop
->groups
;
82 aloop
= gmx_mtop_atomloop_all_init(mtop
);
83 while (gmx_mtop_atomloop_all_next(aloop
, &i
, &atom
))
85 grp
= getGroupType(groups
, egcACC
, i
);
86 if ((grp
< 0) && (grp
>= ngacc
))
88 gmx_incons("Input for acceleration groups wrong");
91 /* This will not work for integrator BD */
92 gstat
[grp
].mA
+= atom
->m
;
93 gstat
[grp
].mB
+= atom
->mB
;
98 void init_ekindata(FILE gmx_unused
*log
, const gmx_mtop_t
*mtop
, const t_grpopts
*opts
,
99 gmx_ekindata_t
*ekind
)
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 || norm2(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
++)
135 #define EKIN_WORK_BUFFER_SIZE 2
136 /* Allocate 2 extra elements on both sides, so in single
138 * EKIN_WORK_BUFFER_SIZE*DIM*DIM*sizeof(real) = 72/144 bytes
139 * buffer on both sides to avoid cache pollution.
141 snew(ekind
->ekin_work_alloc
[thread
], ekind
->ngtc
+2*EKIN_WORK_BUFFER_SIZE
);
142 ekind
->ekin_work
[thread
] = ekind
->ekin_work_alloc
[thread
] + EKIN_WORK_BUFFER_SIZE
;
143 /* Nasty hack so we can have the per-thread accumulation
144 * variable for dekindl in the same thread-local cache lines
145 * as the per-thread accumulation tensors for ekin[fh],
146 * because they are accumulated in the same loop. */
147 ekind
->dekindl_work
[thread
] = &(ekind
->ekin_work
[thread
][ekind
->ngtc
][0][0]);
148 #undef EKIN_WORK_BUFFER_SIZE
150 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
153 ekind
->ngacc
= opts
->ngacc
;
154 snew(ekind
->grpstat
, opts
->ngacc
);
155 init_grpstat(mtop
, opts
->ngacc
, ekind
->grpstat
);
158 void accumulate_u(const t_commrec
*cr
, const t_grpopts
*opts
, gmx_ekindata_t
*ekind
)
160 /* This routine will only be called when it's necessary */
166 for (g
= 0; (g
< opts
->ngacc
); g
++)
168 add_binr(rb
, DIM
, ekind
->grpstat
[g
].u
);
172 for (g
= 0; (g
< opts
->ngacc
); g
++)
174 extract_binr(rb
, DIM
*g
, DIM
, ekind
->grpstat
[g
].u
);
179 void update_ekindata(int start
, int homenr
, gmx_ekindata_t
*ekind
,
180 const t_grpopts
*opts
, const rvec v
[], const t_mdatoms
*md
, real lambda
)
185 /* calculate mean velocities at whole timestep */
186 for (g
= 0; (g
< opts
->ngtc
); g
++)
188 ekind
->tcstat
[g
].T
= 0;
193 for (g
= 0; (g
< opts
->ngacc
); g
++)
195 clear_rvec(ekind
->grpstat
[g
].u
);
199 for (n
= start
; (n
< start
+homenr
); n
++)
205 for (d
= 0; (d
< DIM
); d
++)
207 mv
= md
->massT
[n
]*v
[n
][d
];
208 ekind
->grpstat
[g
].u
[d
] += mv
;
212 for (g
= 0; (g
< opts
->ngacc
); g
++)
214 for (d
= 0; (d
< DIM
); d
++)
216 ekind
->grpstat
[g
].u
[d
] /=
217 (1-lambda
)*ekind
->grpstat
[g
].mA
+ lambda
*ekind
->grpstat
[g
].mB
;
223 real
sum_ekin(const t_grpopts
*opts
, gmx_ekindata_t
*ekind
, real
*dekindlambda
,
224 gmx_bool bEkinAveVel
, gmx_bool bScaleEkin
)
228 t_grp_tcstat
*tcstat
;
237 clear_mat(ekind
->ekin
);
239 for (i
= 0; (i
< ngtc
); i
++)
243 tcstat
= &ekind
->tcstat
[i
];
244 /* Sometimes a group does not have degrees of freedom, e.g.
245 * when it consists of shells and virtual sites, then we just
246 * set the temperatue to 0 and also neglect the kinetic
247 * energy, which should be zero anyway.
256 /* in this case, kinetic energy is from the current velocities already */
257 msmul(tcstat
->ekinf
, tcstat
->ekinscalef_nhc
, tcstat
->ekinf
);
262 /* Calculate the full step Ekin as the average of the half steps */
263 for (j
= 0; (j
< DIM
); j
++)
265 for (m
= 0; (m
< DIM
); m
++)
267 tcstat
->ekinf
[j
][m
] =
268 0.5*(tcstat
->ekinh
[j
][m
]*tcstat
->ekinscaleh_nhc
+ tcstat
->ekinh_old
[j
][m
]);
272 m_add(tcstat
->ekinf
, ekind
->ekin
, ekind
->ekin
);
274 tcstat
->Th
= calc_temp(trace(tcstat
->ekinh
), nd
);
275 tcstat
->T
= calc_temp(trace(tcstat
->ekinf
), nd
);
277 /* after the scaling factors have been multiplied in, we can remove them */
280 tcstat
->ekinscalef_nhc
= 1.0;
284 tcstat
->ekinscaleh_nhc
= 1.0;
303 *dekindlambda
= ekind
->dekindl
;
307 *dekindlambda
= 0.5*(ekind
->dekindl
+ ekind
->dekindl_old
);