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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 /* This file is completely threadsafe - keep it that way! */
43 #include "gromacs/gmxlib/network.h"
44 #include "gromacs/math/functions.h"
45 #include "gromacs/math/invertmatrix.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/math/vecdump.h"
48 #include "gromacs/mdlib/gmx_omp_nthreads.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/md_enums.h"
51 #include "gromacs/pbcutil/pbc.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/gmxassert.h"
55 #include "gromacs/utility/gmxomp.h"
56 #include "gromacs/utility/smalloc.h"
58 t_vcm::t_vcm(const SimulationGroups& groups, const t_inputrec& ir) :
59 integratorConservesMomentum(!EI_RANDOM(ir.eI))
61 mode = (ir.nstcomm > 0) ? ir.comm_mode : ecmNO;
63 timeStep = ir.nstcomm * ir.delta_t;
65 if (mode == ecmANGULAR && ndim < 3)
67 gmx_fatal(FARGS, "Can not have angular comm removal with pbc=%s",
68 c_pbcTypeNames[ir.pbcType].c_str());
73 nr = groups.groups[SimulationAtomGroupType::MassCenterVelocityRemoval].size();
74 /* Allocate one extra for a possible rest group */
76 /* We need vcm->nr+1 elements per thread, but to avoid cache
77 * invalidation we add 2 elements to get a 152 byte separation.
80 if (mode == ecmANGULAR)
89 group_name.resize(size);
92 group_mass.resize(size);
93 group_ndf.resize(size);
94 for (int g = 0; (g < nr); g++)
96 group_ndf[g] = ir.opts.nrdf[g];
98 *groups.groupNames[groups.groups[SimulationAtomGroupType::MassCenterVelocityRemoval][g]];
101 thread_vcm.resize(gmx_omp_nthreads_get(emntDefault) * stride);
104 nFreeze = ir.opts.nFreeze;
109 if (mode == ecmANGULAR)
115 void reportComRemovalInfo(FILE* fp, const t_vcm& vcm)
118 /* Copy pointer to group names and print it. */
119 if (fp && vcm.mode != ecmNO)
121 fprintf(fp, "Center of mass motion removal mode is %s\n", ECOM(vcm.mode));
123 "We have the following groups for center of"
124 " mass motion removal:\n");
126 for (int g = 0; (g < vcm.nr); g++)
129 fprintf(fp, "%3d: %s\n", g, vcm.group_name[g]);
135 static void update_tensor(const rvec x, real m0, tensor I)
139 /* Compute inertia tensor contribution due to this atom */
140 xy = x[XX] * x[YY] * m0;
141 xz = x[XX] * x[ZZ] * m0;
142 yz = x[YY] * x[ZZ] * m0;
143 I[XX][XX] += x[XX] * x[XX] * m0;
144 I[YY][YY] += x[YY] * x[YY] * m0;
145 I[ZZ][ZZ] += x[ZZ] * x[ZZ] * m0;
154 /* Center of mass code for groups */
155 void calc_vcm_grp(const t_mdatoms& md, const rvec x[], const rvec v[], t_vcm* vcm)
157 int nthreads = gmx_omp_nthreads_get(emntDefault);
158 if (vcm->mode != ecmNO)
160 #pragma omp parallel num_threads(nthreads)
162 int t = gmx_omp_get_thread_num();
163 for (int g = 0; g < vcm->size; g++)
165 /* Reset linear momentum */
166 t_vcm_thread* vcm_t = &vcm->thread_vcm[t * vcm->stride + g];
168 clear_rvec(vcm_t->p);
169 if (vcm->mode == ecmANGULAR)
171 /* Reset angular momentum */
172 clear_rvec(vcm_t->j);
173 clear_rvec(vcm_t->x);
178 #pragma omp for schedule(static)
179 for (int i = 0; i < md.homenr; i++)
182 real m0 = md.massT[i];
187 t_vcm_thread* vcm_t = &vcm->thread_vcm[t * vcm->stride + g];
188 /* Calculate linear momentum */
191 for (m = 0; (m < DIM); m++)
193 vcm_t->p[m] += m0 * v[i][m];
196 if (vcm->mode == ecmANGULAR)
198 /* Calculate angular momentum */
200 cprod(x[i], v[i], j0);
202 for (m = 0; (m < DIM); m++)
204 vcm_t->j[m] += m0 * j0[m];
205 vcm_t->x[m] += m0 * x[i][m];
207 /* Update inertia tensor */
208 update_tensor(x[i], m0, vcm_t->i);
212 for (int g = 0; g < vcm->size; g++)
214 /* Reset linear momentum */
215 vcm->group_mass[g] = 0;
216 clear_rvec(vcm->group_p[g]);
217 if (vcm->mode == ecmANGULAR)
219 /* Reset angular momentum */
220 clear_rvec(vcm->group_j[g]);
221 clear_rvec(vcm->group_x[g]);
222 clear_rvec(vcm->group_w[g]);
223 clear_mat(vcm->group_i[g]);
226 for (int t = 0; t < nthreads; t++)
228 t_vcm_thread* vcm_t = &vcm->thread_vcm[t * vcm->stride + g];
229 vcm->group_mass[g] += vcm_t->mass;
230 rvec_inc(vcm->group_p[g], vcm_t->p);
231 if (vcm->mode == ecmANGULAR)
233 rvec_inc(vcm->group_j[g], vcm_t->j);
234 rvec_inc(vcm->group_x[g], vcm_t->x);
235 m_add(vcm_t->i, vcm->group_i[g], vcm->group_i[g]);
242 /*! \brief Remove the COM motion velocity from the velocities
244 * \note This routine should be called from within an OpenMP parallel region.
246 * \tparam numDimensions Correct dimensions 0 to \p numDimensions-1
247 * \param[in] mdatoms The atom property and group information
248 * \param[in,out] v The velocities to correct
249 * \param[in] vcm VCM data
251 template<int numDimensions>
252 static void doStopComMotionLinear(const t_mdatoms& mdatoms, rvec* v, const t_vcm& vcm)
254 const int homenr = mdatoms.homenr;
255 const unsigned short* group_id = mdatoms.cVCM;
257 if (mdatoms.cFREEZE != nullptr)
259 GMX_RELEASE_ASSERT(vcm.nFreeze != nullptr, "Need freeze dimension info with freeze groups");
261 #pragma omp for schedule(static)
262 for (int i = 0; i < homenr; i++)
264 unsigned short vcmGroup = (group_id == nullptr ? 0 : group_id[i]);
265 unsigned short freezeGroup = mdatoms.cFREEZE[i];
266 for (int d = 0; d < numDimensions; d++)
268 if (vcm.nFreeze[freezeGroup][d] == 0)
270 v[i][d] -= vcm.group_v[vcmGroup][d];
275 else if (group_id == nullptr)
277 #pragma omp for schedule(static)
278 for (int i = 0; i < homenr; i++)
280 for (int d = 0; d < numDimensions; d++)
282 v[i][d] -= vcm.group_v[0][d];
288 #pragma omp for schedule(static)
289 for (int i = 0; i < homenr; i++)
291 const int g = group_id[i];
292 for (int d = 0; d < numDimensions; d++)
294 v[i][d] -= vcm.group_v[g][d];
300 /*! \brief Remove the COM motion velocity from the velocities, correct the coordinates assuming constant acceleration
302 * \note This routine should be called from within an OpenMP parallel region.
304 * \tparam numDimensions Correct dimensions 0 to \p numDimensions-1
305 * \param[in] homenr The number of atoms to correct
306 * \param[in] group_id List of VCM group ids, when nullptr is passed all atoms are assumed to be in group 0
307 * \param[in,out] x The coordinates to correct
308 * \param[in,out] v The velocities to correct
309 * \param[in] vcm VCM data
311 template<int numDimensions>
312 static void doStopComMotionAccelerationCorrection(int homenr,
313 const unsigned short* group_id,
314 rvec* gmx_restrict x,
315 rvec* gmx_restrict v,
318 const real xCorrectionFactor = 0.5 * vcm.timeStep;
320 if (group_id == nullptr)
322 #pragma omp for schedule(static)
323 for (int i = 0; i < homenr; i++)
325 for (int d = 0; d < numDimensions; d++)
327 x[i][d] -= vcm.group_v[0][d] * xCorrectionFactor;
328 v[i][d] -= vcm.group_v[0][d];
334 #pragma omp for schedule(static)
335 for (int i = 0; i < homenr; i++)
337 const int g = group_id[i];
338 for (int d = 0; d < numDimensions; d++)
340 x[i][d] -= vcm.group_v[g][d] * xCorrectionFactor;
341 v[i][d] -= vcm.group_v[g][d];
347 static void do_stopcm_grp(const t_mdatoms& mdatoms, rvec x[], rvec v[], const t_vcm& vcm)
349 if (vcm.mode != ecmNO)
351 const int homenr = mdatoms.homenr;
352 const unsigned short* group_id = mdatoms.cVCM;
354 int gmx_unused nth = gmx_omp_nthreads_get(emntDefault);
355 #pragma omp parallel num_threads(nth)
357 if (vcm.mode == ecmLINEAR || vcm.mode == ecmANGULAR
358 || (vcm.mode == ecmLINEAR_ACCELERATION_CORRECTION && x == nullptr))
360 /* Subtract linear momentum for v */
363 case 1: doStopComMotionLinear<1>(mdatoms, v, vcm); break;
364 case 2: doStopComMotionLinear<2>(mdatoms, v, vcm); break;
365 case 3: doStopComMotionLinear<3>(mdatoms, v, vcm); break;
370 GMX_ASSERT(vcm.mode == ecmLINEAR_ACCELERATION_CORRECTION,
371 "When the mode is not linear or angular, it should be acceleration "
373 /* Subtract linear momentum for v and x*/
377 doStopComMotionAccelerationCorrection<1>(homenr, group_id, x, v, vcm);
380 doStopComMotionAccelerationCorrection<2>(homenr, group_id, x, v, vcm);
383 doStopComMotionAccelerationCorrection<3>(homenr, group_id, x, v, vcm);
387 if (vcm.mode == ecmANGULAR)
389 /* Subtract angular momentum */
390 GMX_ASSERT(x, "Need x to compute angular momentum correction");
393 #pragma omp for schedule(static)
394 for (int i = 0; i < homenr; i++)
400 /* Compute the correction to the velocity for each atom */
402 rvec_sub(x[i], vcm.group_x[g], dx);
403 cprod(vcm.group_w[g], dx, dv);
411 static void get_minv(tensor A, tensor B)
417 tmp[XX][XX] = A[YY][YY] + A[ZZ][ZZ];
418 tmp[YY][XX] = -A[XX][YY];
419 tmp[ZZ][XX] = -A[XX][ZZ];
420 tmp[XX][YY] = -A[XX][YY];
421 tmp[YY][YY] = A[XX][XX] + A[ZZ][ZZ];
422 tmp[ZZ][YY] = -A[YY][ZZ];
423 tmp[XX][ZZ] = -A[XX][ZZ];
424 tmp[YY][ZZ] = -A[YY][ZZ];
425 tmp[ZZ][ZZ] = A[XX][XX] + A[YY][YY];
427 /* This is a hack to prevent very large determinants */
428 rfac = (tmp[XX][XX] + tmp[YY][YY] + tmp[ZZ][ZZ]) / 3;
431 gmx_fatal(FARGS, "Can not stop center of mass: maybe 2dimensional system");
434 for (m = 0; (m < DIM); m++)
436 for (n = 0; (n < DIM); n++)
441 gmx::invertMatrix(tmp, B);
442 for (m = 0; (m < DIM); m++)
444 for (n = 0; (n < DIM); n++)
451 /* Processes VCM after reduction over ranks and prints warning with high VMC and fp != nullptr */
452 static void process_and_check_cm_grp(FILE* fp, t_vcm* vcm, real Temp_Max)
455 real ekcm, ekrot, tm, tm_1, Temp_cm;
459 /* First analyse the total results */
460 if (vcm->mode != ecmNO)
462 for (g = 0; (g < vcm->nr); g++)
464 tm = vcm->group_mass[g];
468 svmul(tm_1, vcm->group_p[g], vcm->group_v[g]);
470 /* Else it's zero anyway! */
472 if (vcm->mode == ecmANGULAR)
474 for (g = 0; (g < vcm->nr); g++)
476 tm = vcm->group_mass[g];
481 /* Compute center of mass for this group */
482 for (m = 0; (m < DIM); m++)
484 vcm->group_x[g][m] *= tm_1;
487 /* Subtract the center of mass contribution to the
490 cprod(vcm->group_x[g], vcm->group_v[g], jcm);
491 for (m = 0; (m < DIM); m++)
493 vcm->group_j[g][m] -= tm * jcm[m];
496 /* Subtract the center of mass contribution from the inertia
497 * tensor (this is not as trivial as it seems, but due to
498 * some cancellation we can still do it, even in parallel).
501 update_tensor(vcm->group_x[g], tm, Icm);
502 m_sub(vcm->group_i[g], Icm, vcm->group_i[g]);
504 /* Compute angular velocity, using matrix operation
509 get_minv(vcm->group_i[g], Icm);
510 mvmul(Icm, vcm->group_j[g], vcm->group_w[g]);
512 /* Else it's zero anyway! */
516 for (g = 0; (g < vcm->nr); g++)
519 if (vcm->group_mass[g] != 0 && vcm->group_ndf[g] > 0)
521 for (m = 0; m < vcm->ndim; m++)
523 ekcm += gmx::square(vcm->group_v[g][m]);
525 ekcm *= 0.5 * vcm->group_mass[g];
526 Temp_cm = 2 * ekcm / vcm->group_ndf[g];
528 if ((Temp_cm > Temp_Max) && fp)
530 fprintf(fp, "Large VCM(group %s): %12.5f, %12.5f, %12.5f, Temp-cm: %12.5e\n",
531 vcm->group_name[g], vcm->group_v[g][XX], vcm->group_v[g][YY],
532 vcm->group_v[g][ZZ], Temp_cm);
535 if (vcm->mode == ecmANGULAR)
537 ekrot = 0.5 * iprod(vcm->group_j[g], vcm->group_w[g]);
538 // TODO: Change absolute energy comparison to relative
539 if ((ekrot > 1) && fp && vcm->integratorConservesMomentum)
541 /* if we have an integrator that may not conserve momenta, skip */
542 tm = vcm->group_mass[g];
543 fprintf(fp, "Group %s with mass %12.5e, Ekrot %12.5e Det(I) = %12.5e\n",
544 vcm->group_name[g], tm, ekrot, det(vcm->group_i[g]));
545 fprintf(fp, " COM: %12.5f %12.5f %12.5f\n", vcm->group_x[g][XX],
546 vcm->group_x[g][YY], vcm->group_x[g][ZZ]);
547 fprintf(fp, " P: %12.5f %12.5f %12.5f\n", vcm->group_p[g][XX],
548 vcm->group_p[g][YY], vcm->group_p[g][ZZ]);
549 fprintf(fp, " V: %12.5f %12.5f %12.5f\n", vcm->group_v[g][XX],
550 vcm->group_v[g][YY], vcm->group_v[g][ZZ]);
551 fprintf(fp, " J: %12.5f %12.5f %12.5f\n", vcm->group_j[g][XX],
552 vcm->group_j[g][YY], vcm->group_j[g][ZZ]);
553 fprintf(fp, " w: %12.5f %12.5f %12.5f\n", vcm->group_w[g][XX],
554 vcm->group_w[g][YY], vcm->group_w[g][ZZ]);
555 pr_rvecs(fp, 0, "Inertia tensor", vcm->group_i[g], DIM);
562 void process_and_stopcm_grp(FILE* fplog, t_vcm* vcm, const t_mdatoms& mdatoms, rvec x[], rvec v[])
564 if (vcm->mode != ecmNO)
566 // TODO: Replace fixed temperature of 1 by a system value
567 process_and_check_cm_grp(fplog, vcm, 1);
569 do_stopcm_grp(mdatoms, x, v, *vcm);