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.
38 /* This file is completely threadsafe - keep it that way! */
47 #include "gromacs/math/paddedvector.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/math/veccompare.h"
50 #include "gromacs/mdtypes/awh_history.h"
51 #include "gromacs/mdtypes/df_history.h"
52 #include "gromacs/mdtypes/inputrec.h"
53 #include "gromacs/mdtypes/md_enums.h"
54 #include "gromacs/mdtypes/pull_params.h"
55 #include "gromacs/mdtypes/swaphistory.h"
56 #include "gromacs/pbcutil/boxutilities.h"
57 #include "gromacs/pbcutil/pbc.h"
58 #include "gromacs/utility/compare.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/smalloc.h"
61 #include "gromacs/utility/stringutil.h"
63 #include "checkpointdata.h"
65 /* The source code in this file should be thread-safe.
66 Please keep it that way. */
68 history_t::history_t() :
71 disre_rm3tav(nullptr),
78 ekinstate_t::ekinstate_t() :
89 clear_mat(ekin_total
);
95 * \brief Enum describing the contents ekinstate_t writes to modular checkpoint
97 * When changing the checkpoint content, add a new element just above Count, and adjust the
98 * checkpoint functionality.
100 enum class CheckpointVersion
102 Base
, //!< First version of modular checkpointing
103 Count
//!< Number of entries. Add new versions right above this!
105 constexpr auto c_currentVersion
= CheckpointVersion(int(CheckpointVersion::Count
) - 1);
108 template<gmx::CheckpointDataOperation operation
>
109 void ekinstate_t::doCheckpoint(gmx::CheckpointData
<operation
> checkpointData
)
111 gmx::checkpointVersion(&checkpointData
, "ekinstate_t version", c_currentVersion
);
113 checkpointData
.scalar("bUpToDate", &bUpToDate
);
118 auto numOfTensors
= ekin_n
;
119 checkpointData
.scalar("ekin_n", &numOfTensors
);
120 if (operation
== gmx::CheckpointDataOperation::Read
)
122 // If this isn't matching, we haven't allocated the right amount of data
123 GMX_RELEASE_ASSERT(numOfTensors
== ekin_n
,
124 "ekinstate_t checkpoint reading: Tensor size mismatch.");
126 for (int idx
= 0; idx
< numOfTensors
; ++idx
)
128 checkpointData
.tensor(gmx::formatString("ekinh %d", idx
), ekinh
[idx
]);
129 checkpointData
.tensor(gmx::formatString("ekinf %d", idx
), ekinf
[idx
]);
130 checkpointData
.tensor(gmx::formatString("ekinh_old %d", idx
), ekinh_old
[idx
]);
132 checkpointData
.arrayRef("ekinscalef_nhc", gmx::makeCheckpointArrayRef
<operation
>(ekinscalef_nhc
));
133 checkpointData
.arrayRef("ekinscaleh_nhc", gmx::makeCheckpointArrayRef
<operation
>(ekinscaleh_nhc
));
134 checkpointData
.arrayRef("vscale_nhc", gmx::makeCheckpointArrayRef
<operation
>(vscale_nhc
));
135 checkpointData
.scalar("dekindl", &dekindl
);
136 checkpointData
.scalar("mvcos", &mvcos
);
138 if (operation
== gmx::CheckpointDataOperation::Read
)
140 hasReadEkinState
= true;
144 // Explicit template instantiation
145 template void ekinstate_t::doCheckpoint(gmx::CheckpointData
<gmx::CheckpointDataOperation::Read
> checkpointData
);
146 template void ekinstate_t::doCheckpoint(gmx::CheckpointData
<gmx::CheckpointDataOperation::Write
> checkpointData
);
148 void init_gtc_state(t_state
* state
, int ngtc
, int nnhpres
, int nhchainlength
)
151 state
->nnhpres
= nnhpres
;
152 state
->nhchainlength
= nhchainlength
;
153 state
->nosehoover_xi
.resize(state
->nhchainlength
* state
->ngtc
, 0);
154 state
->nosehoover_vxi
.resize(state
->nhchainlength
* state
->ngtc
, 0);
155 state
->therm_integral
.resize(state
->ngtc
, 0);
156 state
->baros_integral
= 0.0;
157 state
->nhpres_xi
.resize(state
->nhchainlength
* nnhpres
, 0);
158 state
->nhpres_vxi
.resize(state
->nhchainlength
* nnhpres
, 0);
162 /* Checkpoint code relies on this function having no effect if
163 state->natoms is > 0 and passed as natoms. */
164 void state_change_natoms(t_state
* state
, int natoms
)
166 state
->natoms
= natoms
;
168 /* We need padding, since we might use SIMD access, but the
169 * containers here all ensure that. */
170 if (state
->flags
& (1 << estX
))
172 state
->x
.resizeWithPadding(natoms
);
174 if (state
->flags
& (1 << estV
))
176 state
->v
.resizeWithPadding(natoms
);
178 if (state
->flags
& (1 << estCGP
))
180 state
->cg_p
.resizeWithPadding(natoms
);
184 void init_dfhist_state(t_state
* state
, int dfhistNumLambda
)
186 if (dfhistNumLambda
> 0)
188 snew(state
->dfhist
, 1);
189 init_df_history(state
->dfhist
, dfhistNumLambda
);
193 state
->dfhist
= nullptr;
197 void comp_state(const t_state
* st1
, const t_state
* st2
, gmx_bool bRMSD
, real ftol
, real abstol
)
201 fprintf(stdout
, "comparing flags\n");
202 cmp_int(stdout
, "flags", -1, st1
->flags
, st2
->flags
);
203 fprintf(stdout
, "comparing box\n");
204 cmp_rvecs(stdout
, "box", DIM
, st1
->box
, st2
->box
, FALSE
, ftol
, abstol
);
205 fprintf(stdout
, "comparing box_rel\n");
206 cmp_rvecs(stdout
, "box_rel", DIM
, st1
->box_rel
, st2
->box_rel
, FALSE
, ftol
, abstol
);
207 fprintf(stdout
, "comparing boxv\n");
208 cmp_rvecs(stdout
, "boxv", DIM
, st1
->boxv
, st2
->boxv
, FALSE
, ftol
, abstol
);
209 if (st1
->flags
& (1 << estSVIR_PREV
))
211 fprintf(stdout
, "comparing shake vir_prev\n");
212 cmp_rvecs(stdout
, "svir_prev", DIM
, st1
->svir_prev
, st2
->svir_prev
, FALSE
, ftol
, abstol
);
214 if (st1
->flags
& (1 << estFVIR_PREV
))
216 fprintf(stdout
, "comparing force vir_prev\n");
217 cmp_rvecs(stdout
, "fvir_prev", DIM
, st1
->fvir_prev
, st2
->fvir_prev
, FALSE
, ftol
, abstol
);
219 if (st1
->flags
& (1 << estPRES_PREV
))
221 fprintf(stdout
, "comparing prev_pres\n");
222 cmp_rvecs(stdout
, "pres_prev", DIM
, st1
->pres_prev
, st2
->pres_prev
, FALSE
, ftol
, abstol
);
224 cmp_int(stdout
, "ngtc", -1, st1
->ngtc
, st2
->ngtc
);
225 cmp_int(stdout
, "nhchainlength", -1, st1
->nhchainlength
, st2
->nhchainlength
);
226 if (st1
->ngtc
== st2
->ngtc
&& st1
->nhchainlength
== st2
->nhchainlength
)
228 for (i
= 0; i
< st1
->ngtc
; i
++)
230 nc
= i
* st1
->nhchainlength
;
231 for (j
= 0; j
< nc
; j
++)
233 cmp_real(stdout
, "nosehoover_xi", i
, st1
->nosehoover_xi
[nc
+ j
],
234 st2
->nosehoover_xi
[nc
+ j
], ftol
, abstol
);
238 cmp_int(stdout
, "nnhpres", -1, st1
->nnhpres
, st2
->nnhpres
);
239 if (st1
->nnhpres
== st2
->nnhpres
&& st1
->nhchainlength
== st2
->nhchainlength
)
241 for (i
= 0; i
< st1
->nnhpres
; i
++)
243 nc
= i
* st1
->nhchainlength
;
244 for (j
= 0; j
< nc
; j
++)
246 cmp_real(stdout
, "nosehoover_xi", i
, st1
->nhpres_xi
[nc
+ j
], st2
->nhpres_xi
[nc
+ j
],
252 cmp_int(stdout
, "natoms", -1, st1
->natoms
, st2
->natoms
);
253 if (st1
->natoms
== st2
->natoms
)
255 if ((st1
->flags
& (1 << estX
)) && (st2
->flags
& (1 << estX
)))
257 fprintf(stdout
, "comparing x\n");
258 cmp_rvecs(stdout
, "x", st1
->natoms
, st1
->x
.rvec_array(), st2
->x
.rvec_array(), bRMSD
,
261 if ((st1
->flags
& (1 << estV
)) && (st2
->flags
& (1 << estV
)))
263 fprintf(stdout
, "comparing v\n");
264 cmp_rvecs(stdout
, "v", st1
->natoms
, st1
->v
.rvec_array(), st2
->v
.rvec_array(), bRMSD
,
270 rvec
* makeRvecArray(gmx::ArrayRef
<const gmx::RVec
> v
, gmx::index n
)
272 GMX_ASSERT(v
.ssize() >= n
, "We can't copy more elements than the vector size");
278 const rvec
* vPtr
= as_rvec_array(v
.data());
279 for (gmx::index i
= 0; i
< n
; i
++)
281 copy_rvec(vPtr
[i
], dest
[i
]);
308 // It would be nicer to initialize these with {} or {{0}} in the
309 // above initialization list, but uncrustify doesn't understand
311 // TODO Fix this if we switch to clang-format some time.
316 clear_mat(pres_prev
);
317 clear_mat(svir_prev
);
318 clear_mat(fvir_prev
);
321 void set_box_rel(const t_inputrec
* ir
, t_state
* state
)
323 /* Make sure the box obeys the restrictions before we fix the ratios */
324 correct_box(nullptr, 0, state
->box
);
326 clear_mat(state
->box_rel
);
328 if (inputrecPreserveShape(ir
))
330 const int ndim
= ir
->epct
== epctSEMIISOTROPIC
? 2 : 3;
331 do_box_rel(ndim
, ir
->deform
, state
->box_rel
, state
->box
, true);
335 void preserve_box_shape(const t_inputrec
* ir
, matrix box_rel
, matrix box
)
337 if (inputrecPreserveShape(ir
))
339 const int ndim
= ir
->epct
== epctSEMIISOTROPIC
? 2 : 3;
340 do_box_rel(ndim
, ir
->deform
, box_rel
, box
, false);
344 void printLambdaStateToLog(FILE* fplog
, const gmx::ArrayRef
<real
> lambda
, const bool isInitialOutput
)
346 if (fplog
!= nullptr)
348 fprintf(fplog
, "%s vector of lambda components:[ ", isInitialOutput
? "Initial" : "Current");
349 for (const auto& l
: lambda
)
351 fprintf(fplog
, "%10.4f ", l
);
353 fprintf(fplog
, "]\n%s", isInitialOutput
? "" : "\n");
357 void initialize_lambdas(FILE* fplog
, const t_inputrec
& ir
, bool isMaster
, int* fep_state
, gmx::ArrayRef
<real
> lambda
)
359 /* TODO: Clean up initialization of fep_state and lambda in
360 t_state. This function works, but could probably use a logic
361 rewrite to keep all the different types of efep straight. */
363 if ((ir
.efep
== efepNO
) && (!ir
.bSimTemp
))
368 const t_lambda
* fep
= ir
.fepvals
;
371 *fep_state
= fep
->init_fep_state
; /* this might overwrite the checkpoint
372 if checkpoint is set -- a kludge is in for now
376 for (int i
= 0; i
< efptNR
; i
++)
379 /* overwrite lambda state with init_lambda for now for backwards compatibility */
380 if (fep
->init_lambda
>= 0) /* if it's -1, it was never initialized */
382 thisLambda
= fep
->init_lambda
;
386 thisLambda
= fep
->all_lambda
[i
][fep
->init_fep_state
];
390 lambda
[i
] = thisLambda
;
395 /* need to rescale control temperatures to match current state */
396 for (int i
= 0; i
< ir
.opts
.ngtc
; i
++)
398 if (ir
.opts
.ref_t
[i
] > 0)
400 ir
.opts
.ref_t
[i
] = ir
.simtempvals
->temperatures
[fep
->init_fep_state
];
405 /* Send to the log the information on the current lambdas */
406 const bool isInitialOutput
= true;
407 printLambdaStateToLog(fplog
, lambda
, isInitialOutput
);