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