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.
41 * This file contains the definition of the microstate of the simulated system
43 * History of observables that needs to be checkpointed should be stored
44 * in ObservablesHistory.
45 * The state of the mdrun machinery that needs to be checkpointed is also
51 * \ingroup module_mdtypes
54 #ifndef GMX_MDTYPES_STATE_H
55 #define GMX_MDTYPES_STATE_H
61 #include "gromacs/gpu_utils/hostallocator.h"
62 #include "gromacs/math/paddedvector.h"
63 #include "gromacs/math/vectypes.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/utility/arrayref.h"
66 #include "gromacs/utility/basedefinitions.h"
67 #include "gromacs/utility/real.h"
76 //! Convenience alias for until all is moved in the gmx namespace
78 using PaddedHostVector = gmx::PaddedHostVector<T>;
81 * The t_state struct should contain all the (possibly) non-static
82 * information required to define the state of the system.
83 * Currently the random seeds for SD and BD are missing.
86 /* \brief Enum for all entries in \p t_state
88 * These enums are used in flags as (1<<est...).
89 * The order of these enums should not be changed,
90 * since that affects the checkpoint (.cpt) file format.
94 estBOX, estBOX_REL, estBOXV, estPRES_PREV, estNH_XI, estTHERM_INT,
95 estX, estV, estSDX_NOTSUPPORTED, estCGP,
96 estLD_RNG_NOTSUPPORTED, estLD_RNGI_NOTSUPPORTED,
97 estDISRE_INITF, estDISRE_RM3TAV,
98 estORIRE_INITF, estORIRE_DTAV,
99 estSVIR_PREV, estNH_VXI, estVETA, estVOL0, estNHPRES_XI, estNHPRES_VXI, estFVIR_PREV,
100 estFEPSTATE, estMC_RNG_NOTSUPPORTED, estMC_RNGI_NOTSUPPORTED,
101 estBAROS_INT, estPULLCOMPREVSTEP,
105 //! \brief The names of the state entries, defined in src/gmxlib/checkpoint.c
106 extern const char *est_names[estNR];
108 /*! \libinternal \brief History information for NMR distance and orientation restraints
110 * Often this is only used for reporting observables, and thus should not
111 * actually be part of the microstate. But with time-dependent restraining
112 * they are actually part of the (non-Markovian) microstate.
113 * \todo Rename this with a more descriptive name.
120 real disre_initf; //!< The scaling factor for initializing the time av.
121 int ndisrepairs; //!< The number of distance restraints
122 real *disre_rm3tav; //!< The r^-3 time averaged pair distances
123 real orire_initf; //!< The scaling factor for initializing the time av.
124 int norire_Dtav; //!< The number of matrix element in dtav (npair*5)
125 real *orire_Dtav; //!< The time averaged orientation tensors
128 /*! \libinternal \brief Struct used for checkpointing only
130 * This struct would not be required with unlimited precision.
131 * But because of limited precision, the COM motion removal implementation
132 * can cause the kinetic energy in the MD loop to differ by a few bits from
133 * the kinetic energy one would determine from state.v.
140 gmx_bool bUpToDate; //!< Test if all data is up to date
141 int ekin_n; //!< The number of tensors
142 tensor *ekinh; //!< Half step Ekin, size \p ekin_n
143 tensor *ekinf; //!< Full step Ekin, size \p ekin_n
144 tensor *ekinh_old; //!< Half step Ekin of the previous step, size \p ekin_n
145 tensor ekin_total; //!< Total kinetic energy
146 std::vector<double> ekinscalef_nhc; //!< Nose-Hoover Ekin scaling factors for full step Ekin
147 std::vector<double> ekinscaleh_nhc; //!< Nose-Hoover Ekin scaling factors for half step Ekin
148 std::vector<double> vscale_nhc; //!< Nose-Hoover velocity scaling factors
149 real dekindl; //!< dEkin/dlambda, with free-energy
150 real mvcos; //!< Cosine(z) component of the momentum, for viscosity calculations
153 /*! \brief Free-energy sampling history struct
155 * \todo Split out into microstate and observables history.
157 typedef struct df_history_t
159 int nlambda; //!< total number of lambda states - for history
161 gmx_bool bEquil; //!< Have we reached equilibration
162 int *n_at_lam; //!< number of points observed at each lambda
163 real *wl_histo; //!< histogram for WL flatness determination
164 real wl_delta; //!< current wang-landau delta
166 real *sum_weights; //!< weights of the states
167 real *sum_dg; //!< free energies of the states -- not actually used for weighting, but informational
168 real *sum_minvar; //!< corrections to weights for minimum variance
169 real *sum_variance; //!< variances of the states
171 real **accum_p; //!< accumulated bennett weights for n+1
172 real **accum_m; //!< accumulated bennett weights for n-1
173 real **accum_p2; //!< accumulated squared bennett weights for n+1
174 real **accum_m2; //!< accumulated squared bennett weights for n-1
176 real **Tij; //!< transition matrix
177 real **Tij_empirical; //!< Empirical transition matrix
182 /*! \brief The microstate of the system
184 * The global state will contain complete data for all used entries.
185 * The local state with domain decomposition will have partial entries
186 * for which \p stateEntryIsAtomProperty() is true. Some entries that
187 * are used in the global state might not be present in the local state.
188 * \todo Move pure observables history to ObservablesHistory.
196 int natoms; //!< Number of atoms, local + non-local; this is the size of \p x, \p v and \p cg_p, when used
197 int ngtc; //!< The number of temperature coupling groups
198 int nnhpres; //!< The NH-chain length for the MTTK barostat
199 int nhchainlength; //!< The NH-chain length for temperature coupling
200 int flags; //!< Set of bit-flags telling which entries are present, see enum at the top of the file
201 int fep_state; //!< indicates which of the alchemical states we are in
202 std::array<real, efptNR> lambda; //!< Free-energy lambda vector
203 matrix box; //!< Matrix of box vectors
204 matrix box_rel; //!< Relative box vectors to preserve box shape
205 matrix boxv; //!< Box velocities for Parrinello-Rahman P-coupling
206 matrix pres_prev; //!< Pressure of the previous step for pcoupl
207 matrix svir_prev; //!< Shake virial for previous step for pcoupl
208 matrix fvir_prev; //!< Force virial of the previous step for pcoupl
209 std::vector<double> nosehoover_xi; //!< Nose-Hoover coordinates (ngtc)
210 std::vector<double> nosehoover_vxi; //!< Nose-Hoover velocities (ngtc)
211 std::vector<double> nhpres_xi; //!< Pressure Nose-Hoover coordinates
212 std::vector<double> nhpres_vxi; //!< Pressure Nose-Hoover velocities
213 std::vector<double> therm_integral; //!< Work exterted N-H/V-rescale T-coupling (ngtc)
214 double baros_integral; //!< For Berendsen P-coupling conserved quantity
215 real veta; //!< Trotter based isotropic P-coupling
216 real vol0; //!< Initial volume,required for computing MTTK conserved quantity
217 PaddedHostVector<gmx::RVec> x; //!< The coordinates (natoms)
218 PaddedVector<gmx::RVec> v; //!< The velocities (natoms)
219 PaddedVector<gmx::RVec> cg_p; //!< p vector for conjugate gradient minimization
221 ekinstate_t ekinstate; //!< The state of the kinetic energy
223 /* History for special algorithms, should be moved to a history struct */
224 history_t hist; //!< Time history for restraints
225 df_history_t *dfhist; //!< Free-energy history for free energy analysis
226 std::shared_ptr<gmx::AwhHistory> awhHistory; //!< Accelerated weight histogram history
228 int ddp_count; //!< The DD partitioning count for this state
229 int ddp_count_cg_gl; //!< The DD partitioning count for index_gl
230 std::vector<int> cg_gl; //!< The global cg number of the local cgs
232 std::vector<double> pull_com_prev_step; //!< The COM of the previous step of each pull group
236 /* We don't document the structs below, as they don't belong here.
237 * TODO: Move the next two structs out of state.h.
242 std::vector<double> Qinv; /* inverse mass of thermostat -- computed from inputs, but a good place to store */
243 std::vector<double> QPinv; /* inverse mass of thermostat for barostat -- computed from inputs, but a good place to store */
244 double Winv; /* Pressure mass inverse -- computed, not input, but a good place to store. Need to make a matrix later */
245 tensor Winvm; /* inverse pressure mass tensor, computed */
261 //! Resizes the T- and P-coupling state variables
262 void init_gtc_state(t_state *state, int ngtc, int nnhpres, int nhchainlength);
264 //! Change the number of atoms represented by this state, allocating memory as needed.
265 void state_change_natoms(t_state *state, int natoms);
267 //! Allocates memory for free-energy history
268 void init_dfhist_state(t_state *state, int dfhistNumLambda);
270 /*! \brief Compares two states, write the differences to stdout */
271 void comp_state(const t_state *st1, const t_state *st2, gmx_bool bRMSD, real ftol, real abstol);
273 /*! \brief Allocates an rvec pointer and copy the contents of v to it */
274 rvec *makeRvecArray(gmx::ArrayRef<const gmx::RVec> v,
277 /*! \brief Determine the relative box components
279 * Set box_rel e.g. used in mdrun state, used to preserve the box shape
280 * \param[in] ir Input record
281 * \param[inout] state State
283 void set_box_rel(const t_inputrec *ir, t_state *state);
285 /*! \brief Make sure the relative box shape remains the same
287 * This function ensures that the relative box dimensions are
288 * preserved, which otherwise might diffuse away due to rounding
289 * errors in pressure coupling or the deform option.
291 * \param[in] ir Input record
292 * \param[in] box_rel Relative box dimensions
293 * \param[inout] box The corrected actual box dimensions
295 void preserve_box_shape(const t_inputrec *ir, matrix box_rel, matrix box);
297 /*! \brief Returns an arrayRef to the positions in \p state when \p state!=null
299 * When \p state=nullptr, returns an empty arrayRef.
301 * \note The size returned is the number of atoms, without padding.
303 * \param[in] state The state, can be nullptr
305 static inline gmx::ArrayRef<const gmx::RVec>
306 positionsFromStatePointer(const t_state *state)
310 return gmx::makeConstArrayRef(state->x).subArray(0, state->natoms);
314 return gmx::EmptyArrayRef();