6c791cdba20e4fb26a442f12117d20422eb86295
[alexxy/gromacs.git] / src / gromacs / mdtypes / state.cpp
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 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
40
41 #include "state.h"
42
43 #include <cstring>
44
45 #include <algorithm>
46
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"
62
63 #include "checkpointdata.h"
64
65 /* The source code in this file should be thread-safe.
66       Please keep it that way. */
67
68 history_t::history_t() :
69     disre_initf(0),
70     ndisrepairs(0),
71     disre_rm3tav(nullptr),
72     orire_initf(0),
73     norire_Dtav(0),
74     orire_Dtav(nullptr)
75 {
76 }
77
78 ekinstate_t::ekinstate_t() :
79     bUpToDate(FALSE),
80     ekin_n(0),
81     ekinh(nullptr),
82     ekinf(nullptr),
83     ekinh_old(nullptr),
84     ekin_total(),
85
86     dekindl(0),
87     mvcos(0)
88 {
89     clear_mat(ekin_total);
90 }
91
92 namespace
93 {
94 /*!
95  * \brief Enum describing the contents ekinstate_t writes to modular checkpoint
96  *
97  * When changing the checkpoint content, add a new element just above Count, and adjust the
98  * checkpoint functionality.
99  */
100 enum class CheckpointVersion
101 {
102     Base, //!< First version of modular checkpointing
103     Count //!< Number of entries. Add new versions right above this!
104 };
105 constexpr auto c_currentVersion = CheckpointVersion(int(CheckpointVersion::Count) - 1);
106 } // namespace
107
108 template<gmx::CheckpointDataOperation operation>
109 void ekinstate_t::doCheckpoint(gmx::CheckpointData<operation> checkpointData)
110 {
111     auto version = c_currentVersion;
112     checkpointData.enumScalar("version", &version);
113     if (version != c_currentVersion)
114     {
115         throw gmx::FileIOError("ekinstate_t checkpoint version mismatch.");
116     }
117
118     checkpointData.scalar("bUpToDate", &bUpToDate);
119     if (!bUpToDate)
120     {
121         return;
122     }
123     auto numOfTensors = ekin_n;
124     checkpointData.scalar("ekin_n", &numOfTensors);
125     if (operation == gmx::CheckpointDataOperation::Read)
126     {
127         // If this isn't matching, we haven't allocated the right amount of data
128         GMX_RELEASE_ASSERT(numOfTensors == ekin_n,
129                            "ekinstate_t checkpoint reading: Tensor size mismatch.");
130     }
131     for (int idx = 0; idx < numOfTensors; ++idx)
132     {
133         checkpointData.tensor(gmx::formatString("ekinh %d", idx), ekinh[idx]);
134         checkpointData.tensor(gmx::formatString("ekinf %d", idx), ekinf[idx]);
135         checkpointData.tensor(gmx::formatString("ekinh_old %d", idx), ekinh_old[idx]);
136     }
137     checkpointData.arrayRef("ekinscalef_nhc", gmx::makeCheckpointArrayRef<operation>(ekinscalef_nhc));
138     checkpointData.arrayRef("ekinscaleh_nhc", gmx::makeCheckpointArrayRef<operation>(ekinscaleh_nhc));
139     checkpointData.arrayRef("vscale_nhc", gmx::makeCheckpointArrayRef<operation>(vscale_nhc));
140     checkpointData.scalar("dekindl", &dekindl);
141     checkpointData.scalar("mvcos", &mvcos);
142
143     if (operation == gmx::CheckpointDataOperation::Read)
144     {
145         hasReadEkinState = true;
146     }
147 }
148
149 // Explicit template instantiation
150 template void ekinstate_t::doCheckpoint(gmx::CheckpointData<gmx::CheckpointDataOperation::Read> checkpointData);
151 template void ekinstate_t::doCheckpoint(gmx::CheckpointData<gmx::CheckpointDataOperation::Write> checkpointData);
152
153 void init_gtc_state(t_state* state, int ngtc, int nnhpres, int nhchainlength)
154 {
155     state->ngtc          = ngtc;
156     state->nnhpres       = nnhpres;
157     state->nhchainlength = nhchainlength;
158     state->nosehoover_xi.resize(state->nhchainlength * state->ngtc, 0);
159     state->nosehoover_vxi.resize(state->nhchainlength * state->ngtc, 0);
160     state->therm_integral.resize(state->ngtc, 0);
161     state->baros_integral = 0.0;
162     state->nhpres_xi.resize(state->nhchainlength * nnhpres, 0);
163     state->nhpres_vxi.resize(state->nhchainlength * nnhpres, 0);
164 }
165
166
167 /* Checkpoint code relies on this function having no effect if
168    state->natoms is > 0 and passed as natoms. */
169 void state_change_natoms(t_state* state, int natoms)
170 {
171     state->natoms = natoms;
172
173     /* We need padding, since we might use SIMD access, but the
174      * containers here all ensure that. */
175     if (state->flags & (1 << estX))
176     {
177         state->x.resizeWithPadding(natoms);
178     }
179     if (state->flags & (1 << estV))
180     {
181         state->v.resizeWithPadding(natoms);
182     }
183     if (state->flags & (1 << estCGP))
184     {
185         state->cg_p.resizeWithPadding(natoms);
186     }
187 }
188
189 void init_dfhist_state(t_state* state, int dfhistNumLambda)
190 {
191     if (dfhistNumLambda > 0)
192     {
193         snew(state->dfhist, 1);
194         init_df_history(state->dfhist, dfhistNumLambda);
195     }
196     else
197     {
198         state->dfhist = nullptr;
199     }
200 }
201
202 void comp_state(const t_state* st1, const t_state* st2, gmx_bool bRMSD, real ftol, real abstol)
203 {
204     int i, j, nc;
205
206     fprintf(stdout, "comparing flags\n");
207     cmp_int(stdout, "flags", -1, st1->flags, st2->flags);
208     fprintf(stdout, "comparing box\n");
209     cmp_rvecs(stdout, "box", DIM, st1->box, st2->box, FALSE, ftol, abstol);
210     fprintf(stdout, "comparing box_rel\n");
211     cmp_rvecs(stdout, "box_rel", DIM, st1->box_rel, st2->box_rel, FALSE, ftol, abstol);
212     fprintf(stdout, "comparing boxv\n");
213     cmp_rvecs(stdout, "boxv", DIM, st1->boxv, st2->boxv, FALSE, ftol, abstol);
214     if (st1->flags & (1 << estSVIR_PREV))
215     {
216         fprintf(stdout, "comparing shake vir_prev\n");
217         cmp_rvecs(stdout, "svir_prev", DIM, st1->svir_prev, st2->svir_prev, FALSE, ftol, abstol);
218     }
219     if (st1->flags & (1 << estFVIR_PREV))
220     {
221         fprintf(stdout, "comparing force vir_prev\n");
222         cmp_rvecs(stdout, "fvir_prev", DIM, st1->fvir_prev, st2->fvir_prev, FALSE, ftol, abstol);
223     }
224     if (st1->flags & (1 << estPRES_PREV))
225     {
226         fprintf(stdout, "comparing prev_pres\n");
227         cmp_rvecs(stdout, "pres_prev", DIM, st1->pres_prev, st2->pres_prev, FALSE, ftol, abstol);
228     }
229     cmp_int(stdout, "ngtc", -1, st1->ngtc, st2->ngtc);
230     cmp_int(stdout, "nhchainlength", -1, st1->nhchainlength, st2->nhchainlength);
231     if (st1->ngtc == st2->ngtc && st1->nhchainlength == st2->nhchainlength)
232     {
233         for (i = 0; i < st1->ngtc; i++)
234         {
235             nc = i * st1->nhchainlength;
236             for (j = 0; j < nc; j++)
237             {
238                 cmp_real(stdout, "nosehoover_xi", i, st1->nosehoover_xi[nc + j],
239                          st2->nosehoover_xi[nc + j], ftol, abstol);
240             }
241         }
242     }
243     cmp_int(stdout, "nnhpres", -1, st1->nnhpres, st2->nnhpres);
244     if (st1->nnhpres == st2->nnhpres && st1->nhchainlength == st2->nhchainlength)
245     {
246         for (i = 0; i < st1->nnhpres; i++)
247         {
248             nc = i * st1->nhchainlength;
249             for (j = 0; j < nc; j++)
250             {
251                 cmp_real(stdout, "nosehoover_xi", i, st1->nhpres_xi[nc + j], st2->nhpres_xi[nc + j],
252                          ftol, abstol);
253             }
254         }
255     }
256
257     cmp_int(stdout, "natoms", -1, st1->natoms, st2->natoms);
258     if (st1->natoms == st2->natoms)
259     {
260         if ((st1->flags & (1 << estX)) && (st2->flags & (1 << estX)))
261         {
262             fprintf(stdout, "comparing x\n");
263             cmp_rvecs(stdout, "x", st1->natoms, st1->x.rvec_array(), st2->x.rvec_array(), bRMSD,
264                       ftol, abstol);
265         }
266         if ((st1->flags & (1 << estV)) && (st2->flags & (1 << estV)))
267         {
268             fprintf(stdout, "comparing v\n");
269             cmp_rvecs(stdout, "v", st1->natoms, st1->v.rvec_array(), st2->v.rvec_array(), bRMSD,
270                       ftol, abstol);
271         }
272     }
273 }
274
275 rvec* makeRvecArray(gmx::ArrayRef<const gmx::RVec> v, gmx::index n)
276 {
277     GMX_ASSERT(v.ssize() >= n, "We can't copy more elements than the vector size");
278
279     rvec* dest;
280
281     snew(dest, n);
282
283     const rvec* vPtr = as_rvec_array(v.data());
284     for (gmx::index i = 0; i < n; i++)
285     {
286         copy_rvec(vPtr[i], dest[i]);
287     }
288
289     return dest;
290 }
291
292 t_state::t_state() :
293     natoms(0),
294     ngtc(0),
295     nnhpres(0),
296     nhchainlength(0),
297     flags(0),
298     fep_state(0),
299     lambda(),
300
301     baros_integral(0),
302     veta(0),
303     vol0(0),
304
305     ekinstate(),
306     hist(),
307     dfhist(nullptr),
308     awhHistory(nullptr),
309     ddp_count(0),
310     ddp_count_cg_gl(0)
311
312 {
313     // It would be nicer to initialize these with {} or {{0}} in the
314     // above initialization list, but uncrustify doesn't understand
315     // that.
316     // TODO Fix this if we switch to clang-format some time.
317     lambda = { { 0 } };
318     clear_mat(box);
319     clear_mat(box_rel);
320     clear_mat(boxv);
321     clear_mat(pres_prev);
322     clear_mat(svir_prev);
323     clear_mat(fvir_prev);
324 }
325
326 void set_box_rel(const t_inputrec* ir, t_state* state)
327 {
328     /* Make sure the box obeys the restrictions before we fix the ratios */
329     correct_box(nullptr, 0, state->box);
330
331     clear_mat(state->box_rel);
332
333     if (inputrecPreserveShape(ir))
334     {
335         const int ndim = ir->epct == epctSEMIISOTROPIC ? 2 : 3;
336         do_box_rel(ndim, ir->deform, state->box_rel, state->box, true);
337     }
338 }
339
340 void preserve_box_shape(const t_inputrec* ir, matrix box_rel, matrix box)
341 {
342     if (inputrecPreserveShape(ir))
343     {
344         const int ndim = ir->epct == epctSEMIISOTROPIC ? 2 : 3;
345         do_box_rel(ndim, ir->deform, box_rel, box, false);
346     }
347 }
348
349 void printLambdaStateToLog(FILE* fplog, const gmx::ArrayRef<real> lambda, const bool isInitialOutput)
350 {
351     if (fplog != nullptr)
352     {
353         fprintf(fplog, "%s vector of lambda components:[ ", isInitialOutput ? "Initial" : "Current");
354         for (const auto& l : lambda)
355         {
356             fprintf(fplog, "%10.4f ", l);
357         }
358         fprintf(fplog, "]\n%s", isInitialOutput ? "" : "\n");
359     }
360 }
361
362 void initialize_lambdas(FILE*               fplog,
363                         const t_inputrec&   ir,
364                         bool                isMaster,
365                         int*                fep_state,
366                         gmx::ArrayRef<real> lambda,
367                         double*             lam0)
368 {
369     /* TODO: Clean up initialization of fep_state and lambda in
370        t_state.  This function works, but could probably use a logic
371        rewrite to keep all the different types of efep straight. */
372
373     if ((ir.efep == efepNO) && (!ir.bSimTemp))
374     {
375         return;
376     }
377
378     const t_lambda* fep = ir.fepvals;
379     if (isMaster)
380     {
381         *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
382                                              if checkpoint is set -- a kludge is in for now
383                                              to prevent this.*/
384     }
385
386     for (int i = 0; i < efptNR; i++)
387     {
388         double thisLambda;
389         /* overwrite lambda state with init_lambda for now for backwards compatibility */
390         if (fep->init_lambda >= 0) /* if it's -1, it was never initialized */
391         {
392             thisLambda = fep->init_lambda;
393         }
394         else
395         {
396             thisLambda = fep->all_lambda[i][fep->init_fep_state];
397         }
398         if (isMaster)
399         {
400             lambda[i] = thisLambda;
401         }
402         if (lam0 != nullptr)
403         {
404             lam0[i] = thisLambda;
405         }
406     }
407     if (ir.bSimTemp)
408     {
409         /* need to rescale control temperatures to match current state */
410         for (int i = 0; i < ir.opts.ngtc; i++)
411         {
412             if (ir.opts.ref_t[i] > 0)
413             {
414                 ir.opts.ref_t[i] = ir.simtempvals->temperatures[fep->init_fep_state];
415             }
416         }
417     }
418
419     /* Send to the log the information on the current lambdas */
420     const bool isInitialOutput = true;
421     printLambdaStateToLog(fplog, lambda, isInitialOutput);
422 }