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 #ifndef GMX_MDLIB_MD_SUPPORT_H
38 #define GMX_MDLIB_MD_SUPPORT_H
40 #include "gromacs/mdlib/vcm.h"
41 #include "gromacs/timing/wallcycle.h"
43 struct gmx_ekindata_t;
44 struct gmx_enerdata_t;
45 struct gmx_global_stat;
46 struct gmx_multisim_t;
47 struct gmx_signalling_t;
63 class SimulationSignaller;
66 /* Define a number of flags to better control the information
67 * passed to compute_globals in md.c and global_stat.
70 /* we are computing the kinetic energy from average velocities */
71 #define CGLO_EKINAVEVEL (1u << 2u)
72 /* we are removing the center of mass momenta */
73 #define CGLO_STOPCM (1u << 3u)
74 /* bGStat is defined in do_md */
75 #define CGLO_GSTAT (1u << 4u)
76 /* Sum the energy terms in global computation */
77 #define CGLO_ENERGY (1u << 6u)
78 /* Sum the kinetic energy terms in global computation */
79 #define CGLO_TEMPERATURE (1u << 7u)
80 /* Sum the kinetic energy terms in global computation */
81 #define CGLO_PRESSURE (1u << 8u)
82 /* Sum the constraint term in global computation */
83 #define CGLO_CONSTRAINT (1u << 9u)
84 /* Reading ekin from the trajectory */
85 #define CGLO_READEKIN (1u << 10u)
86 /* we need to reset the ekin rescaling factor here */
87 #define CGLO_SCALEEKIN (1u << 11u)
88 /* After a new DD partitioning, we need to set a flag to schedule
89 * global reduction of the total number of bonded interactions that
90 * will be computed, to check none are missing. */
91 #define CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS (1u << 12u)
94 /*! \brief Return the number of steps that will take place between
95 * intra-simulation communications, given the constraints of the
97 int computeGlobalCommunicationPeriod(const gmx::MDLogger& mdlog, t_inputrec* ir, const t_commrec* cr);
99 /*! \brief Return true if the \p value is equal across the set of multi-simulations
101 * \todo This duplicates some of check_multi_int. Consolidate. */
102 bool multisim_int_all_are_equal(const gmx_multisim_t* ms, int64_t value);
104 void rerun_parallel_comm(t_commrec* cr, t_trxframe* fr, gmx_bool* bLastStep);
106 //! \brief Allocate and initialize node-local state entries
107 void set_state_entries(t_state* state, const t_inputrec* ir);
109 /* Set the lambda values in the global state from a frame read with rerun */
110 void setCurrentLambdasRerun(int64_t step,
111 const t_lambda* fepvals,
112 const t_trxframe* rerun_fr,
114 t_state* globalState);
116 /* Set the lambda values at each step of mdrun when they change */
117 void setCurrentLambdasLocal(int64_t step,
118 const t_lambda* fepvals,
120 gmx::ArrayRef<real> lambda,
121 int currentFEPState);
123 int multisim_min(const gmx_multisim_t* ms, int nmin, int n);
124 /* Set an appropriate value for n across the whole multi-simulation */
127 /* Compute global variables during integration
129 * Coordinates x are needed for kinetic energy calculation with cosine accelation
130 * and for COM removal with rotational and acceleration correction modes.
131 * Velocities v are needed for kinetic energy calculation and for COM removal.
133 void compute_globals(gmx_global_stat* gstat,
135 const t_inputrec* ir,
137 gmx_ekindata_t* ekind,
142 const t_mdatoms* mdatoms,
145 gmx_wallcycle_t wcycle,
146 gmx_enerdata_t* enerd,
152 gmx::Constraints* constr,
153 gmx::SimulationSignaller* signalCoordinator,
154 const matrix lastbox,
155 int* totalNumberOfBondedInteractions,
156 gmx_bool* bSumEkinhOld,