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,2021, 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 #ifndef GMX_MDLIB_VCM_H
39 #define GMX_MDLIB_VCM_H
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/mdtypes/md_enums.h"
47 #include "gromacs/utility/basedefinitions.h"
48 #include "gromacs/utility/real.h"
50 struct SimulationGroups;
78 //! Size of group arrays
80 //! Stride for thread data
82 //! One of the enums above
83 ComRemovalAlgorithm mode = ComRemovalAlgorithm::Linear;
84 //! The number of dimensions for corr.
86 //! The time step for COMM removal
88 //! Number of degrees of freedom
89 std::vector<real> group_ndf;
91 std::vector<real> group_mass;
92 //! Linear momentum per group
93 std::vector<gmx::RVec> group_p;
94 //! Linear velocity per group
95 std::vector<gmx::RVec> group_v;
96 //! Center of mass per group
97 std::vector<gmx::RVec> group_x;
98 //! Angular momentum per group
99 std::vector<gmx::RVec> group_j;
100 //! Angular velocity (omega)
101 std::vector<gmx::RVec> group_w;
102 //! Moment of inertia per group
103 tensor* group_i = nullptr;
104 //! These two are copies to pointers in
105 std::vector<char*> group_name;
106 //! Tells whether dimensions are frozen per freeze group
107 ivec* nFreeze = nullptr;
108 //! Temporary data per thread and group
109 std::vector<t_vcm_thread> thread_vcm;
111 //! Tell whether the integrator conserves momentum
112 bool integratorConservesMomentum = false;
114 t_vcm(const SimulationGroups& groups, const t_inputrec& ir);
118 /* print COM removal info to log */
119 void reportComRemovalInfo(FILE* fp, const t_vcm& vcm);
122 /* Do a per group center of mass things */
123 void calc_vcm_grp(const t_mdatoms& md,
124 gmx::ArrayRef<const gmx::RVec> x,
125 gmx::ArrayRef<const gmx::RVec> v,
128 /* Set the COM velocity to zero and potentially correct the COM position.
130 * Processes the kinetic energy reduced over MPI before removing COM motion.
131 * With mode linear, nullptr can be passed for x.
132 * With acceleration correction nullptr should be passed for x at initialization
133 * and a pointer to the coordinates at normal MD steps.
134 * When fplog != nullptr, a warning is printed to fplog with high COM velocity.
136 void process_and_stopcm_grp(FILE* fplog,
138 const t_mdatoms& mdatoms,
139 gmx::ArrayRef<gmx::RVec> x,
140 gmx::ArrayRef<gmx::RVec> v);