Simplified uniform GPU selection in CMake
[gromacs.git] / src / gromacs / mdtypes / state.h
blob016711df7e3137df3fd11bd18d7be0e83e997353
1 /*
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.
39 /*! \file
41 * \brief
42 * This file contains the definition of the microstate of the simulated system
44 * History of observables that needs to be checkpointed should be stored
45 * in ObservablesHistory.
46 * The state of the mdrun machinery that needs to be checkpointed is also
47 * stored elsewhere.
49 * \author Berk Hess
51 * \inlibraryapi
52 * \ingroup module_mdtypes
55 #ifndef GMX_MDTYPES_STATE_H
56 #define GMX_MDTYPES_STATE_H
58 #include <array>
59 #include <memory>
60 #include <vector>
62 #include "gromacs/gpu_utils/hostallocator.h"
63 #include "gromacs/math/paddedvector.h"
64 #include "gromacs/math/vectypes.h"
65 #include "gromacs/mdtypes/md_enums.h"
66 #include "gromacs/utility/arrayref.h"
67 #include "gromacs/utility/basedefinitions.h"
68 #include "gromacs/utility/real.h"
70 struct t_inputrec;
72 namespace gmx
74 struct AwhHistory;
77 //! Convenience alias for until all is moved in the gmx namespace
78 template<class T>
79 using PaddedHostVector = gmx::PaddedHostVector<T>;
82 * The t_state struct should contain all the (possibly) non-static
83 * information required to define the state of the system.
84 * Currently the random seeds for SD and BD are missing.
87 /* \brief Enum for all entries in \p t_state
89 * These enums are used in flags as (1<<est...).
90 * The order of these enums should not be changed,
91 * since that affects the checkpoint (.cpt) file format.
93 enum
95 estLAMBDA,
96 estBOX,
97 estBOX_REL,
98 estBOXV,
99 estPRES_PREV,
100 estNH_XI,
101 estTHERM_INT,
102 estX,
103 estV,
104 estSDX_NOTSUPPORTED,
105 estCGP,
106 estLD_RNG_NOTSUPPORTED,
107 estLD_RNGI_NOTSUPPORTED,
108 estDISRE_INITF,
109 estDISRE_RM3TAV,
110 estORIRE_INITF,
111 estORIRE_DTAV,
112 estSVIR_PREV,
113 estNH_VXI,
114 estVETA,
115 estVOL0,
116 estNHPRES_XI,
117 estNHPRES_VXI,
118 estFVIR_PREV,
119 estFEPSTATE,
120 estMC_RNG_NOTSUPPORTED,
121 estMC_RNGI_NOTSUPPORTED,
122 estBAROS_INT,
123 estPULLCOMPREVSTEP,
124 estNR
127 //! \brief The names of the state entries, defined in src/gmxlib/checkpoint.c
128 extern const char* est_names[estNR];
130 /*! \libinternal \brief History information for NMR distance and orientation restraints
132 * Often this is only used for reporting observables, and thus should not
133 * actually be part of the microstate. But with time-dependent restraining
134 * they are actually part of the (non-Markovian) microstate.
135 * \todo Rename this with a more descriptive name.
137 class history_t
139 public:
140 history_t();
142 real disre_initf; //!< The scaling factor for initializing the time av.
143 int ndisrepairs; //!< The number of distance restraints
144 real* disre_rm3tav; //!< The r^-3 time averaged pair distances
145 real orire_initf; //!< The scaling factor for initializing the time av.
146 int norire_Dtav; //!< The number of matrix element in dtav (npair*5)
147 real* orire_Dtav; //!< The time averaged orientation tensors
150 /*! \libinternal \brief Struct used for checkpointing only
152 * This struct would not be required with unlimited precision.
153 * But because of limited precision, the COM motion removal implementation
154 * can cause the kinetic energy in the MD loop to differ by a few bits from
155 * the kinetic energy one would determine from state.v.
157 class ekinstate_t
159 public:
160 ekinstate_t();
162 gmx_bool bUpToDate; //!< Test if all data is up to date
163 int ekin_n; //!< The number of tensors
164 tensor* ekinh; //!< Half step Ekin, size \p ekin_n
165 tensor* ekinf; //!< Full step Ekin, size \p ekin_n
166 tensor* ekinh_old; //!< Half step Ekin of the previous step, size \p ekin_n
167 tensor ekin_total; //!< Total kinetic energy
168 std::vector<double> ekinscalef_nhc; //!< Nose-Hoover Ekin scaling factors for full step Ekin
169 std::vector<double> ekinscaleh_nhc; //!< Nose-Hoover Ekin scaling factors for half step Ekin
170 std::vector<double> vscale_nhc; //!< Nose-Hoover velocity scaling factors
171 real dekindl; //!< dEkin/dlambda, with free-energy
172 real mvcos; //!< Cosine(z) component of the momentum, for viscosity calculations
173 /*! \brief Whether KE terms have been read from the checkpoint.
175 * Only used for managing whether the call to compute_globals
176 * before we enter the MD loop should compute these quantities
177 * fresh, or not. */
178 bool hasReadEkinState;
181 /*! \brief Free-energy sampling history struct
183 * \todo Split out into microstate and observables history.
185 typedef struct df_history_t
187 int nlambda; //!< total number of lambda states - for history
189 gmx_bool bEquil; //!< Have we reached equilibration
190 int* n_at_lam; //!< number of points observed at each lambda
191 real* wl_histo; //!< histogram for WL flatness determination
192 real wl_delta; //!< current wang-landau delta
194 real* sum_weights; //!< weights of the states
195 real* sum_dg; //!< free energies of the states -- not actually used for weighting, but informational
196 real* sum_minvar; //!< corrections to weights for minimum variance
197 real* sum_variance; //!< variances of the states
199 real** accum_p; //!< accumulated bennett weights for n+1
200 real** accum_m; //!< accumulated bennett weights for n-1
201 real** accum_p2; //!< accumulated squared bennett weights for n+1
202 real** accum_m2; //!< accumulated squared bennett weights for n-1
204 real** Tij; //!< transition matrix
205 real** Tij_empirical; //!< Empirical transition matrix
207 } df_history_t;
210 /*! \brief The microstate of the system
212 * The global state will contain complete data for all used entries.
213 * The local state with domain decomposition will have partial entries
214 * for which \p stateEntryIsAtomProperty() is true. Some entries that
215 * are used in the global state might not be present in the local state.
216 * \todo Move pure observables history to ObservablesHistory.
218 class t_state
220 public:
221 t_state();
223 // All things public
224 int natoms; //!< Number of atoms, local + non-local; this is the size of \p x, \p v and \p cg_p, when used
225 int ngtc; //!< The number of temperature coupling groups
226 int nnhpres; //!< The NH-chain length for the MTTK barostat
227 int nhchainlength; //!< The NH-chain length for temperature coupling
228 int flags; //!< Set of bit-flags telling which entries are present, see enum at the top of the file
229 int fep_state; //!< indicates which of the alchemical states we are in
230 std::array<real, efptNR> lambda; //!< Free-energy lambda vector
231 matrix box; //!< Matrix of box vectors
232 matrix box_rel; //!< Relative box vectors to preserve box shape
233 matrix boxv; //!< Box velocities for Parrinello-Rahman P-coupling
234 matrix pres_prev; //!< Pressure of the previous step for pcoupl
235 matrix svir_prev; //!< Shake virial for previous step for pcoupl
236 matrix fvir_prev; //!< Force virial of the previous step for pcoupl
237 std::vector<double> nosehoover_xi; //!< Nose-Hoover coordinates (ngtc)
238 std::vector<double> nosehoover_vxi; //!< Nose-Hoover velocities (ngtc)
239 std::vector<double> nhpres_xi; //!< Pressure Nose-Hoover coordinates
240 std::vector<double> nhpres_vxi; //!< Pressure Nose-Hoover velocities
241 std::vector<double> therm_integral; //!< Work exterted N-H/V-rescale T-coupling (ngtc)
242 double baros_integral; //!< For Berendsen P-coupling conserved quantity
243 real veta; //!< Trotter based isotropic P-coupling
244 real vol0; //!< Initial volume,required for computing MTTK conserved quantity
245 PaddedHostVector<gmx::RVec> x; //!< The coordinates (natoms)
246 PaddedHostVector<gmx::RVec> v; //!< The velocities (natoms)
247 PaddedHostVector<gmx::RVec> cg_p; //!< p vector for conjugate gradient minimization
249 ekinstate_t ekinstate; //!< The state of the kinetic energy
251 /* History for special algorithms, should be moved to a history struct */
252 history_t hist; //!< Time history for restraints
253 df_history_t* dfhist; //!< Free-energy history for free energy analysis
254 std::shared_ptr<gmx::AwhHistory> awhHistory; //!< Accelerated weight histogram history
256 int ddp_count; //!< The DD partitioning count for this state
257 int ddp_count_cg_gl; //!< The DD partitioning count for index_gl
258 std::vector<int> cg_gl; //!< The global cg number of the local cgs
260 std::vector<double> pull_com_prev_step; //!< The COM of the previous step of each pull group
263 #ifndef DOXYGEN
264 /* We don't document the structs below, as they don't belong here.
265 * TODO: Move the next two structs out of state.h.
268 struct t_extmass
270 std::vector<double> Qinv; /* inverse mass of thermostat -- computed from inputs, but a good place to store */
271 std::vector<double> QPinv; /* inverse mass of thermostat for barostat -- computed from inputs, but a good place to store */
272 double Winv; /* Pressure mass inverse -- computed, not input, but a good place to store. Need to make a matrix later */
273 tensor Winvm; /* inverse pressure mass tensor, computed */
277 typedef struct
279 real veta;
280 double rscale;
281 double vscale;
282 double rvscale;
283 double alpha;
284 double* vscale_nhc;
285 } t_vetavars;
287 #endif // DOXYGEN
289 //! Resizes the T- and P-coupling state variables
290 void init_gtc_state(t_state* state, int ngtc, int nnhpres, int nhchainlength);
292 //! Change the number of atoms represented by this state, allocating memory as needed.
293 void state_change_natoms(t_state* state, int natoms);
295 //! Allocates memory for free-energy history
296 void init_dfhist_state(t_state* state, int dfhistNumLambda);
298 /*! \brief Compares two states, write the differences to stdout */
299 void comp_state(const t_state* st1, const t_state* st2, gmx_bool bRMSD, real ftol, real abstol);
301 /*! \brief Allocates an rvec pointer and copy the contents of v to it */
302 rvec* makeRvecArray(gmx::ArrayRef<const gmx::RVec> v, gmx::index n);
304 /*! \brief Determine the relative box components
306 * Set box_rel e.g. used in mdrun state, used to preserve the box shape
307 * \param[in] ir Input record
308 * \param[inout] state State
310 void set_box_rel(const t_inputrec* ir, t_state* state);
312 /*! \brief Make sure the relative box shape remains the same
314 * This function ensures that the relative box dimensions are
315 * preserved, which otherwise might diffuse away due to rounding
316 * errors in pressure coupling or the deform option.
318 * \param[in] ir Input record
319 * \param[in] box_rel Relative box dimensions
320 * \param[inout] box The corrected actual box dimensions
322 void preserve_box_shape(const t_inputrec* ir, matrix box_rel, matrix box);
324 /*! \brief Returns an arrayRef to the positions in \p state when \p state!=null
326 * When \p state=nullptr, returns an empty arrayRef.
328 * \note The size returned is the number of atoms, without padding.
330 * \param[in] state The state, can be nullptr
332 static inline gmx::ArrayRef<const gmx::RVec> positionsFromStatePointer(const t_state* state)
334 if (state)
336 return gmx::makeConstArrayRef(state->x).subArray(0, state->natoms);
338 else
340 return {};
344 /*! \brief Fills fep_state, lambda, and lam0 if needed
346 * If FEP or simulated tempering is in use:
348 * fills non-null lam0 with the initial lambda values, and
349 * on master rank fills fep_state and lambda.
351 * Reports the initial lambda state to the log file. */
352 void initialize_lambdas(FILE* fplog,
353 const t_inputrec& ir,
354 bool isMaster,
355 int* fep_state,
356 gmx::ArrayRef<real> lambda,
357 double* lam0);
359 #endif