801174e5e55f312f5f50c5bb489a3336350ef257
[alexxy/gromacs.git] / src / gromacs / fileio / tpxio.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,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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 /* This file is completely threadsafe - keep it that way! */
40
41 #include "tpxio.h"
42
43 #include <cstdio>
44 #include <cstdlib>
45 #include <cstring>
46
47 #include <algorithm>
48 #include <memory>
49 #include <vector>
50
51 #include "gromacs/fileio/filetypes.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/gmxfio_xdr.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdtypes/awh_history.h"
57 #include "gromacs/mdtypes/awh_params.h"
58 #include "gromacs/mdtypes/inputrec.h"
59 #include "gromacs/mdtypes/md_enums.h"
60 #include "gromacs/mdtypes/pull_params.h"
61 #include "gromacs/mdtypes/state.h"
62 #include "gromacs/pbcutil/boxutilities.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/topology/block.h"
65 #include "gromacs/topology/ifunc.h"
66 #include "gromacs/topology/mtop_util.h"
67 #include "gromacs/topology/symtab.h"
68 #include "gromacs/topology/topology.h"
69 #include "gromacs/utility/arraysize.h"
70 #include "gromacs/utility/baseversion.h"
71 #include "gromacs/utility/cstringutil.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/futil.h"
74 #include "gromacs/utility/gmxassert.h"
75 #include "gromacs/utility/inmemoryserializer.h"
76 #include "gromacs/utility/keyvaluetreebuilder.h"
77 #include "gromacs/utility/keyvaluetreeserializer.h"
78 #include "gromacs/utility/smalloc.h"
79 #include "gromacs/utility/snprintf.h"
80 #include "gromacs/utility/txtdump.h"
81
82 #define TPX_TAG_RELEASE "release"
83
84 /*! \brief Tag string for the file format written to run input files
85  * written by this version of the code.
86  *
87  * Change this if you want to change the run input format in a feature
88  * branch. This ensures that there will not be different run input
89  * formats around which cannot be distinguished, while not causing
90  * problems rebasing the feature branch onto upstream changes. When
91  * merging with mainstream GROMACS, set this tag string back to
92  * TPX_TAG_RELEASE, and instead add an element to tpxv.
93  */
94 static const char* tpx_tag = TPX_TAG_RELEASE;
95
96 /*! \brief Enum of values that describe the contents of a tpr file
97  * whose format matches a version number
98  *
99  * The enum helps the code be more self-documenting and ensure merges
100  * do not silently resolve when two patches make the same bump. When
101  * adding new functionality, add a new element just above tpxv_Count
102  * in this enumeration, and write code below that does the right thing
103  * according to the value of file_version.
104  */
105 enum tpxv
106 {
107     tpxv_ComputationalElectrophysiology =
108             96, /**< support for ion/water position swaps (computational electrophysiology) */
109     tpxv_Use64BitRandomSeed, /**< change ld_seed from int to int64_t */
110     tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
111     tpxv_InteractiveMolecularDynamics, /**< interactive molecular dynamics (IMD) */
112     tpxv_RemoveObsoleteParameters1,    /**< remove optimize_fft, dihre_fc, nstcheckpoint */
113     tpxv_PullCoordTypeGeom,            /**< add pull type and geometry per group and flat-bottom */
114     tpxv_PullGeomDirRel,               /**< add pull geometry direction-relative */
115     tpxv_IntermolecularBondeds, /**< permit inter-molecular bonded interactions in the topology */
116     tpxv_CompElWithSwapLayerOffset, /**< added parameters for improved CompEl setups */
117     tpxv_CompElPolyatomicIonsAndMultipleIonTypes, /**< CompEl now can handle polyatomic ions and more than two types of ions */
118     tpxv_RemoveAdress,                            /**< removed support for AdResS */
119     tpxv_PullCoordNGroup,               /**< add ngroup to pull coord */
120     tpxv_RemoveTwinRange,               /**< removed support for twin-range interactions */
121     tpxv_ReplacePullPrintCOM12,         /**< Replaced print-com-1, 2 with pull-print-com */
122     tpxv_PullExternalPotential,         /**< Added pull type external potential */
123     tpxv_GenericParamsForElectricField, /**< Introduced KeyValueTree and moved electric field parameters */
124     tpxv_AcceleratedWeightHistogram, /**< sampling with accelerated weight histogram method (AWH) */
125     tpxv_RemoveImplicitSolvation,    /**< removed support for implicit solvation */
126     tpxv_PullPrevStepCOMAsReference, /**< Enabled using the COM of the pull group of the last frame as reference for PBC */
127     tpxv_MimicQMMM,   /**< Introduced support for MiMiC QM/MM interface */
128     tpxv_PullAverage, /**< Added possibility to output average pull force and position */
129     tpxv_GenericInternalParameters, /**< Added internal parameters for mdrun modules*/
130     tpxv_VSite2FD,                  /**< Added 2FD type virtual site */
131     tpxv_AddSizeField, /**< Added field with information about the size of the serialized tpr file in bytes, excluding the header */
132     tpxv_Count         /**< the total number of tpxv versions */
133 };
134
135 /*! \brief Version number of the file format written to run input
136  * files by this version of the code.
137  *
138  * The tpx_version increases whenever the file format in the main
139  * development branch changes, due to an extension of the tpxv enum above.
140  * Backward compatibility for reading old run input files is maintained
141  * by checking this version number against that of the file and then using
142  * the correct code path.
143  *
144  * When developing a feature branch that needs to change the run input
145  * file format, change tpx_tag instead. */
146 static const int tpx_version = tpxv_Count - 1;
147
148
149 /* This number should only be increased when you edit the TOPOLOGY section
150  * or the HEADER of the tpx format.
151  * This way we can maintain forward compatibility too for all analysis tools
152  * and/or external programs that only need to know the atom/residue names,
153  * charges, and bond connectivity.
154  *
155  * It first appeared in tpx version 26, when I also moved the inputrecord
156  * to the end of the tpx file, so we can just skip it if we only
157  * want the topology.
158  *
159  * In particular, it must be increased when adding new elements to
160  * ftupd, so that old code can read new .tpr files.
161  *
162  * Updated for added field that contains the number of bytes of the tpr body, excluding the header.
163  */
164 static const int tpx_generation = 27;
165
166 /* This number should be the most recent backwards incompatible version
167  * I.e., if this number is 9, we cannot read tpx version 9 with this code.
168  */
169 static const int tpx_incompatible_version = 57; // GMX4.0 has version 58
170
171
172 /* Struct used to maintain tpx compatibility when function types are added */
173 typedef struct
174 {
175     int fvnr;  /* file version number in which the function type first appeared */
176     int ftype; /* function type */
177 } t_ftupd;
178
179 /*
180  * TODO The following three lines make little sense, please clarify if
181  * you've had to work out how ftupd works.
182  *
183  * The entries should be ordered in:
184  * 1. ascending function type number
185  * 2. ascending file version number
186  *
187  * Because we support reading of old .tpr file versions (even when
188  * mdrun can no longer run the simulation), we need to be able to read
189  * obsolete t_interaction_function types. Any data read from such
190  * fields is discarded. Their names have _NOLONGERUSED appended to
191  * them to make things clear.
192  */
193 static const t_ftupd ftupd[] = {
194     { 70, F_RESTRBONDS },
195     { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
196     { 76, F_LINEAR_ANGLES },
197     { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
198     { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
199     { 65, F_CMAP },
200     { 60, F_GB12_NOLONGERUSED },
201     { 61, F_GB13_NOLONGERUSED },
202     { 61, F_GB14_NOLONGERUSED },
203     { 72, F_GBPOL_NOLONGERUSED },
204     { 72, F_NPSOLVATION_NOLONGERUSED },
205     { 93, F_LJ_RECIP },
206     { 90, F_FBPOSRES },
207     { 69, F_VTEMP_NOLONGERUSED },
208     { 66, F_PDISPCORR },
209     { 76, F_ANHARM_POL },
210     { 79, F_DVDL_COUL },
211     {
212             79,
213             F_DVDL_VDW,
214     },
215     {
216             79,
217             F_DVDL_BONDED,
218     },
219     { 79, F_DVDL_RESTRAINT },
220     { 79, F_DVDL_TEMPERATURE },
221     { tpxv_GenericInternalParameters, F_DENSITYFITTING },
222     { tpxv_VSite2FD, F_VSITE2FD },
223 };
224 #define NFTUPD asize(ftupd)
225
226 /* Needed for backward compatibility */
227 #define MAXNODES 256
228
229 /**************************************************************
230  *
231  * Now the higer level routines that do io of the structures and arrays
232  *
233  **************************************************************/
234 static void do_pullgrp_tpx_pre95(gmx::ISerializer* serializer, t_pull_group* pgrp, t_pull_coord* pcrd)
235 {
236     rvec tmp;
237
238     serializer->doInt(&pgrp->nat);
239     if (serializer->reading())
240     {
241         snew(pgrp->ind, pgrp->nat);
242     }
243     serializer->doIntArray(pgrp->ind, pgrp->nat);
244     serializer->doInt(&pgrp->nweight);
245     if (serializer->reading())
246     {
247         snew(pgrp->weight, pgrp->nweight);
248     }
249     serializer->doRealArray(pgrp->weight, pgrp->nweight);
250     serializer->doInt(&pgrp->pbcatom);
251     serializer->doRvec(&pcrd->vec);
252     clear_rvec(pcrd->origin);
253     serializer->doRvec(&tmp);
254     pcrd->init = tmp[0];
255     serializer->doReal(&pcrd->rate);
256     serializer->doReal(&pcrd->k);
257     serializer->doReal(&pcrd->kB);
258 }
259
260 static void do_pull_group(gmx::ISerializer* serializer, t_pull_group* pgrp)
261 {
262     serializer->doInt(&pgrp->nat);
263     if (serializer->reading())
264     {
265         snew(pgrp->ind, pgrp->nat);
266     }
267     serializer->doIntArray(pgrp->ind, pgrp->nat);
268     serializer->doInt(&pgrp->nweight);
269     if (serializer->reading())
270     {
271         snew(pgrp->weight, pgrp->nweight);
272     }
273     serializer->doRealArray(pgrp->weight, pgrp->nweight);
274     serializer->doInt(&pgrp->pbcatom);
275 }
276
277 static void do_pull_coord(gmx::ISerializer* serializer,
278                           t_pull_coord*     pcrd,
279                           int               file_version,
280                           int               ePullOld,
281                           int               eGeomOld,
282                           ivec              dimOld)
283 {
284     if (file_version >= tpxv_PullCoordNGroup)
285     {
286         serializer->doInt(&pcrd->eType);
287         if (file_version >= tpxv_PullExternalPotential)
288         {
289             if (pcrd->eType == epullEXTERNAL)
290             {
291                 std::string buf;
292                 if (serializer->reading())
293                 {
294                     serializer->doString(&buf);
295                     pcrd->externalPotentialProvider = gmx_strdup(buf.c_str());
296                 }
297                 else
298                 {
299                     buf = pcrd->externalPotentialProvider;
300                     serializer->doString(&buf);
301                 }
302             }
303             else
304             {
305                 pcrd->externalPotentialProvider = nullptr;
306             }
307         }
308         else
309         {
310             if (serializer->reading())
311             {
312                 pcrd->externalPotentialProvider = nullptr;
313             }
314         }
315         /* Note that we try to support adding new geometries without
316          * changing the tpx version. This requires checks when printing the
317          * geometry string and a check and fatal_error in init_pull.
318          */
319         serializer->doInt(&pcrd->eGeom);
320         serializer->doInt(&pcrd->ngroup);
321         if (pcrd->ngroup <= c_pullCoordNgroupMax)
322         {
323             serializer->doIntArray(pcrd->group, pcrd->ngroup);
324         }
325         else
326         {
327             /* More groups in file than supported, this must be a new geometry
328              * that is not supported by our current code. Since we will not
329              * use the groups for this coord (checks in the pull and WHAM code
330              * ensure this), we can ignore the groups and set ngroup=0.
331              */
332             int* dum;
333             snew(dum, pcrd->ngroup);
334             serializer->doIntArray(dum, pcrd->ngroup);
335             sfree(dum);
336
337             pcrd->ngroup = 0;
338         }
339         serializer->doIvec(&pcrd->dim);
340     }
341     else
342     {
343         pcrd->ngroup = 2;
344         serializer->doInt(&pcrd->group[0]);
345         serializer->doInt(&pcrd->group[1]);
346         if (file_version >= tpxv_PullCoordTypeGeom)
347         {
348             pcrd->ngroup = (pcrd->eGeom == epullgDIRRELATIVE ? 4 : 2);
349             serializer->doInt(&pcrd->eType);
350             serializer->doInt(&pcrd->eGeom);
351             if (pcrd->ngroup == 4)
352             {
353                 serializer->doInt(&pcrd->group[2]);
354                 serializer->doInt(&pcrd->group[3]);
355             }
356             serializer->doIvec(&pcrd->dim);
357         }
358         else
359         {
360             pcrd->eType = ePullOld;
361             pcrd->eGeom = eGeomOld;
362             copy_ivec(dimOld, pcrd->dim);
363         }
364     }
365     serializer->doRvec(&pcrd->origin);
366     serializer->doRvec(&pcrd->vec);
367     if (file_version >= tpxv_PullCoordTypeGeom)
368     {
369         serializer->doBool(&pcrd->bStart);
370     }
371     else
372     {
373         /* This parameter is only printed, but not actually used by mdrun */
374         pcrd->bStart = FALSE;
375     }
376     serializer->doReal(&pcrd->init);
377     serializer->doReal(&pcrd->rate);
378     serializer->doReal(&pcrd->k);
379     serializer->doReal(&pcrd->kB);
380 }
381
382 static void do_expandedvals(gmx::ISerializer* serializer, t_expanded* expand, t_lambda* fepvals, int file_version)
383 {
384     int n_lambda = fepvals->n_lambda;
385
386     /* reset the lambda calculation window */
387     fepvals->lambda_start_n = 0;
388     fepvals->lambda_stop_n  = n_lambda;
389     if (file_version >= 79)
390     {
391         if (n_lambda > 0)
392         {
393             if (serializer->reading())
394             {
395                 snew(expand->init_lambda_weights, n_lambda);
396             }
397             serializer->doRealArray(expand->init_lambda_weights, n_lambda);
398             serializer->doBool(&expand->bInit_weights);
399         }
400
401         serializer->doInt(&expand->nstexpanded);
402         serializer->doInt(&expand->elmcmove);
403         serializer->doInt(&expand->elamstats);
404         serializer->doInt(&expand->lmc_repeats);
405         serializer->doInt(&expand->gibbsdeltalam);
406         serializer->doInt(&expand->lmc_forced_nstart);
407         serializer->doInt(&expand->lmc_seed);
408         serializer->doReal(&expand->mc_temp);
409         serializer->doBool(&expand->bSymmetrizedTMatrix);
410         serializer->doInt(&expand->nstTij);
411         serializer->doInt(&expand->minvarmin);
412         serializer->doInt(&expand->c_range);
413         serializer->doReal(&expand->wl_scale);
414         serializer->doReal(&expand->wl_ratio);
415         serializer->doReal(&expand->init_wl_delta);
416         serializer->doBool(&expand->bWLoneovert);
417         serializer->doInt(&expand->elmceq);
418         serializer->doInt(&expand->equil_steps);
419         serializer->doInt(&expand->equil_samples);
420         serializer->doInt(&expand->equil_n_at_lam);
421         serializer->doReal(&expand->equil_wl_delta);
422         serializer->doReal(&expand->equil_ratio);
423     }
424 }
425
426 static void do_simtempvals(gmx::ISerializer* serializer, t_simtemp* simtemp, int n_lambda, int file_version)
427 {
428     if (file_version >= 79)
429     {
430         serializer->doInt(&simtemp->eSimTempScale);
431         serializer->doReal(&simtemp->simtemp_high);
432         serializer->doReal(&simtemp->simtemp_low);
433         if (n_lambda > 0)
434         {
435             if (serializer->reading())
436             {
437                 snew(simtemp->temperatures, n_lambda);
438             }
439             serializer->doRealArray(simtemp->temperatures, n_lambda);
440         }
441     }
442 }
443
444 static void do_imd(gmx::ISerializer* serializer, t_IMD* imd)
445 {
446     serializer->doInt(&imd->nat);
447     if (serializer->reading())
448     {
449         snew(imd->ind, imd->nat);
450     }
451     serializer->doIntArray(imd->ind, imd->nat);
452 }
453
454 static void do_fepvals(gmx::ISerializer* serializer, t_lambda* fepvals, int file_version)
455 {
456     /* i is defined in the ndo_double macro; use g to iterate. */
457     int  g;
458     real rdum;
459
460     /* free energy values */
461
462     if (file_version >= 79)
463     {
464         serializer->doInt(&fepvals->init_fep_state);
465         serializer->doDouble(&fepvals->init_lambda);
466         serializer->doDouble(&fepvals->delta_lambda);
467     }
468     else if (file_version >= 59)
469     {
470         serializer->doDouble(&fepvals->init_lambda);
471         serializer->doDouble(&fepvals->delta_lambda);
472     }
473     else
474     {
475         serializer->doReal(&rdum);
476         fepvals->init_lambda = rdum;
477         serializer->doReal(&rdum);
478         fepvals->delta_lambda = rdum;
479     }
480     if (file_version >= 79)
481     {
482         serializer->doInt(&fepvals->n_lambda);
483         if (serializer->reading())
484         {
485             snew(fepvals->all_lambda, efptNR);
486         }
487         for (g = 0; g < efptNR; g++)
488         {
489             if (fepvals->n_lambda > 0)
490             {
491                 if (serializer->reading())
492                 {
493                     snew(fepvals->all_lambda[g], fepvals->n_lambda);
494                 }
495                 serializer->doDoubleArray(fepvals->all_lambda[g], fepvals->n_lambda);
496                 serializer->doBoolArray(fepvals->separate_dvdl, efptNR);
497             }
498             else if (fepvals->init_lambda >= 0)
499             {
500                 fepvals->separate_dvdl[efptFEP] = TRUE;
501             }
502         }
503     }
504     else if (file_version >= 64)
505     {
506         serializer->doInt(&fepvals->n_lambda);
507         if (serializer->reading())
508         {
509             int g;
510
511             snew(fepvals->all_lambda, efptNR);
512             /* still allocate the all_lambda array's contents. */
513             for (g = 0; g < efptNR; g++)
514             {
515                 if (fepvals->n_lambda > 0)
516                 {
517                     snew(fepvals->all_lambda[g], fepvals->n_lambda);
518                 }
519             }
520         }
521         serializer->doDoubleArray(fepvals->all_lambda[efptFEP], fepvals->n_lambda);
522         if (fepvals->init_lambda >= 0)
523         {
524             int g, h;
525
526             fepvals->separate_dvdl[efptFEP] = TRUE;
527
528             if (serializer->reading())
529             {
530                 /* copy the contents of the efptFEP lambda component to all
531                    the other components */
532                 for (g = 0; g < efptNR; g++)
533                 {
534                     for (h = 0; h < fepvals->n_lambda; h++)
535                     {
536                         if (g != efptFEP)
537                         {
538                             fepvals->all_lambda[g][h] = fepvals->all_lambda[efptFEP][h];
539                         }
540                     }
541                 }
542             }
543         }
544     }
545     else
546     {
547         fepvals->n_lambda   = 0;
548         fepvals->all_lambda = nullptr;
549         if (fepvals->init_lambda >= 0)
550         {
551             fepvals->separate_dvdl[efptFEP] = TRUE;
552         }
553     }
554     serializer->doReal(&fepvals->sc_alpha);
555     serializer->doInt(&fepvals->sc_power);
556     if (file_version >= 79)
557     {
558         serializer->doReal(&fepvals->sc_r_power);
559     }
560     else
561     {
562         fepvals->sc_r_power = 6.0;
563     }
564     serializer->doReal(&fepvals->sc_sigma);
565     if (serializer->reading())
566     {
567         if (file_version >= 71)
568         {
569             fepvals->sc_sigma_min = fepvals->sc_sigma;
570         }
571         else
572         {
573             fepvals->sc_sigma_min = 0;
574         }
575     }
576     if (file_version >= 79)
577     {
578         serializer->doBool(&fepvals->bScCoul);
579     }
580     else
581     {
582         fepvals->bScCoul = TRUE;
583     }
584     if (file_version >= 64)
585     {
586         serializer->doInt(&fepvals->nstdhdl);
587     }
588     else
589     {
590         fepvals->nstdhdl = 1;
591     }
592
593     if (file_version >= 73)
594     {
595         serializer->doInt(&fepvals->separate_dhdl_file);
596         serializer->doInt(&fepvals->dhdl_derivatives);
597     }
598     else
599     {
600         fepvals->separate_dhdl_file = esepdhdlfileYES;
601         fepvals->dhdl_derivatives   = edhdlderivativesYES;
602     }
603     if (file_version >= 71)
604     {
605         serializer->doInt(&fepvals->dh_hist_size);
606         serializer->doDouble(&fepvals->dh_hist_spacing);
607     }
608     else
609     {
610         fepvals->dh_hist_size    = 0;
611         fepvals->dh_hist_spacing = 0.1;
612     }
613     if (file_version >= 79)
614     {
615         serializer->doInt(&fepvals->edHdLPrintEnergy);
616     }
617     else
618     {
619         fepvals->edHdLPrintEnergy = edHdLPrintEnergyNO;
620     }
621
622     /* handle lambda_neighbors */
623     if ((file_version >= 83 && file_version < 90) || file_version >= 92)
624     {
625         serializer->doInt(&fepvals->lambda_neighbors);
626         if ((fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0)
627             && (fepvals->init_lambda < 0))
628         {
629             fepvals->lambda_start_n = (fepvals->init_fep_state - fepvals->lambda_neighbors);
630             fepvals->lambda_stop_n  = (fepvals->init_fep_state + fepvals->lambda_neighbors + 1);
631             if (fepvals->lambda_start_n < 0)
632             {
633                 fepvals->lambda_start_n = 0;
634             }
635             if (fepvals->lambda_stop_n >= fepvals->n_lambda)
636             {
637                 fepvals->lambda_stop_n = fepvals->n_lambda;
638             }
639         }
640         else
641         {
642             fepvals->lambda_start_n = 0;
643             fepvals->lambda_stop_n  = fepvals->n_lambda;
644         }
645     }
646     else
647     {
648         fepvals->lambda_start_n = 0;
649         fepvals->lambda_stop_n  = fepvals->n_lambda;
650     }
651 }
652
653 static void do_awhBias(gmx::ISerializer* serializer, gmx::AwhBiasParams* awhBiasParams, int gmx_unused file_version)
654 {
655     serializer->doInt(&awhBiasParams->eTarget);
656     serializer->doDouble(&awhBiasParams->targetBetaScaling);
657     serializer->doDouble(&awhBiasParams->targetCutoff);
658     serializer->doInt(&awhBiasParams->eGrowth);
659     serializer->doInt(&awhBiasParams->bUserData);
660     serializer->doDouble(&awhBiasParams->errorInitial);
661     serializer->doInt(&awhBiasParams->ndim);
662     serializer->doInt(&awhBiasParams->shareGroup);
663     serializer->doBool(&awhBiasParams->equilibrateHistogram);
664
665     if (serializer->reading())
666     {
667         snew(awhBiasParams->dimParams, awhBiasParams->ndim);
668     }
669
670     for (int d = 0; d < awhBiasParams->ndim; d++)
671     {
672         gmx::AwhDimParams* dimParams = &awhBiasParams->dimParams[d];
673
674         serializer->doInt(&dimParams->eCoordProvider);
675         serializer->doInt(&dimParams->coordIndex);
676         serializer->doDouble(&dimParams->origin);
677         serializer->doDouble(&dimParams->end);
678         serializer->doDouble(&dimParams->period);
679         serializer->doDouble(&dimParams->forceConstant);
680         serializer->doDouble(&dimParams->diffusion);
681         serializer->doDouble(&dimParams->coordValueInit);
682         serializer->doDouble(&dimParams->coverDiameter);
683     }
684 }
685
686 static void do_awh(gmx::ISerializer* serializer, gmx::AwhParams* awhParams, int gmx_unused file_version)
687 {
688     serializer->doInt(&awhParams->numBias);
689     serializer->doInt(&awhParams->nstOut);
690     serializer->doInt64(&awhParams->seed);
691     serializer->doInt(&awhParams->nstSampleCoord);
692     serializer->doInt(&awhParams->numSamplesUpdateFreeEnergy);
693     serializer->doInt(&awhParams->ePotential);
694     serializer->doBool(&awhParams->shareBiasMultisim);
695
696     if (awhParams->numBias > 0)
697     {
698         if (serializer->reading())
699         {
700             snew(awhParams->awhBiasParams, awhParams->numBias);
701         }
702
703         for (int k = 0; k < awhParams->numBias; k++)
704         {
705             do_awhBias(serializer, &awhParams->awhBiasParams[k], file_version);
706         }
707     }
708 }
709
710 static void do_pull(gmx::ISerializer* serializer, pull_params_t* pull, int file_version, int ePullOld)
711 {
712     int  eGeomOld = -1;
713     ivec dimOld;
714     int  g;
715
716     if (file_version >= 95)
717     {
718         serializer->doInt(&pull->ngroup);
719     }
720     serializer->doInt(&pull->ncoord);
721     if (file_version < 95)
722     {
723         pull->ngroup = pull->ncoord + 1;
724     }
725     if (file_version < tpxv_PullCoordTypeGeom)
726     {
727         real dum;
728
729         serializer->doInt(&eGeomOld);
730         serializer->doIvec(&dimOld);
731         /* The inner cylinder radius, now removed */
732         serializer->doReal(&dum);
733     }
734     serializer->doReal(&pull->cylinder_r);
735     serializer->doReal(&pull->constr_tol);
736     if (file_version >= 95)
737     {
738         serializer->doBool(&pull->bPrintCOM);
739         /* With file_version < 95 this value is set below */
740     }
741     if (file_version >= tpxv_ReplacePullPrintCOM12)
742     {
743         serializer->doBool(&pull->bPrintRefValue);
744         serializer->doBool(&pull->bPrintComp);
745     }
746     else if (file_version >= tpxv_PullCoordTypeGeom)
747     {
748         int idum;
749         serializer->doInt(&idum); /* used to be bPrintCOM2 */
750         serializer->doBool(&pull->bPrintRefValue);
751         serializer->doBool(&pull->bPrintComp);
752     }
753     else
754     {
755         pull->bPrintRefValue = FALSE;
756         pull->bPrintComp     = TRUE;
757     }
758     serializer->doInt(&pull->nstxout);
759     serializer->doInt(&pull->nstfout);
760     if (file_version >= tpxv_PullPrevStepCOMAsReference)
761     {
762         serializer->doBool(&pull->bSetPbcRefToPrevStepCOM);
763     }
764     else
765     {
766         pull->bSetPbcRefToPrevStepCOM = FALSE;
767     }
768     if (serializer->reading())
769     {
770         snew(pull->group, pull->ngroup);
771         snew(pull->coord, pull->ncoord);
772     }
773     if (file_version < 95)
774     {
775         /* epullgPOS for position pulling, before epullgDIRPBC was removed */
776         if (eGeomOld == epullgDIRPBC)
777         {
778             gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
779         }
780         if (eGeomOld > epullgDIRPBC)
781         {
782             eGeomOld -= 1;
783         }
784
785         for (g = 0; g < pull->ngroup; g++)
786         {
787             /* We read and ignore a pull coordinate for group 0 */
788             do_pullgrp_tpx_pre95(serializer, &pull->group[g], &pull->coord[std::max(g - 1, 0)]);
789             if (g > 0)
790             {
791                 pull->coord[g - 1].group[0] = 0;
792                 pull->coord[g - 1].group[1] = g;
793             }
794         }
795
796         pull->bPrintCOM = (pull->group[0].nat > 0);
797     }
798     else
799     {
800         for (g = 0; g < pull->ngroup; g++)
801         {
802             do_pull_group(serializer, &pull->group[g]);
803         }
804         for (g = 0; g < pull->ncoord; g++)
805         {
806             do_pull_coord(serializer, &pull->coord[g], file_version, ePullOld, eGeomOld, dimOld);
807         }
808     }
809     if (file_version >= tpxv_PullAverage)
810     {
811         gmx_bool v;
812
813         v = pull->bXOutAverage;
814         serializer->doBool(&v);
815         pull->bXOutAverage = v;
816         v                  = pull->bFOutAverage;
817         serializer->doBool(&v);
818         pull->bFOutAverage = v;
819     }
820 }
821
822
823 static void do_rotgrp(gmx::ISerializer* serializer, t_rotgrp* rotg)
824 {
825     serializer->doInt(&rotg->eType);
826     serializer->doInt(&rotg->bMassW);
827     serializer->doInt(&rotg->nat);
828     if (serializer->reading())
829     {
830         snew(rotg->ind, rotg->nat);
831     }
832     serializer->doIntArray(rotg->ind, rotg->nat);
833     if (serializer->reading())
834     {
835         snew(rotg->x_ref, rotg->nat);
836     }
837     serializer->doRvecArray(rotg->x_ref, rotg->nat);
838     serializer->doRvec(&rotg->inputVec);
839     serializer->doRvec(&rotg->pivot);
840     serializer->doReal(&rotg->rate);
841     serializer->doReal(&rotg->k);
842     serializer->doReal(&rotg->slab_dist);
843     serializer->doReal(&rotg->min_gaussian);
844     serializer->doReal(&rotg->eps);
845     serializer->doInt(&rotg->eFittype);
846     serializer->doInt(&rotg->PotAngle_nstep);
847     serializer->doReal(&rotg->PotAngle_step);
848 }
849
850 static void do_rot(gmx::ISerializer* serializer, t_rot* rot)
851 {
852     int g;
853
854     serializer->doInt(&rot->ngrp);
855     serializer->doInt(&rot->nstrout);
856     serializer->doInt(&rot->nstsout);
857     if (serializer->reading())
858     {
859         snew(rot->grp, rot->ngrp);
860     }
861     for (g = 0; g < rot->ngrp; g++)
862     {
863         do_rotgrp(serializer, &rot->grp[g]);
864     }
865 }
866
867
868 static void do_swapgroup(gmx::ISerializer* serializer, t_swapGroup* g)
869 {
870
871     /* Name of the group or molecule */
872     std::string buf;
873     if (serializer->reading())
874     {
875         serializer->doString(&buf);
876         g->molname = gmx_strdup(buf.c_str());
877     }
878     else
879     {
880         buf = g->molname;
881         serializer->doString(&buf);
882     }
883
884     /* Number of atoms in the group */
885     serializer->doInt(&g->nat);
886
887     /* The group's atom indices */
888     if (serializer->reading())
889     {
890         snew(g->ind, g->nat);
891     }
892     serializer->doIntArray(g->ind, g->nat);
893
894     /* Requested counts for compartments A and B */
895     serializer->doIntArray(g->nmolReq, eCompNR);
896 }
897
898 static void do_swapcoords_tpx(gmx::ISerializer* serializer, t_swapcoords* swap, int file_version)
899 {
900     /* Enums for better readability of the code */
901     enum
902     {
903         eCompA = 0,
904         eCompB
905     };
906     enum
907     {
908         eChannel0 = 0,
909         eChannel1
910     };
911
912
913     if (file_version >= tpxv_CompElPolyatomicIonsAndMultipleIonTypes)
914     {
915         /* The total number of swap groups is the sum of the fixed groups
916          * (split0, split1, solvent), and the user-defined groups (2+ types of ions)
917          */
918         serializer->doInt(&swap->ngrp);
919         if (serializer->reading())
920         {
921             snew(swap->grp, swap->ngrp);
922         }
923         for (int ig = 0; ig < swap->ngrp; ig++)
924         {
925             do_swapgroup(serializer, &swap->grp[ig]);
926         }
927         serializer->doBool(&swap->massw_split[eChannel0]);
928         serializer->doBool(&swap->massw_split[eChannel1]);
929         serializer->doInt(&swap->nstswap);
930         serializer->doInt(&swap->nAverage);
931         serializer->doReal(&swap->threshold);
932         serializer->doReal(&swap->cyl0r);
933         serializer->doReal(&swap->cyl0u);
934         serializer->doReal(&swap->cyl0l);
935         serializer->doReal(&swap->cyl1r);
936         serializer->doReal(&swap->cyl1u);
937         serializer->doReal(&swap->cyl1l);
938     }
939     else
940     {
941         /*** Support reading older CompEl .tpr files ***/
942
943         /* In the original CompEl .tpr files, we always have 5 groups: */
944         swap->ngrp = 5;
945         snew(swap->grp, swap->ngrp);
946
947         swap->grp[eGrpSplit0].molname  = gmx_strdup("split0");  // group 0: split0
948         swap->grp[eGrpSplit1].molname  = gmx_strdup("split1");  // group 1: split1
949         swap->grp[eGrpSolvent].molname = gmx_strdup("solvent"); // group 2: solvent
950         swap->grp[3].molname           = gmx_strdup("anions");  // group 3: anions
951         swap->grp[4].molname           = gmx_strdup("cations"); // group 4: cations
952
953         serializer->doInt(&swap->grp[3].nat);
954         serializer->doInt(&swap->grp[eGrpSolvent].nat);
955         serializer->doInt(&swap->grp[eGrpSplit0].nat);
956         serializer->doBool(&swap->massw_split[eChannel0]);
957         serializer->doInt(&swap->grp[eGrpSplit1].nat);
958         serializer->doBool(&swap->massw_split[eChannel1]);
959         serializer->doInt(&swap->nstswap);
960         serializer->doInt(&swap->nAverage);
961         serializer->doReal(&swap->threshold);
962         serializer->doReal(&swap->cyl0r);
963         serializer->doReal(&swap->cyl0u);
964         serializer->doReal(&swap->cyl0l);
965         serializer->doReal(&swap->cyl1r);
966         serializer->doReal(&swap->cyl1u);
967         serializer->doReal(&swap->cyl1l);
968
969         // The order[] array keeps compatibility with older .tpr files
970         // by reading in the groups in the classic order
971         {
972             const int order[4] = { 3, eGrpSolvent, eGrpSplit0, eGrpSplit1 };
973
974             for (int ig = 0; ig < 4; ig++)
975             {
976                 int g = order[ig];
977                 snew(swap->grp[g].ind, swap->grp[g].nat);
978                 serializer->doIntArray(swap->grp[g].ind, swap->grp[g].nat);
979             }
980         }
981
982         for (int j = eCompA; j <= eCompB; j++)
983         {
984             serializer->doInt(&swap->grp[3].nmolReq[j]); // group 3 = anions
985             serializer->doInt(&swap->grp[4].nmolReq[j]); // group 4 = cations
986         }
987     } /* End support reading older CompEl .tpr files */
988
989     if (file_version >= tpxv_CompElWithSwapLayerOffset)
990     {
991         serializer->doReal(&swap->bulkOffset[eCompA]);
992         serializer->doReal(&swap->bulkOffset[eCompB]);
993     }
994 }
995
996 static void do_legacy_efield(gmx::ISerializer* serializer, gmx::KeyValueTreeObjectBuilder* root)
997 {
998     const char* const dimName[] = { "x", "y", "z" };
999
1000     auto appliedForcesObj = root->addObject("applied-forces");
1001     auto efieldObj        = appliedForcesObj.addObject("electric-field");
1002     // The content of the tpr file for this feature has
1003     // been the same since gromacs 4.0 that was used for
1004     // developing.
1005     for (int j = 0; j < DIM; ++j)
1006     {
1007         int n, nt;
1008         serializer->doInt(&n);
1009         serializer->doInt(&nt);
1010         std::vector<real> aa(n + 1), phi(nt + 1), at(nt + 1), phit(nt + 1);
1011         serializer->doRealArray(aa.data(), n);
1012         serializer->doRealArray(phi.data(), n);
1013         serializer->doRealArray(at.data(), nt);
1014         serializer->doRealArray(phit.data(), nt);
1015         if (n > 0)
1016         {
1017             if (n > 1 || nt > 1)
1018             {
1019                 gmx_fatal(FARGS,
1020                           "Can not handle tpr files with more than one electric field term per "
1021                           "direction.");
1022             }
1023             auto dimObj = efieldObj.addObject(dimName[j]);
1024             dimObj.addValue<real>("E0", aa[0]);
1025             dimObj.addValue<real>("omega", at[0]);
1026             dimObj.addValue<real>("t0", phi[0]);
1027             dimObj.addValue<real>("sigma", phit[0]);
1028         }
1029     }
1030 }
1031
1032
1033 static void do_inputrec(gmx::ISerializer* serializer, t_inputrec* ir, int file_version)
1034 {
1035     int      i, j, k, idum = 0;
1036     real     rdum;
1037     gmx_bool bdum = false;
1038
1039     if (file_version != tpx_version)
1040     {
1041         /* Give a warning about features that are not accessible */
1042         fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n", file_version, tpx_version);
1043     }
1044
1045     if (file_version == 0)
1046     {
1047         return;
1048     }
1049
1050     gmx::KeyValueTreeBuilder       paramsBuilder;
1051     gmx::KeyValueTreeObjectBuilder paramsObj = paramsBuilder.rootObject();
1052
1053     /* Basic inputrec stuff */
1054     serializer->doInt(&ir->eI);
1055     if (file_version >= 62)
1056     {
1057         serializer->doInt64(&ir->nsteps);
1058     }
1059     else
1060     {
1061         serializer->doInt(&idum);
1062         ir->nsteps = idum;
1063     }
1064
1065     if (file_version >= 62)
1066     {
1067         serializer->doInt64(&ir->init_step);
1068     }
1069     else
1070     {
1071         serializer->doInt(&idum);
1072         ir->init_step = idum;
1073     }
1074
1075     serializer->doInt(&ir->simulation_part);
1076
1077     if (file_version >= 67)
1078     {
1079         serializer->doInt(&ir->nstcalcenergy);
1080     }
1081     else
1082     {
1083         ir->nstcalcenergy = 1;
1084     }
1085     if (file_version >= 81)
1086     {
1087         serializer->doInt(&ir->cutoff_scheme);
1088         if (file_version < 94)
1089         {
1090             ir->cutoff_scheme = 1 - ir->cutoff_scheme;
1091         }
1092     }
1093     else
1094     {
1095         ir->cutoff_scheme = ecutsGROUP;
1096     }
1097     serializer->doInt(&idum); /* used to be ns_type; not used anymore */
1098     serializer->doInt(&ir->nstlist);
1099     serializer->doInt(&idum); /* used to be ndelta; not used anymore */
1100
1101     serializer->doReal(&ir->rtpi);
1102
1103     serializer->doInt(&ir->nstcomm);
1104     serializer->doInt(&ir->comm_mode);
1105
1106     /* ignore nstcheckpoint */
1107     if (file_version < tpxv_RemoveObsoleteParameters1)
1108     {
1109         serializer->doInt(&idum);
1110     }
1111
1112     serializer->doInt(&ir->nstcgsteep);
1113
1114     serializer->doInt(&ir->nbfgscorr);
1115
1116     serializer->doInt(&ir->nstlog);
1117     serializer->doInt(&ir->nstxout);
1118     serializer->doInt(&ir->nstvout);
1119     serializer->doInt(&ir->nstfout);
1120     serializer->doInt(&ir->nstenergy);
1121     serializer->doInt(&ir->nstxout_compressed);
1122     if (file_version >= 59)
1123     {
1124         serializer->doDouble(&ir->init_t);
1125         serializer->doDouble(&ir->delta_t);
1126     }
1127     else
1128     {
1129         serializer->doReal(&rdum);
1130         ir->init_t = rdum;
1131         serializer->doReal(&rdum);
1132         ir->delta_t = rdum;
1133     }
1134     serializer->doReal(&ir->x_compression_precision);
1135     if (file_version >= 81)
1136     {
1137         serializer->doReal(&ir->verletbuf_tol);
1138     }
1139     else
1140     {
1141         ir->verletbuf_tol = 0;
1142     }
1143     serializer->doReal(&ir->rlist);
1144     if (file_version >= 67 && file_version < tpxv_RemoveTwinRange)
1145     {
1146         if (serializer->reading())
1147         {
1148             // Reading such a file version could invoke the twin-range
1149             // scheme, about which mdrun should give a fatal error.
1150             real dummy_rlistlong = -1;
1151             serializer->doReal(&dummy_rlistlong);
1152
1153             if (ir->rlist > 0 && (dummy_rlistlong == 0 || dummy_rlistlong > ir->rlist))
1154             {
1155                 // Get mdrun to issue an error (regardless of
1156                 // ir->cutoff_scheme).
1157                 ir->useTwinRange = true;
1158             }
1159             else
1160             {
1161                 // grompp used to set rlistlong actively. Users were
1162                 // probably also confused and set rlistlong == rlist.
1163                 // However, in all remaining cases, it is safe to let
1164                 // mdrun proceed normally.
1165                 ir->useTwinRange = false;
1166             }
1167         }
1168     }
1169     else
1170     {
1171         // No need to read or write anything
1172         ir->useTwinRange = false;
1173     }
1174     if (file_version >= 82 && file_version != 90)
1175     {
1176         // Multiple time-stepping is no longer enabled, but the old
1177         // support required the twin-range scheme, for which mdrun
1178         // already emits a fatal error.
1179         int dummy_nstcalclr = -1;
1180         serializer->doInt(&dummy_nstcalclr);
1181     }
1182     serializer->doInt(&ir->coulombtype);
1183     if (file_version >= 81)
1184     {
1185         serializer->doInt(&ir->coulomb_modifier);
1186     }
1187     else
1188     {
1189         ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1190     }
1191     serializer->doReal(&ir->rcoulomb_switch);
1192     serializer->doReal(&ir->rcoulomb);
1193     serializer->doInt(&ir->vdwtype);
1194     if (file_version >= 81)
1195     {
1196         serializer->doInt(&ir->vdw_modifier);
1197     }
1198     else
1199     {
1200         ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1201     }
1202     serializer->doReal(&ir->rvdw_switch);
1203     serializer->doReal(&ir->rvdw);
1204     serializer->doInt(&ir->eDispCorr);
1205     serializer->doReal(&ir->epsilon_r);
1206     serializer->doReal(&ir->epsilon_rf);
1207     serializer->doReal(&ir->tabext);
1208
1209     // This permits reading a .tpr file that used implicit solvent,
1210     // and later permitting mdrun to refuse to run it.
1211     if (serializer->reading())
1212     {
1213         if (file_version < tpxv_RemoveImplicitSolvation)
1214         {
1215             serializer->doInt(&idum);
1216             serializer->doInt(&idum);
1217             serializer->doReal(&rdum);
1218             serializer->doReal(&rdum);
1219             serializer->doInt(&idum);
1220             ir->implicit_solvent = (idum > 0);
1221         }
1222         else
1223         {
1224             ir->implicit_solvent = false;
1225         }
1226         if (file_version < tpxv_RemoveImplicitSolvation)
1227         {
1228             serializer->doReal(&rdum);
1229             serializer->doReal(&rdum);
1230             serializer->doReal(&rdum);
1231             serializer->doReal(&rdum);
1232             if (file_version >= 60)
1233             {
1234                 serializer->doReal(&rdum);
1235                 serializer->doInt(&idum);
1236             }
1237             serializer->doReal(&rdum);
1238         }
1239     }
1240
1241     if (file_version >= 81)
1242     {
1243         serializer->doReal(&ir->fourier_spacing);
1244     }
1245     else
1246     {
1247         ir->fourier_spacing = 0.0;
1248     }
1249     serializer->doInt(&ir->nkx);
1250     serializer->doInt(&ir->nky);
1251     serializer->doInt(&ir->nkz);
1252     serializer->doInt(&ir->pme_order);
1253     serializer->doReal(&ir->ewald_rtol);
1254
1255     if (file_version >= 93)
1256     {
1257         serializer->doReal(&ir->ewald_rtol_lj);
1258     }
1259     else
1260     {
1261         ir->ewald_rtol_lj = ir->ewald_rtol;
1262     }
1263     serializer->doInt(&ir->ewald_geometry);
1264     serializer->doReal(&ir->epsilon_surface);
1265
1266     /* ignore bOptFFT */
1267     if (file_version < tpxv_RemoveObsoleteParameters1)
1268     {
1269         serializer->doBool(&bdum);
1270     }
1271
1272     if (file_version >= 93)
1273     {
1274         serializer->doInt(&ir->ljpme_combination_rule);
1275     }
1276     serializer->doBool(&ir->bContinuation);
1277     serializer->doInt(&ir->etc);
1278     /* before version 18, ir->etc was a gmx_bool (ir->btc),
1279      * but the values 0 and 1 still mean no and
1280      * berendsen temperature coupling, respectively.
1281      */
1282     if (file_version >= 79)
1283     {
1284         serializer->doBool(&ir->bPrintNHChains);
1285     }
1286     if (file_version >= 71)
1287     {
1288         serializer->doInt(&ir->nsttcouple);
1289     }
1290     else
1291     {
1292         ir->nsttcouple = ir->nstcalcenergy;
1293     }
1294     serializer->doInt(&ir->epc);
1295     serializer->doInt(&ir->epct);
1296     if (file_version >= 71)
1297     {
1298         serializer->doInt(&ir->nstpcouple);
1299     }
1300     else
1301     {
1302         ir->nstpcouple = ir->nstcalcenergy;
1303     }
1304     serializer->doReal(&ir->tau_p);
1305     serializer->doRvec(&ir->ref_p[XX]);
1306     serializer->doRvec(&ir->ref_p[YY]);
1307     serializer->doRvec(&ir->ref_p[ZZ]);
1308     serializer->doRvec(&ir->compress[XX]);
1309     serializer->doRvec(&ir->compress[YY]);
1310     serializer->doRvec(&ir->compress[ZZ]);
1311     serializer->doInt(&ir->refcoord_scaling);
1312     serializer->doRvec(&ir->posres_com);
1313     serializer->doRvec(&ir->posres_comB);
1314
1315     if (file_version < 79)
1316     {
1317         serializer->doInt(&ir->andersen_seed);
1318     }
1319     else
1320     {
1321         ir->andersen_seed = 0;
1322     }
1323
1324     serializer->doReal(&ir->shake_tol);
1325
1326     serializer->doInt(&ir->efep);
1327     do_fepvals(serializer, ir->fepvals, file_version);
1328
1329     if (file_version >= 79)
1330     {
1331         serializer->doBool(&ir->bSimTemp);
1332         if (ir->bSimTemp)
1333         {
1334             ir->bSimTemp = TRUE;
1335         }
1336     }
1337     else
1338     {
1339         ir->bSimTemp = FALSE;
1340     }
1341     if (ir->bSimTemp)
1342     {
1343         do_simtempvals(serializer, ir->simtempvals, ir->fepvals->n_lambda, file_version);
1344     }
1345
1346     if (file_version >= 79)
1347     {
1348         serializer->doBool(&ir->bExpanded);
1349         if (ir->bExpanded)
1350         {
1351             ir->bExpanded = TRUE;
1352         }
1353         else
1354         {
1355             ir->bExpanded = FALSE;
1356         }
1357     }
1358     if (ir->bExpanded)
1359     {
1360         do_expandedvals(serializer, ir->expandedvals, ir->fepvals, file_version);
1361     }
1362
1363     serializer->doInt(&ir->eDisre);
1364     serializer->doInt(&ir->eDisreWeighting);
1365     serializer->doBool(&ir->bDisreMixed);
1366     serializer->doReal(&ir->dr_fc);
1367     serializer->doReal(&ir->dr_tau);
1368     serializer->doInt(&ir->nstdisreout);
1369     serializer->doReal(&ir->orires_fc);
1370     serializer->doReal(&ir->orires_tau);
1371     serializer->doInt(&ir->nstorireout);
1372
1373     /* ignore dihre_fc */
1374     if (file_version < 79)
1375     {
1376         serializer->doReal(&rdum);
1377     }
1378
1379     serializer->doReal(&ir->em_stepsize);
1380     serializer->doReal(&ir->em_tol);
1381     serializer->doBool(&ir->bShakeSOR);
1382     serializer->doInt(&ir->niter);
1383     serializer->doReal(&ir->fc_stepsize);
1384     serializer->doInt(&ir->eConstrAlg);
1385     serializer->doInt(&ir->nProjOrder);
1386     serializer->doReal(&ir->LincsWarnAngle);
1387     serializer->doInt(&ir->nLincsIter);
1388     serializer->doReal(&ir->bd_fric);
1389     if (file_version >= tpxv_Use64BitRandomSeed)
1390     {
1391         serializer->doInt64(&ir->ld_seed);
1392     }
1393     else
1394     {
1395         serializer->doInt(&idum);
1396         ir->ld_seed = idum;
1397     }
1398
1399     for (i = 0; i < DIM; i++)
1400     {
1401         serializer->doRvec(&ir->deform[i]);
1402     }
1403     serializer->doReal(&ir->cos_accel);
1404
1405     serializer->doInt(&ir->userint1);
1406     serializer->doInt(&ir->userint2);
1407     serializer->doInt(&ir->userint3);
1408     serializer->doInt(&ir->userint4);
1409     serializer->doReal(&ir->userreal1);
1410     serializer->doReal(&ir->userreal2);
1411     serializer->doReal(&ir->userreal3);
1412     serializer->doReal(&ir->userreal4);
1413
1414     /* AdResS is removed, but we need to be able to read old files,
1415        and let mdrun refuse to run them */
1416     if (file_version >= 77 && file_version < tpxv_RemoveAdress)
1417     {
1418         serializer->doBool(&ir->bAdress);
1419         if (ir->bAdress)
1420         {
1421             int  idum, numThermoForceGroups, numEnergyGroups;
1422             real rdum;
1423             rvec rvecdum;
1424             serializer->doInt(&idum);
1425             serializer->doReal(&rdum);
1426             serializer->doReal(&rdum);
1427             serializer->doReal(&rdum);
1428             serializer->doInt(&idum);
1429             serializer->doInt(&idum);
1430             serializer->doRvec(&rvecdum);
1431             serializer->doInt(&numThermoForceGroups);
1432             serializer->doReal(&rdum);
1433             serializer->doInt(&numEnergyGroups);
1434             serializer->doInt(&idum);
1435
1436             if (numThermoForceGroups > 0)
1437             {
1438                 std::vector<int> idumn(numThermoForceGroups);
1439                 serializer->doIntArray(idumn.data(), idumn.size());
1440             }
1441             if (numEnergyGroups > 0)
1442             {
1443                 std::vector<int> idumn(numEnergyGroups);
1444                 serializer->doIntArray(idumn.data(), idumn.size());
1445             }
1446         }
1447     }
1448     else
1449     {
1450         ir->bAdress = FALSE;
1451     }
1452
1453     /* pull stuff */
1454     {
1455         int ePullOld = 0;
1456
1457         if (file_version >= tpxv_PullCoordTypeGeom)
1458         {
1459             serializer->doBool(&ir->bPull);
1460         }
1461         else
1462         {
1463             serializer->doInt(&ePullOld);
1464             ir->bPull = (ePullOld > 0);
1465             /* We removed the first ePull=ePullNo for the enum */
1466             ePullOld -= 1;
1467         }
1468         if (ir->bPull)
1469         {
1470             if (serializer->reading())
1471             {
1472                 snew(ir->pull, 1);
1473             }
1474             do_pull(serializer, ir->pull, file_version, ePullOld);
1475         }
1476     }
1477
1478     if (file_version >= tpxv_AcceleratedWeightHistogram)
1479     {
1480         serializer->doBool(&ir->bDoAwh);
1481
1482         if (ir->bDoAwh)
1483         {
1484             if (serializer->reading())
1485             {
1486                 snew(ir->awhParams, 1);
1487             }
1488             do_awh(serializer, ir->awhParams, file_version);
1489         }
1490     }
1491     else
1492     {
1493         ir->bDoAwh = FALSE;
1494     }
1495
1496     /* Enforced rotation */
1497     if (file_version >= 74)
1498     {
1499         serializer->doBool(&ir->bRot);
1500         if (ir->bRot)
1501         {
1502             if (serializer->reading())
1503             {
1504                 snew(ir->rot, 1);
1505             }
1506             do_rot(serializer, ir->rot);
1507         }
1508     }
1509     else
1510     {
1511         ir->bRot = FALSE;
1512     }
1513
1514     /* Interactive molecular dynamics */
1515     if (file_version >= tpxv_InteractiveMolecularDynamics)
1516     {
1517         serializer->doBool(&ir->bIMD);
1518         if (ir->bIMD)
1519         {
1520             if (serializer->reading())
1521             {
1522                 snew(ir->imd, 1);
1523             }
1524             do_imd(serializer, ir->imd);
1525         }
1526     }
1527     else
1528     {
1529         /* We don't support IMD sessions for old .tpr files */
1530         ir->bIMD = FALSE;
1531     }
1532
1533     /* grpopts stuff */
1534     serializer->doInt(&ir->opts.ngtc);
1535     if (file_version >= 69)
1536     {
1537         serializer->doInt(&ir->opts.nhchainlength);
1538     }
1539     else
1540     {
1541         ir->opts.nhchainlength = 1;
1542     }
1543     serializer->doInt(&ir->opts.ngacc);
1544     serializer->doInt(&ir->opts.ngfrz);
1545     serializer->doInt(&ir->opts.ngener);
1546
1547     if (serializer->reading())
1548     {
1549         snew(ir->opts.nrdf, ir->opts.ngtc);
1550         snew(ir->opts.ref_t, ir->opts.ngtc);
1551         snew(ir->opts.annealing, ir->opts.ngtc);
1552         snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1553         snew(ir->opts.anneal_time, ir->opts.ngtc);
1554         snew(ir->opts.anneal_temp, ir->opts.ngtc);
1555         snew(ir->opts.tau_t, ir->opts.ngtc);
1556         snew(ir->opts.nFreeze, ir->opts.ngfrz);
1557         snew(ir->opts.acc, ir->opts.ngacc);
1558         snew(ir->opts.egp_flags, ir->opts.ngener * ir->opts.ngener);
1559     }
1560     if (ir->opts.ngtc > 0)
1561     {
1562         serializer->doRealArray(ir->opts.nrdf, ir->opts.ngtc);
1563         serializer->doRealArray(ir->opts.ref_t, ir->opts.ngtc);
1564         serializer->doRealArray(ir->opts.tau_t, ir->opts.ngtc);
1565     }
1566     if (ir->opts.ngfrz > 0)
1567     {
1568         serializer->doIvecArray(ir->opts.nFreeze, ir->opts.ngfrz);
1569     }
1570     if (ir->opts.ngacc > 0)
1571     {
1572         serializer->doRvecArray(ir->opts.acc, ir->opts.ngacc);
1573     }
1574     serializer->doIntArray(ir->opts.egp_flags, ir->opts.ngener * ir->opts.ngener);
1575
1576     /* First read the lists with annealing and npoints for each group */
1577     serializer->doIntArray(ir->opts.annealing, ir->opts.ngtc);
1578     serializer->doIntArray(ir->opts.anneal_npoints, ir->opts.ngtc);
1579     for (j = 0; j < (ir->opts.ngtc); j++)
1580     {
1581         k = ir->opts.anneal_npoints[j];
1582         if (serializer->reading())
1583         {
1584             snew(ir->opts.anneal_time[j], k);
1585             snew(ir->opts.anneal_temp[j], k);
1586         }
1587         serializer->doRealArray(ir->opts.anneal_time[j], k);
1588         serializer->doRealArray(ir->opts.anneal_temp[j], k);
1589     }
1590     /* Walls */
1591     {
1592         serializer->doInt(&ir->nwall);
1593         serializer->doInt(&ir->wall_type);
1594         serializer->doReal(&ir->wall_r_linpot);
1595         serializer->doInt(&ir->wall_atomtype[0]);
1596         serializer->doInt(&ir->wall_atomtype[1]);
1597         serializer->doReal(&ir->wall_density[0]);
1598         serializer->doReal(&ir->wall_density[1]);
1599         serializer->doReal(&ir->wall_ewald_zfac);
1600     }
1601
1602     /* Cosine stuff for electric fields */
1603     if (file_version < tpxv_GenericParamsForElectricField)
1604     {
1605         do_legacy_efield(serializer, &paramsObj);
1606     }
1607
1608     /* Swap ions */
1609     if (file_version >= tpxv_ComputationalElectrophysiology)
1610     {
1611         serializer->doInt(&ir->eSwapCoords);
1612         if (ir->eSwapCoords != eswapNO)
1613         {
1614             if (serializer->reading())
1615             {
1616                 snew(ir->swap, 1);
1617             }
1618             do_swapcoords_tpx(serializer, ir->swap, file_version);
1619         }
1620     }
1621
1622     /* QMMM stuff */
1623     {
1624         serializer->doBool(&ir->bQMMM);
1625         serializer->doInt(&ir->QMMMscheme);
1626         serializer->doReal(&ir->scalefactor);
1627         serializer->doInt(&ir->opts.ngQM);
1628         if (serializer->reading())
1629         {
1630             snew(ir->opts.QMmethod, ir->opts.ngQM);
1631             snew(ir->opts.QMbasis, ir->opts.ngQM);
1632             snew(ir->opts.QMcharge, ir->opts.ngQM);
1633             snew(ir->opts.QMmult, ir->opts.ngQM);
1634             snew(ir->opts.bSH, ir->opts.ngQM);
1635             snew(ir->opts.CASorbitals, ir->opts.ngQM);
1636             snew(ir->opts.CASelectrons, ir->opts.ngQM);
1637             snew(ir->opts.SAon, ir->opts.ngQM);
1638             snew(ir->opts.SAoff, ir->opts.ngQM);
1639             snew(ir->opts.SAsteps, ir->opts.ngQM);
1640         }
1641         if (ir->opts.ngQM > 0 && ir->bQMMM)
1642         {
1643             serializer->doIntArray(ir->opts.QMmethod, ir->opts.ngQM);
1644             serializer->doIntArray(ir->opts.QMbasis, ir->opts.ngQM);
1645             serializer->doIntArray(ir->opts.QMcharge, ir->opts.ngQM);
1646             serializer->doIntArray(ir->opts.QMmult, ir->opts.ngQM);
1647             serializer->doBoolArray(ir->opts.bSH, ir->opts.ngQM);
1648             serializer->doIntArray(ir->opts.CASorbitals, ir->opts.ngQM);
1649             serializer->doIntArray(ir->opts.CASelectrons, ir->opts.ngQM);
1650             serializer->doRealArray(ir->opts.SAon, ir->opts.ngQM);
1651             serializer->doRealArray(ir->opts.SAoff, ir->opts.ngQM);
1652             serializer->doIntArray(ir->opts.SAsteps, ir->opts.ngQM);
1653             /* We leave in dummy i/o for removed parameters to avoid
1654              * changing the tpr format for every QMMM change.
1655              */
1656             std::vector<int> dummy(ir->opts.ngQM, 0);
1657             serializer->doIntArray(dummy.data(), ir->opts.ngQM);
1658             serializer->doIntArray(dummy.data(), ir->opts.ngQM);
1659         }
1660         /* end of QMMM stuff */
1661     }
1662
1663     if (file_version >= tpxv_GenericParamsForElectricField)
1664     {
1665         if (serializer->reading())
1666         {
1667             paramsObj.mergeObject(gmx::deserializeKeyValueTree(serializer));
1668         }
1669         else
1670         {
1671             GMX_RELEASE_ASSERT(ir->params != nullptr,
1672                                "Parameters should be present when writing inputrec");
1673             gmx::serializeKeyValueTree(*ir->params, serializer);
1674         }
1675     }
1676     if (serializer->reading())
1677     {
1678         ir->params = new gmx::KeyValueTreeObject(paramsBuilder.build());
1679         // Initialize internal parameters to an empty kvt for all tpr versions
1680         ir->internalParameters = std::make_unique<gmx::KeyValueTreeObject>();
1681     }
1682
1683     if (file_version >= tpxv_GenericInternalParameters)
1684     {
1685         if (serializer->reading())
1686         {
1687             ir->internalParameters =
1688                     std::make_unique<gmx::KeyValueTreeObject>(gmx::deserializeKeyValueTree(serializer));
1689         }
1690         else
1691         {
1692             GMX_RELEASE_ASSERT(ir->internalParameters != nullptr,
1693                                "Parameters should be present when writing inputrec");
1694             gmx::serializeKeyValueTree(*ir->internalParameters, serializer);
1695         }
1696     }
1697 }
1698
1699
1700 static void do_harm(gmx::ISerializer* serializer, t_iparams* iparams)
1701 {
1702     serializer->doReal(&iparams->harmonic.rA);
1703     serializer->doReal(&iparams->harmonic.krA);
1704     serializer->doReal(&iparams->harmonic.rB);
1705     serializer->doReal(&iparams->harmonic.krB);
1706 }
1707
1708 static void do_iparams(gmx::ISerializer* serializer, t_functype ftype, t_iparams* iparams, int file_version)
1709 {
1710     int  idum;
1711     real rdum;
1712
1713     switch (ftype)
1714     {
1715         case F_ANGLES:
1716         case F_G96ANGLES:
1717         case F_BONDS:
1718         case F_G96BONDS:
1719         case F_HARMONIC:
1720         case F_IDIHS:
1721             do_harm(serializer, iparams);
1722             if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && serializer->reading())
1723             {
1724                 /* Correct incorrect storage of parameters */
1725                 iparams->pdihs.phiB = iparams->pdihs.phiA;
1726                 iparams->pdihs.cpB  = iparams->pdihs.cpA;
1727             }
1728             break;
1729         case F_RESTRANGLES:
1730             serializer->doReal(&iparams->harmonic.rA);
1731             serializer->doReal(&iparams->harmonic.krA);
1732             break;
1733         case F_LINEAR_ANGLES:
1734             serializer->doReal(&iparams->linangle.klinA);
1735             serializer->doReal(&iparams->linangle.aA);
1736             serializer->doReal(&iparams->linangle.klinB);
1737             serializer->doReal(&iparams->linangle.aB);
1738             break;
1739         case F_FENEBONDS:
1740             serializer->doReal(&iparams->fene.bm);
1741             serializer->doReal(&iparams->fene.kb);
1742             break;
1743
1744         case F_RESTRBONDS:
1745             serializer->doReal(&iparams->restraint.lowA);
1746             serializer->doReal(&iparams->restraint.up1A);
1747             serializer->doReal(&iparams->restraint.up2A);
1748             serializer->doReal(&iparams->restraint.kA);
1749             serializer->doReal(&iparams->restraint.lowB);
1750             serializer->doReal(&iparams->restraint.up1B);
1751             serializer->doReal(&iparams->restraint.up2B);
1752             serializer->doReal(&iparams->restraint.kB);
1753             break;
1754         case F_TABBONDS:
1755         case F_TABBONDSNC:
1756         case F_TABANGLES:
1757         case F_TABDIHS:
1758             serializer->doReal(&iparams->tab.kA);
1759             serializer->doInt(&iparams->tab.table);
1760             serializer->doReal(&iparams->tab.kB);
1761             break;
1762         case F_CROSS_BOND_BONDS:
1763             serializer->doReal(&iparams->cross_bb.r1e);
1764             serializer->doReal(&iparams->cross_bb.r2e);
1765             serializer->doReal(&iparams->cross_bb.krr);
1766             break;
1767         case F_CROSS_BOND_ANGLES:
1768             serializer->doReal(&iparams->cross_ba.r1e);
1769             serializer->doReal(&iparams->cross_ba.r2e);
1770             serializer->doReal(&iparams->cross_ba.r3e);
1771             serializer->doReal(&iparams->cross_ba.krt);
1772             break;
1773         case F_UREY_BRADLEY:
1774             serializer->doReal(&iparams->u_b.thetaA);
1775             serializer->doReal(&iparams->u_b.kthetaA);
1776             serializer->doReal(&iparams->u_b.r13A);
1777             serializer->doReal(&iparams->u_b.kUBA);
1778             if (file_version >= 79)
1779             {
1780                 serializer->doReal(&iparams->u_b.thetaB);
1781                 serializer->doReal(&iparams->u_b.kthetaB);
1782                 serializer->doReal(&iparams->u_b.r13B);
1783                 serializer->doReal(&iparams->u_b.kUBB);
1784             }
1785             else
1786             {
1787                 iparams->u_b.thetaB  = iparams->u_b.thetaA;
1788                 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1789                 iparams->u_b.r13B    = iparams->u_b.r13A;
1790                 iparams->u_b.kUBB    = iparams->u_b.kUBA;
1791             }
1792             break;
1793         case F_QUARTIC_ANGLES:
1794             serializer->doReal(&iparams->qangle.theta);
1795             serializer->doRealArray(iparams->qangle.c, 5);
1796             break;
1797         case F_BHAM:
1798             serializer->doReal(&iparams->bham.a);
1799             serializer->doReal(&iparams->bham.b);
1800             serializer->doReal(&iparams->bham.c);
1801             break;
1802         case F_MORSE:
1803             serializer->doReal(&iparams->morse.b0A);
1804             serializer->doReal(&iparams->morse.cbA);
1805             serializer->doReal(&iparams->morse.betaA);
1806             if (file_version >= 79)
1807             {
1808                 serializer->doReal(&iparams->morse.b0B);
1809                 serializer->doReal(&iparams->morse.cbB);
1810                 serializer->doReal(&iparams->morse.betaB);
1811             }
1812             else
1813             {
1814                 iparams->morse.b0B   = iparams->morse.b0A;
1815                 iparams->morse.cbB   = iparams->morse.cbA;
1816                 iparams->morse.betaB = iparams->morse.betaA;
1817             }
1818             break;
1819         case F_CUBICBONDS:
1820             serializer->doReal(&iparams->cubic.b0);
1821             serializer->doReal(&iparams->cubic.kb);
1822             serializer->doReal(&iparams->cubic.kcub);
1823             break;
1824         case F_CONNBONDS: break;
1825         case F_POLARIZATION: serializer->doReal(&iparams->polarize.alpha); break;
1826         case F_ANHARM_POL:
1827             serializer->doReal(&iparams->anharm_polarize.alpha);
1828             serializer->doReal(&iparams->anharm_polarize.drcut);
1829             serializer->doReal(&iparams->anharm_polarize.khyp);
1830             break;
1831         case F_WATER_POL:
1832             serializer->doReal(&iparams->wpol.al_x);
1833             serializer->doReal(&iparams->wpol.al_y);
1834             serializer->doReal(&iparams->wpol.al_z);
1835             serializer->doReal(&iparams->wpol.rOH);
1836             serializer->doReal(&iparams->wpol.rHH);
1837             serializer->doReal(&iparams->wpol.rOD);
1838             break;
1839         case F_THOLE_POL:
1840             serializer->doReal(&iparams->thole.a);
1841             serializer->doReal(&iparams->thole.alpha1);
1842             serializer->doReal(&iparams->thole.alpha2);
1843             serializer->doReal(&iparams->thole.rfac);
1844             break;
1845         case F_LJ:
1846             serializer->doReal(&iparams->lj.c6);
1847             serializer->doReal(&iparams->lj.c12);
1848             break;
1849         case F_LJ14:
1850             serializer->doReal(&iparams->lj14.c6A);
1851             serializer->doReal(&iparams->lj14.c12A);
1852             serializer->doReal(&iparams->lj14.c6B);
1853             serializer->doReal(&iparams->lj14.c12B);
1854             break;
1855         case F_LJC14_Q:
1856             serializer->doReal(&iparams->ljc14.fqq);
1857             serializer->doReal(&iparams->ljc14.qi);
1858             serializer->doReal(&iparams->ljc14.qj);
1859             serializer->doReal(&iparams->ljc14.c6);
1860             serializer->doReal(&iparams->ljc14.c12);
1861             break;
1862         case F_LJC_PAIRS_NB:
1863             serializer->doReal(&iparams->ljcnb.qi);
1864             serializer->doReal(&iparams->ljcnb.qj);
1865             serializer->doReal(&iparams->ljcnb.c6);
1866             serializer->doReal(&iparams->ljcnb.c12);
1867             break;
1868         case F_PDIHS:
1869         case F_PIDIHS:
1870         case F_ANGRES:
1871         case F_ANGRESZ:
1872             serializer->doReal(&iparams->pdihs.phiA);
1873             serializer->doReal(&iparams->pdihs.cpA);
1874             serializer->doReal(&iparams->pdihs.phiB);
1875             serializer->doReal(&iparams->pdihs.cpB);
1876             serializer->doInt(&iparams->pdihs.mult);
1877             break;
1878         case F_RESTRDIHS:
1879             serializer->doReal(&iparams->pdihs.phiA);
1880             serializer->doReal(&iparams->pdihs.cpA);
1881             break;
1882         case F_DISRES:
1883             serializer->doInt(&iparams->disres.label);
1884             serializer->doInt(&iparams->disres.type);
1885             serializer->doReal(&iparams->disres.low);
1886             serializer->doReal(&iparams->disres.up1);
1887             serializer->doReal(&iparams->disres.up2);
1888             serializer->doReal(&iparams->disres.kfac);
1889             break;
1890         case F_ORIRES:
1891             serializer->doInt(&iparams->orires.ex);
1892             serializer->doInt(&iparams->orires.label);
1893             serializer->doInt(&iparams->orires.power);
1894             serializer->doReal(&iparams->orires.c);
1895             serializer->doReal(&iparams->orires.obs);
1896             serializer->doReal(&iparams->orires.kfac);
1897             break;
1898         case F_DIHRES:
1899             if (file_version < 82)
1900             {
1901                 serializer->doInt(&idum);
1902                 serializer->doInt(&idum);
1903             }
1904             serializer->doReal(&iparams->dihres.phiA);
1905             serializer->doReal(&iparams->dihres.dphiA);
1906             serializer->doReal(&iparams->dihres.kfacA);
1907             if (file_version >= 82)
1908             {
1909                 serializer->doReal(&iparams->dihres.phiB);
1910                 serializer->doReal(&iparams->dihres.dphiB);
1911                 serializer->doReal(&iparams->dihres.kfacB);
1912             }
1913             else
1914             {
1915                 iparams->dihres.phiB  = iparams->dihres.phiA;
1916                 iparams->dihres.dphiB = iparams->dihres.dphiA;
1917                 iparams->dihres.kfacB = iparams->dihres.kfacA;
1918             }
1919             break;
1920         case F_POSRES:
1921             serializer->doRvec(&iparams->posres.pos0A);
1922             serializer->doRvec(&iparams->posres.fcA);
1923             serializer->doRvec(&iparams->posres.pos0B);
1924             serializer->doRvec(&iparams->posres.fcB);
1925             break;
1926         case F_FBPOSRES:
1927             serializer->doInt(&iparams->fbposres.geom);
1928             serializer->doRvec(&iparams->fbposres.pos0);
1929             serializer->doReal(&iparams->fbposres.r);
1930             serializer->doReal(&iparams->fbposres.k);
1931             break;
1932         case F_CBTDIHS: serializer->doRealArray(iparams->cbtdihs.cbtcA, NR_CBTDIHS); break;
1933         case F_RBDIHS:
1934             serializer->doRealArray(iparams->rbdihs.rbcA, NR_RBDIHS);
1935             serializer->doRealArray(iparams->rbdihs.rbcB, NR_RBDIHS);
1936             break;
1937         case F_FOURDIHS:
1938             /* Fourier dihedrals are internally represented
1939              * as Ryckaert-Bellemans since those are faster to compute.
1940              */
1941             serializer->doRealArray(iparams->rbdihs.rbcA, NR_RBDIHS);
1942             serializer->doRealArray(iparams->rbdihs.rbcB, NR_RBDIHS);
1943             break;
1944         case F_CONSTR:
1945         case F_CONSTRNC:
1946             serializer->doReal(&iparams->constr.dA);
1947             serializer->doReal(&iparams->constr.dB);
1948             break;
1949         case F_SETTLE:
1950             serializer->doReal(&iparams->settle.doh);
1951             serializer->doReal(&iparams->settle.dhh);
1952             break;
1953         case F_VSITE2:
1954         case F_VSITE2FD: serializer->doReal(&iparams->vsite.a); break;
1955         case F_VSITE3:
1956         case F_VSITE3FD:
1957         case F_VSITE3FAD:
1958             serializer->doReal(&iparams->vsite.a);
1959             serializer->doReal(&iparams->vsite.b);
1960             break;
1961         case F_VSITE3OUT:
1962         case F_VSITE4FD:
1963         case F_VSITE4FDN:
1964             serializer->doReal(&iparams->vsite.a);
1965             serializer->doReal(&iparams->vsite.b);
1966             serializer->doReal(&iparams->vsite.c);
1967             break;
1968         case F_VSITEN:
1969             serializer->doInt(&iparams->vsiten.n);
1970             serializer->doReal(&iparams->vsiten.a);
1971             break;
1972         case F_GB12_NOLONGERUSED:
1973         case F_GB13_NOLONGERUSED:
1974         case F_GB14_NOLONGERUSED:
1975             // Implicit solvent parameters can still be read, but never used
1976             if (serializer->reading())
1977             {
1978                 if (file_version < 68)
1979                 {
1980                     serializer->doReal(&rdum);
1981                     serializer->doReal(&rdum);
1982                     serializer->doReal(&rdum);
1983                     serializer->doReal(&rdum);
1984                 }
1985                 if (file_version < tpxv_RemoveImplicitSolvation)
1986                 {
1987                     serializer->doReal(&rdum);
1988                     serializer->doReal(&rdum);
1989                     serializer->doReal(&rdum);
1990                     serializer->doReal(&rdum);
1991                     serializer->doReal(&rdum);
1992                 }
1993             }
1994             break;
1995         case F_CMAP:
1996             serializer->doInt(&iparams->cmap.cmapA);
1997             serializer->doInt(&iparams->cmap.cmapB);
1998             break;
1999         default:
2000             gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d", ftype,
2001                       interaction_function[ftype].name, __FILE__, __LINE__);
2002     }
2003 }
2004
2005 static void do_ilist(gmx::ISerializer* serializer, InteractionList* ilist)
2006 {
2007     int nr = ilist->size();
2008     serializer->doInt(&nr);
2009     if (serializer->reading())
2010     {
2011         ilist->iatoms.resize(nr);
2012     }
2013     serializer->doIntArray(ilist->iatoms.data(), ilist->size());
2014 }
2015
2016 static void do_ffparams(gmx::ISerializer* serializer, gmx_ffparams_t* ffparams, int file_version)
2017 {
2018     serializer->doInt(&ffparams->atnr);
2019     int numTypes = ffparams->numTypes();
2020     serializer->doInt(&numTypes);
2021     if (serializer->reading())
2022     {
2023         ffparams->functype.resize(numTypes);
2024         ffparams->iparams.resize(numTypes);
2025     }
2026     /* Read/write all the function types */
2027     serializer->doIntArray(ffparams->functype.data(), ffparams->functype.size());
2028
2029     if (file_version >= 66)
2030     {
2031         serializer->doDouble(&ffparams->reppow);
2032     }
2033     else
2034     {
2035         ffparams->reppow = 12.0;
2036     }
2037
2038     serializer->doReal(&ffparams->fudgeQQ);
2039
2040     /* Check whether all these function types are supported by the code.
2041      * In practice the code is backwards compatible, which means that the
2042      * numbering may have to be altered from old numbering to new numbering
2043      */
2044     for (int i = 0; i < ffparams->numTypes(); i++)
2045     {
2046         if (serializer->reading())
2047         {
2048             /* Loop over file versions */
2049             for (int k = 0; k < NFTUPD; k++)
2050             {
2051                 /* Compare the read file_version to the update table */
2052                 if ((file_version < ftupd[k].fvnr) && (ffparams->functype[i] >= ftupd[k].ftype))
2053                 {
2054                     ffparams->functype[i] += 1;
2055                 }
2056             }
2057         }
2058
2059         do_iparams(serializer, ffparams->functype[i], &ffparams->iparams[i], file_version);
2060     }
2061 }
2062
2063 static void add_settle_atoms(InteractionList* ilist)
2064 {
2065     int i;
2066
2067     /* Settle used to only store the first atom: add the other two */
2068     ilist->iatoms.resize(2 * ilist->size());
2069     for (i = ilist->size() / 4 - 1; i >= 0; i--)
2070     {
2071         ilist->iatoms[4 * i + 0] = ilist->iatoms[2 * i + 0];
2072         ilist->iatoms[4 * i + 1] = ilist->iatoms[2 * i + 1];
2073         ilist->iatoms[4 * i + 2] = ilist->iatoms[2 * i + 1] + 1;
2074         ilist->iatoms[4 * i + 3] = ilist->iatoms[2 * i + 1] + 2;
2075     }
2076 }
2077
2078 static void do_ilists(gmx::ISerializer* serializer, InteractionLists* ilists, int file_version)
2079 {
2080     GMX_RELEASE_ASSERT(ilists, "Need a valid ilists object");
2081     GMX_RELEASE_ASSERT(ilists->size() == F_NRE,
2082                        "The code needs to be in sync with InteractionLists");
2083
2084     for (int j = 0; j < F_NRE; j++)
2085     {
2086         InteractionList& ilist  = (*ilists)[j];
2087         gmx_bool         bClear = FALSE;
2088         if (serializer->reading())
2089         {
2090             for (int k = 0; k < NFTUPD; k++)
2091             {
2092                 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2093                 {
2094                     bClear = TRUE;
2095                 }
2096             }
2097         }
2098         if (bClear)
2099         {
2100             ilist.iatoms.clear();
2101         }
2102         else
2103         {
2104             do_ilist(serializer, &ilist);
2105             if (file_version < 78 && j == F_SETTLE && ilist.size() > 0)
2106             {
2107                 add_settle_atoms(&ilist);
2108             }
2109         }
2110     }
2111 }
2112
2113 static void do_block(gmx::ISerializer* serializer, t_block* block)
2114 {
2115     serializer->doInt(&block->nr);
2116     if (serializer->reading())
2117     {
2118         if ((block->nalloc_index > 0) && (nullptr != block->index))
2119         {
2120             sfree(block->index);
2121         }
2122         block->nalloc_index = block->nr + 1;
2123         snew(block->index, block->nalloc_index);
2124     }
2125     serializer->doIntArray(block->index, block->nr + 1);
2126 }
2127
2128 static void doListOfLists(gmx::ISerializer* serializer, gmx::ListOfLists<int>* listOfLists)
2129 {
2130     int numLists = listOfLists->ssize();
2131     serializer->doInt(&numLists);
2132     int numElements = listOfLists->elementsView().ssize();
2133     serializer->doInt(&numElements);
2134     if (serializer->reading())
2135     {
2136         std::vector<int> listRanges(numLists + 1);
2137         serializer->doIntArray(listRanges.data(), numLists + 1);
2138         std::vector<int> elements(numElements);
2139         serializer->doIntArray(elements.data(), numElements);
2140         *listOfLists = gmx::ListOfLists<int>(std::move(listRanges), std::move(elements));
2141     }
2142     else
2143     {
2144         serializer->doIntArray(const_cast<int*>(listOfLists->listRangesView().data()), numLists + 1);
2145         serializer->doIntArray(const_cast<int*>(listOfLists->elementsView().data()), numElements);
2146     }
2147 }
2148
2149 /* This is a primitive routine to make it possible to translate atomic numbers
2150  * to element names when reading TPR files, without making the Gromacs library
2151  * directory a dependency on mdrun (which is the case if we need elements.dat).
2152  */
2153 static const char* atomicnumber_to_element(int atomicnumber)
2154 {
2155     const char* p;
2156
2157     /* This does not have to be complete, so we only include elements likely
2158      * to occur in PDB files.
2159      */
2160     switch (atomicnumber)
2161     {
2162         case 1: p = "H"; break;
2163         case 5: p = "B"; break;
2164         case 6: p = "C"; break;
2165         case 7: p = "N"; break;
2166         case 8: p = "O"; break;
2167         case 9: p = "F"; break;
2168         case 11: p = "Na"; break;
2169         case 12: p = "Mg"; break;
2170         case 15: p = "P"; break;
2171         case 16: p = "S"; break;
2172         case 17: p = "Cl"; break;
2173         case 18: p = "Ar"; break;
2174         case 19: p = "K"; break;
2175         case 20: p = "Ca"; break;
2176         case 25: p = "Mn"; break;
2177         case 26: p = "Fe"; break;
2178         case 28: p = "Ni"; break;
2179         case 29: p = "Cu"; break;
2180         case 30: p = "Zn"; break;
2181         case 35: p = "Br"; break;
2182         case 47: p = "Ag"; break;
2183         default: p = ""; break;
2184     }
2185     return p;
2186 }
2187
2188
2189 static void do_atom(gmx::ISerializer* serializer, t_atom* atom)
2190 {
2191     serializer->doReal(&atom->m);
2192     serializer->doReal(&atom->q);
2193     serializer->doReal(&atom->mB);
2194     serializer->doReal(&atom->qB);
2195     serializer->doUShort(&atom->type);
2196     serializer->doUShort(&atom->typeB);
2197     serializer->doInt(&atom->ptype);
2198     serializer->doInt(&atom->resind);
2199     serializer->doInt(&atom->atomnumber);
2200     if (serializer->reading())
2201     {
2202         /* Set element string from atomic number if present.
2203          * This routine returns an empty string if the name is not found.
2204          */
2205         std::strncpy(atom->elem, atomicnumber_to_element(atom->atomnumber), 4);
2206         /* avoid warnings about potentially unterminated string */
2207         atom->elem[3] = '\0';
2208     }
2209 }
2210
2211 static void do_grps(gmx::ISerializer* serializer, gmx::ArrayRef<AtomGroupIndices> grps)
2212 {
2213     for (auto& group : grps)
2214     {
2215         int size = group.size();
2216         serializer->doInt(&size);
2217         if (serializer->reading())
2218         {
2219             group.resize(size);
2220         }
2221         serializer->doIntArray(group.data(), size);
2222     }
2223 }
2224
2225 static void do_symstr(gmx::ISerializer* serializer, char*** nm, t_symtab* symtab)
2226 {
2227     int ls;
2228
2229     if (serializer->reading())
2230     {
2231         serializer->doInt(&ls);
2232         *nm = get_symtab_handle(symtab, ls);
2233     }
2234     else
2235     {
2236         ls = lookup_symtab(symtab, *nm);
2237         serializer->doInt(&ls);
2238     }
2239 }
2240
2241 static void do_strstr(gmx::ISerializer* serializer, int nstr, char*** nm, t_symtab* symtab)
2242 {
2243     int j;
2244
2245     for (j = 0; (j < nstr); j++)
2246     {
2247         do_symstr(serializer, &(nm[j]), symtab);
2248     }
2249 }
2250
2251 static void do_resinfo(gmx::ISerializer* serializer, int n, t_resinfo* ri, t_symtab* symtab, int file_version)
2252 {
2253     int j;
2254
2255     for (j = 0; (j < n); j++)
2256     {
2257         do_symstr(serializer, &(ri[j].name), symtab);
2258         if (file_version >= 63)
2259         {
2260             serializer->doInt(&ri[j].nr);
2261             serializer->doUChar(&ri[j].ic);
2262         }
2263         else
2264         {
2265             ri[j].nr = j + 1;
2266             ri[j].ic = ' ';
2267         }
2268     }
2269 }
2270
2271 static void do_atoms(gmx::ISerializer* serializer, t_atoms* atoms, t_symtab* symtab, int file_version)
2272 {
2273     int i;
2274
2275     serializer->doInt(&atoms->nr);
2276     serializer->doInt(&atoms->nres);
2277     if (serializer->reading())
2278     {
2279         /* Since we have always written all t_atom properties in the tpr file
2280          * (at least for all backward compatible versions), we don't store
2281          * but simple set the booleans here.
2282          */
2283         atoms->haveMass    = TRUE;
2284         atoms->haveCharge  = TRUE;
2285         atoms->haveType    = TRUE;
2286         atoms->haveBState  = TRUE;
2287         atoms->havePdbInfo = FALSE;
2288
2289         snew(atoms->atom, atoms->nr);
2290         snew(atoms->atomname, atoms->nr);
2291         snew(atoms->atomtype, atoms->nr);
2292         snew(atoms->atomtypeB, atoms->nr);
2293         snew(atoms->resinfo, atoms->nres);
2294         atoms->pdbinfo = nullptr;
2295     }
2296     else
2297     {
2298         GMX_RELEASE_ASSERT(atoms->haveMass && atoms->haveCharge && atoms->haveType && atoms->haveBState,
2299                            "Mass, charge, atomtype and B-state parameters should be present in "
2300                            "t_atoms when writing a tpr file");
2301     }
2302     for (i = 0; (i < atoms->nr); i++)
2303     {
2304         do_atom(serializer, &atoms->atom[i]);
2305     }
2306     do_strstr(serializer, atoms->nr, atoms->atomname, symtab);
2307     do_strstr(serializer, atoms->nr, atoms->atomtype, symtab);
2308     do_strstr(serializer, atoms->nr, atoms->atomtypeB, symtab);
2309
2310     do_resinfo(serializer, atoms->nres, atoms->resinfo, symtab, file_version);
2311 }
2312
2313 static void do_groups(gmx::ISerializer* serializer, SimulationGroups* groups, t_symtab* symtab)
2314 {
2315     do_grps(serializer, groups->groups);
2316     int numberOfGroupNames = groups->groupNames.size();
2317     serializer->doInt(&numberOfGroupNames);
2318     if (serializer->reading())
2319     {
2320         groups->groupNames.resize(numberOfGroupNames);
2321     }
2322     do_strstr(serializer, numberOfGroupNames, groups->groupNames.data(), symtab);
2323     for (auto group : gmx::keysOf(groups->groupNumbers))
2324     {
2325         int numberOfGroupNumbers = groups->numberOfGroupNumbers(group);
2326         serializer->doInt(&numberOfGroupNumbers);
2327         if (numberOfGroupNumbers != 0)
2328         {
2329             if (serializer->reading())
2330             {
2331                 groups->groupNumbers[group].resize(numberOfGroupNumbers);
2332             }
2333             serializer->doUCharArray(groups->groupNumbers[group].data(), numberOfGroupNumbers);
2334         }
2335     }
2336 }
2337
2338 static void do_atomtypes(gmx::ISerializer* serializer, t_atomtypes* atomtypes, int file_version)
2339 {
2340     int j;
2341
2342     serializer->doInt(&atomtypes->nr);
2343     j = atomtypes->nr;
2344     if (serializer->reading())
2345     {
2346         snew(atomtypes->atomnumber, j);
2347     }
2348     if (serializer->reading() && file_version < tpxv_RemoveImplicitSolvation)
2349     {
2350         std::vector<real> dummy(atomtypes->nr, 0);
2351         serializer->doRealArray(dummy.data(), dummy.size());
2352         serializer->doRealArray(dummy.data(), dummy.size());
2353         serializer->doRealArray(dummy.data(), dummy.size());
2354     }
2355     serializer->doIntArray(atomtypes->atomnumber, j);
2356
2357     if (serializer->reading() && file_version >= 60 && file_version < tpxv_RemoveImplicitSolvation)
2358     {
2359         std::vector<real> dummy(atomtypes->nr, 0);
2360         serializer->doRealArray(dummy.data(), dummy.size());
2361         serializer->doRealArray(dummy.data(), dummy.size());
2362     }
2363 }
2364
2365 static void do_symtab(gmx::ISerializer* serializer, t_symtab* symtab)
2366 {
2367     int       i, nr;
2368     t_symbuf* symbuf;
2369
2370     serializer->doInt(&symtab->nr);
2371     nr = symtab->nr;
2372     if (serializer->reading())
2373     {
2374         snew(symtab->symbuf, 1);
2375         symbuf          = symtab->symbuf;
2376         symbuf->bufsize = nr;
2377         snew(symbuf->buf, nr);
2378         for (i = 0; (i < nr); i++)
2379         {
2380             std::string buf;
2381             serializer->doString(&buf);
2382             symbuf->buf[i] = gmx_strdup(buf.c_str());
2383         }
2384     }
2385     else
2386     {
2387         symbuf = symtab->symbuf;
2388         while (symbuf != nullptr)
2389         {
2390             for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2391             {
2392                 std::string buf = symbuf->buf[i];
2393                 serializer->doString(&buf);
2394             }
2395             nr -= i;
2396             symbuf = symbuf->next;
2397         }
2398         if (nr != 0)
2399         {
2400             gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2401         }
2402     }
2403 }
2404
2405 static void do_cmap(gmx::ISerializer* serializer, gmx_cmap_t* cmap_grid)
2406 {
2407
2408     int ngrid = cmap_grid->cmapdata.size();
2409     serializer->doInt(&ngrid);
2410     serializer->doInt(&cmap_grid->grid_spacing);
2411
2412     int gs    = cmap_grid->grid_spacing;
2413     int nelem = gs * gs;
2414
2415     if (serializer->reading())
2416     {
2417         cmap_grid->cmapdata.resize(ngrid);
2418
2419         for (int i = 0; i < ngrid; i++)
2420         {
2421             cmap_grid->cmapdata[i].cmap.resize(4 * nelem);
2422         }
2423     }
2424
2425     for (int i = 0; i < ngrid; i++)
2426     {
2427         for (int j = 0; j < nelem; j++)
2428         {
2429             serializer->doReal(&cmap_grid->cmapdata[i].cmap[j * 4]);
2430             serializer->doReal(&cmap_grid->cmapdata[i].cmap[j * 4 + 1]);
2431             serializer->doReal(&cmap_grid->cmapdata[i].cmap[j * 4 + 2]);
2432             serializer->doReal(&cmap_grid->cmapdata[i].cmap[j * 4 + 3]);
2433         }
2434     }
2435 }
2436
2437
2438 static void do_moltype(gmx::ISerializer* serializer, gmx_moltype_t* molt, t_symtab* symtab, int file_version)
2439 {
2440     do_symstr(serializer, &(molt->name), symtab);
2441
2442     do_atoms(serializer, &molt->atoms, symtab, file_version);
2443
2444     do_ilists(serializer, &molt->ilist, file_version);
2445
2446     /* TODO: Remove the obsolete charge group index from the file */
2447     t_block cgs;
2448     cgs.nr           = molt->atoms.nr;
2449     cgs.nalloc_index = cgs.nr + 1;
2450     snew(cgs.index, cgs.nalloc_index);
2451     for (int i = 0; i < cgs.nr + 1; i++)
2452     {
2453         cgs.index[i] = i;
2454     }
2455     do_block(serializer, &cgs);
2456     sfree(cgs.index);
2457
2458     /* This used to be in the atoms struct */
2459     doListOfLists(serializer, &molt->excls);
2460 }
2461
2462 static void do_molblock(gmx::ISerializer* serializer, gmx_molblock_t* molb, int numAtomsPerMolecule)
2463 {
2464     serializer->doInt(&molb->type);
2465     serializer->doInt(&molb->nmol);
2466     /* To maintain forward topology reading compatibility, we store #atoms.
2467      * TODO: Change this to conditional reading of a dummy int when we
2468      *       increase tpx_generation.
2469      */
2470     serializer->doInt(&numAtomsPerMolecule);
2471     /* Position restraint coordinates */
2472     int numPosres_xA = molb->posres_xA.size();
2473     serializer->doInt(&numPosres_xA);
2474     if (numPosres_xA > 0)
2475     {
2476         if (serializer->reading())
2477         {
2478             molb->posres_xA.resize(numPosres_xA);
2479         }
2480         serializer->doRvecArray(as_rvec_array(molb->posres_xA.data()), numPosres_xA);
2481     }
2482     int numPosres_xB = molb->posres_xB.size();
2483     serializer->doInt(&numPosres_xB);
2484     if (numPosres_xB > 0)
2485     {
2486         if (serializer->reading())
2487         {
2488             molb->posres_xB.resize(numPosres_xB);
2489         }
2490         serializer->doRvecArray(as_rvec_array(molb->posres_xB.data()), numPosres_xB);
2491     }
2492 }
2493
2494 static void set_disres_npair(gmx_mtop_t* mtop)
2495 {
2496     gmx_mtop_ilistloop_t iloop;
2497     int                  nmol;
2498
2499     gmx::ArrayRef<t_iparams> ip = mtop->ffparams.iparams;
2500
2501     iloop = gmx_mtop_ilistloop_init(mtop);
2502     while (const InteractionLists* ilist = gmx_mtop_ilistloop_next(iloop, &nmol))
2503     {
2504         const InteractionList& il = (*ilist)[F_DISRES];
2505
2506         if (il.size() > 0)
2507         {
2508             gmx::ArrayRef<const int> a     = il.iatoms;
2509             int                      npair = 0;
2510             for (int i = 0; i < il.size(); i += 3)
2511             {
2512                 npair++;
2513                 if (i + 3 == il.size() || ip[a[i]].disres.label != ip[a[i + 3]].disres.label)
2514                 {
2515                     ip[a[i]].disres.npair = npair;
2516                     npair                 = 0;
2517                 }
2518             }
2519         }
2520     }
2521 }
2522
2523 static void do_mtop(gmx::ISerializer* serializer, gmx_mtop_t* mtop, int file_version)
2524 {
2525     do_symtab(serializer, &(mtop->symtab));
2526
2527     do_symstr(serializer, &(mtop->name), &(mtop->symtab));
2528
2529     do_ffparams(serializer, &mtop->ffparams, file_version);
2530
2531     int nmoltype = mtop->moltype.size();
2532     serializer->doInt(&nmoltype);
2533     if (serializer->reading())
2534     {
2535         mtop->moltype.resize(nmoltype);
2536     }
2537     for (gmx_moltype_t& moltype : mtop->moltype)
2538     {
2539         do_moltype(serializer, &moltype, &mtop->symtab, file_version);
2540     }
2541
2542     int nmolblock = mtop->molblock.size();
2543     serializer->doInt(&nmolblock);
2544     if (serializer->reading())
2545     {
2546         mtop->molblock.resize(nmolblock);
2547     }
2548     for (gmx_molblock_t& molblock : mtop->molblock)
2549     {
2550         int numAtomsPerMolecule = (serializer->reading() ? 0 : mtop->moltype[molblock.type].atoms.nr);
2551         do_molblock(serializer, &molblock, numAtomsPerMolecule);
2552     }
2553     serializer->doInt(&mtop->natoms);
2554
2555     if (file_version >= tpxv_IntermolecularBondeds)
2556     {
2557         serializer->doBool(&mtop->bIntermolecularInteractions);
2558         if (mtop->bIntermolecularInteractions)
2559         {
2560             if (serializer->reading())
2561             {
2562                 mtop->intermolecular_ilist = std::make_unique<InteractionLists>();
2563             }
2564             do_ilists(serializer, mtop->intermolecular_ilist.get(), file_version);
2565         }
2566     }
2567     else
2568     {
2569         mtop->bIntermolecularInteractions = FALSE;
2570     }
2571
2572     do_atomtypes(serializer, &(mtop->atomtypes), file_version);
2573
2574     if (file_version >= 65)
2575     {
2576         do_cmap(serializer, &mtop->ffparams.cmap_grid);
2577     }
2578     else
2579     {
2580         mtop->ffparams.cmap_grid.grid_spacing = 0;
2581         mtop->ffparams.cmap_grid.cmapdata.clear();
2582     }
2583
2584     do_groups(serializer, &mtop->groups, &(mtop->symtab));
2585
2586     mtop->haveMoleculeIndices = true;
2587
2588     if (serializer->reading())
2589     {
2590         close_symtab(&(mtop->symtab));
2591     }
2592 }
2593
2594 /*! \brief
2595  * Read the first part of the TPR file to find general system information.
2596  *
2597  * If \p TopOnlyOK is true then we can read even future versions
2598  * of tpx files, provided the \p fileGeneration hasn't changed.
2599  * If it is false, we need the \p ir too, and bail out
2600  * if the file is newer than the program.
2601  *
2602  * The version and generation of the topology (see top of this file)
2603  * are returned in the two last arguments, if those arguments are non-nullptr.
2604  *
2605  * If possible, we will read the \p ir even when \p TopOnlyOK is true.
2606  *
2607  * \param[in,out] serializer The serializer used to handle header processing.
2608  * \param[in,out] tpx File header datastructure.
2609  * \param[in]     filename The name of the file being read/written
2610  * \param[in,out] fio File handle.
2611  * \param[in] TopOnlyOK If not reading \p ir is fine or not.
2612  */
2613 static void do_tpxheader(gmx::FileIOXdrSerializer* serializer,
2614                          TpxFileHeader*            tpx,
2615                          const char*               filename,
2616                          t_fileio*                 fio,
2617                          bool                      TopOnlyOK)
2618 {
2619     int  precision;
2620     int  idum = 0;
2621     real rdum = 0;
2622
2623     /* XDR binary topology file */
2624     precision = sizeof(real);
2625     std::string buf;
2626     std::string fileTag;
2627     if (serializer->reading())
2628     {
2629         serializer->doString(&buf);
2630         if (std::strncmp(buf.c_str(), "VERSION", 7) != 0)
2631         {
2632             gmx_fatal(
2633                     FARGS,
2634                     "Can not read file %s,\n"
2635                     "             this file is from a GROMACS version which is older than 2.0\n"
2636                     "             Make a new one with grompp or use a gro or pdb file, if possible",
2637                     filename);
2638         }
2639         // We need to know the precision used to write the TPR file, to match it
2640         // to the precision of the currently running binary. If the precisions match
2641         // there is no problem, but mismatching precision needs to be accounted for
2642         // by reading into temporary variables of the correct precision instead
2643         // of the desired target datastructures.
2644         serializer->doInt(&precision);
2645         tpx->isDouble = (precision == sizeof(double));
2646         if ((precision != sizeof(float)) && !tpx->isDouble)
2647         {
2648             gmx_fatal(FARGS,
2649                       "Unknown precision in file %s: real is %d bytes "
2650                       "instead of %zu or %zu",
2651                       filename, precision, sizeof(float), sizeof(double));
2652         }
2653         gmx_fio_setprecision(fio, tpx->isDouble);
2654         fprintf(stderr, "Reading file %s, %s (%s precision)\n", filename, buf.c_str(),
2655                 tpx->isDouble ? "double" : "single");
2656     }
2657     else
2658     {
2659         buf = gmx::formatString("VERSION %s", gmx_version());
2660         serializer->doString(&buf);
2661         gmx_fio_setprecision(fio, tpx->isDouble);
2662         serializer->doInt(&precision);
2663         fileTag = gmx::formatString("%s", tpx_tag);
2664     }
2665
2666     /* Check versions! */
2667     serializer->doInt(&tpx->fileVersion);
2668
2669     /* This is for backward compatibility with development versions 77-79
2670      * where the tag was, mistakenly, placed before the generation,
2671      * which would cause a segv instead of a proper error message
2672      * when reading the topology only from tpx with <77 code.
2673      */
2674     if (tpx->fileVersion >= 77 && tpx->fileVersion <= 79)
2675     {
2676         serializer->doString(&fileTag);
2677     }
2678
2679     serializer->doInt(&tpx->fileGeneration);
2680
2681     if (tpx->fileVersion >= 81)
2682     {
2683         serializer->doString(&fileTag);
2684     }
2685     if (serializer->reading())
2686     {
2687         if (tpx->fileVersion < 77)
2688         {
2689             /* Versions before 77 don't have the tag, set it to release */
2690             fileTag = gmx::formatString("%s", TPX_TAG_RELEASE);
2691         }
2692
2693         if (fileTag != tpx_tag)
2694         {
2695             fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n", fileTag.c_str(), tpx_tag);
2696
2697             /* We only support reading tpx files with the same tag as the code
2698              * or tpx files with the release tag and with lower version number.
2699              */
2700             if (fileTag != TPX_TAG_RELEASE && tpx->fileVersion < tpx_version)
2701             {
2702                 gmx_fatal(FARGS,
2703                           "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' "
2704                           "with program for tpx version %d, tag '%s'",
2705                           filename, tpx->fileVersion, fileTag.c_str(), tpx_version, tpx_tag);
2706             }
2707         }
2708     }
2709
2710     if ((tpx->fileVersion <= tpx_incompatible_version)
2711         || ((tpx->fileVersion > tpx_version) && !TopOnlyOK) || (tpx->fileGeneration > tpx_generation)
2712         || tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
2713     {
2714         gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program", filename,
2715                   tpx->fileVersion, tpx_version);
2716     }
2717
2718     serializer->doInt(&tpx->natoms);
2719     serializer->doInt(&tpx->ngtc);
2720
2721     if (tpx->fileVersion < 62)
2722     {
2723         serializer->doInt(&idum);
2724         serializer->doReal(&rdum);
2725     }
2726     if (tpx->fileVersion >= 79)
2727     {
2728         serializer->doInt(&tpx->fep_state);
2729     }
2730     serializer->doReal(&tpx->lambda);
2731     serializer->doBool(&tpx->bIr);
2732     serializer->doBool(&tpx->bTop);
2733     serializer->doBool(&tpx->bX);
2734     serializer->doBool(&tpx->bV);
2735     serializer->doBool(&tpx->bF);
2736     serializer->doBool(&tpx->bBox);
2737
2738     if (tpx->fileVersion >= tpxv_AddSizeField && tpx->fileGeneration >= 27)
2739     {
2740         if (!serializer->reading())
2741         {
2742             GMX_RELEASE_ASSERT(tpx->sizeOfTprBody != 0,
2743                                "Not possible to write new file with zero TPR body size");
2744         }
2745         serializer->doInt64(&tpx->sizeOfTprBody);
2746     }
2747
2748     if ((tpx->fileGeneration > tpx_generation))
2749     {
2750         /* This can only happen if TopOnlyOK=TRUE */
2751         tpx->bIr = FALSE;
2752     }
2753 }
2754
2755 #define do_test(serializer, b, p)                            \
2756     if ((serializer)->reading() && ((p) != nullptr) && !(b)) \
2757     gmx_fatal(FARGS, "No %s in input file", #p)
2758
2759 /*! \brief
2760  * Process the first part of the TPR into the state datastructure.
2761  *
2762  * Due to the structure of the legacy code, it is necessary
2763  * to split up the state reading into two parts, with the
2764  * box and legacy temperature coupling processed before the
2765  * topology datastructures.
2766  *
2767  * See the documentation for do_tpx_body for the correct order of
2768  * the operations for reading a tpr file.
2769  *
2770  * \param[in] serializer Abstract serializer used to read/write data.
2771  * \param[in] tpx The file header data.
2772  * \param[in, out] state Global state data.
2773  */
2774 static void do_tpx_state_first(gmx::ISerializer* serializer, TpxFileHeader* tpx, t_state* state)
2775 {
2776     if (serializer->reading())
2777     {
2778         state->flags = 0;
2779         init_gtc_state(state, tpx->ngtc, 0, 0);
2780     }
2781     do_test(serializer, tpx->bBox, state->box);
2782     if (tpx->bBox)
2783     {
2784         serializer->doRvecArray(state->box, DIM);
2785         if (tpx->fileVersion >= 51)
2786         {
2787             serializer->doRvecArray(state->box_rel, DIM);
2788         }
2789         else
2790         {
2791             /* We initialize box_rel after reading the inputrec */
2792             clear_mat(state->box_rel);
2793         }
2794         serializer->doRvecArray(state->boxv, DIM);
2795         if (tpx->fileVersion < 56)
2796         {
2797             matrix mdum;
2798             serializer->doRvecArray(mdum, DIM);
2799         }
2800     }
2801
2802     if (state->ngtc > 0)
2803     {
2804         real* dumv;
2805         snew(dumv, state->ngtc);
2806         if (tpx->fileVersion < 69)
2807         {
2808             serializer->doRealArray(dumv, state->ngtc);
2809         }
2810         /* These used to be the Berendsen tcoupl_lambda's */
2811         serializer->doRealArray(dumv, state->ngtc);
2812         sfree(dumv);
2813     }
2814 }
2815
2816 /*! \brief
2817  * Process global topology data.
2818  *
2819  * See the documentation for do_tpx_body for the correct order of
2820  * the operations for reading a tpr file.
2821  *
2822  * \param[in] serializer Abstract serializer  used to read/write data.
2823  * \param[in] tpx The file header data.
2824  * \param[in,out] mtop Global topology.
2825  */
2826 static void do_tpx_mtop(gmx::ISerializer* serializer, TpxFileHeader* tpx, gmx_mtop_t* mtop)
2827 {
2828     do_test(serializer, tpx->bTop, mtop);
2829     if (tpx->bTop)
2830     {
2831         if (mtop)
2832         {
2833             do_mtop(serializer, mtop, tpx->fileVersion);
2834             set_disres_npair(mtop);
2835             gmx_mtop_finalize(mtop);
2836         }
2837         else
2838         {
2839             gmx_mtop_t dum_top;
2840             do_mtop(serializer, &dum_top, tpx->fileVersion);
2841         }
2842     }
2843 }
2844 /*! \brief
2845  * Process coordinate vectors for state data.
2846  *
2847  * Main part of state gets processed here.
2848  *
2849  * See the documentation for do_tpx_body for the correct order of
2850  * the operations for reading a tpr file.
2851  *
2852  * \param[in] serializer Abstract serializer used to read/write data.
2853  * \param[in] tpx The file header data.
2854  * \param[in,out] state Global state data.
2855  * \param[in,out] x Individual coordinates for processing, deprecated.
2856  * \param[in,out] v Individual velocities for processing, deprecated.
2857  */
2858 static void do_tpx_state_second(gmx::ISerializer* serializer, TpxFileHeader* tpx, t_state* state, rvec* x, rvec* v)
2859 {
2860     if (!serializer->reading())
2861     {
2862         GMX_RELEASE_ASSERT(
2863                 x == nullptr && v == nullptr,
2864                 "Passing separate x and v pointers to do_tpx() is not supported when writing");
2865     }
2866     else
2867     {
2868         GMX_RELEASE_ASSERT(!(x == nullptr && v != nullptr),
2869                            "Passing x==NULL and v!=NULL is not supported");
2870     }
2871
2872     if (serializer->reading())
2873     {
2874         if (x == nullptr)
2875         {
2876             // v is also nullptr by the above assertion, so we may
2877             // need to make memory in state for storing the contents
2878             // of the tpx file.
2879             if (tpx->bX)
2880             {
2881                 state->flags |= (1 << estX);
2882             }
2883             if (tpx->bV)
2884             {
2885                 state->flags |= (1 << estV);
2886             }
2887             state_change_natoms(state, tpx->natoms);
2888         }
2889     }
2890
2891     if (x == nullptr)
2892     {
2893         x = state->x.rvec_array();
2894         v = state->v.rvec_array();
2895     }
2896     do_test(serializer, tpx->bX, x);
2897     if (tpx->bX)
2898     {
2899         if (serializer->reading())
2900         {
2901             state->flags |= (1 << estX);
2902         }
2903         serializer->doRvecArray(x, tpx->natoms);
2904     }
2905
2906     do_test(serializer, tpx->bV, v);
2907     if (tpx->bV)
2908     {
2909         if (serializer->reading())
2910         {
2911             state->flags |= (1 << estV);
2912         }
2913         if (!v)
2914         {
2915             std::vector<gmx::RVec> dummyVelocities(tpx->natoms);
2916             serializer->doRvecArray(as_rvec_array(dummyVelocities.data()), tpx->natoms);
2917         }
2918         else
2919         {
2920             serializer->doRvecArray(v, tpx->natoms);
2921         }
2922     }
2923
2924     // No need to run do_test when the last argument is NULL
2925     if (tpx->bF)
2926     {
2927         std::vector<gmx::RVec> dummyForces(state->natoms);
2928         serializer->doRvecArray(as_rvec_array(dummyForces.data()), tpx->natoms);
2929     }
2930 }
2931 /*! \brief
2932  * Process simulation parameters.
2933  *
2934  * See the documentation for do_tpx_body for the correct order of
2935  * the operations for reading a tpr file.
2936  *
2937  * \param[in] serializer Abstract serializer used to read/write data.
2938  * \param[in] tpx The file header data.
2939  * \param[in,out] ir Datastructure with simulation parameters.
2940  */
2941 static int do_tpx_ir(gmx::ISerializer* serializer, TpxFileHeader* tpx, t_inputrec* ir)
2942 {
2943     int      ePBC;
2944     gmx_bool bPeriodicMols;
2945
2946     /* Starting with tpx version 26, we have the inputrec
2947      * at the end of the file, so we can ignore it
2948      * if the file is never than the software (but still the
2949      * same generation - see comments at the top of this file.
2950      *
2951      *
2952      */
2953     ePBC          = -1;
2954     bPeriodicMols = FALSE;
2955
2956     do_test(serializer, tpx->bIr, ir);
2957     if (tpx->bIr)
2958     {
2959         if (tpx->fileVersion >= 53)
2960         {
2961             /* Removed the pbc info from do_inputrec, since we always want it */
2962             if (!serializer->reading())
2963             {
2964                 ePBC          = ir->ePBC;
2965                 bPeriodicMols = ir->bPeriodicMols;
2966             }
2967             serializer->doInt(&ePBC);
2968             serializer->doBool(&bPeriodicMols);
2969         }
2970         if (tpx->fileGeneration <= tpx_generation && ir)
2971         {
2972             do_inputrec(serializer, ir, tpx->fileVersion);
2973             if (tpx->fileVersion < 53)
2974             {
2975                 ePBC          = ir->ePBC;
2976                 bPeriodicMols = ir->bPeriodicMols;
2977             }
2978         }
2979         if (serializer->reading() && ir && tpx->fileVersion >= 53)
2980         {
2981             /* We need to do this after do_inputrec, since that initializes ir */
2982             ir->ePBC          = ePBC;
2983             ir->bPeriodicMols = bPeriodicMols;
2984         }
2985     }
2986     return ePBC;
2987 }
2988
2989 /*! \brief
2990  * Correct and finalize read information.
2991  *
2992  * If \p state is nullptr, skip the parts dependent on it.
2993  *
2994  * See the documentation for do_tpx_body for the correct order of
2995  * the operations for reading a tpr file.
2996  *
2997  * \param[in] tpx The file header used to check version numbers.
2998  * \param[out] ir Input rec that needs correction.
2999  * \param[out] state State needing correction.
3000  * \param[out] mtop Topology to finalize.
3001  */
3002 static void do_tpx_finalize(TpxFileHeader* tpx, t_inputrec* ir, t_state* state, gmx_mtop_t* mtop)
3003 {
3004     if (tpx->fileVersion < 51 && state)
3005     {
3006         set_box_rel(ir, state);
3007     }
3008     if (tpx->bIr && ir)
3009     {
3010         if (state && state->ngtc == 0)
3011         {
3012             /* Reading old version without tcoupl state data: set it */
3013             init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3014         }
3015         if (tpx->bTop && mtop)
3016         {
3017             if (tpx->fileVersion < 57)
3018             {
3019                 if (mtop->moltype[0].ilist[F_DISRES].size() > 0)
3020                 {
3021                     ir->eDisre = edrSimple;
3022                 }
3023                 else
3024                 {
3025                     ir->eDisre = edrNone;
3026                 }
3027             }
3028         }
3029     }
3030 }
3031
3032 /*! \brief
3033  * Process TPR data for file reading/writing.
3034  *
3035  * The TPR file gets processed in in four stages due to the organization
3036  * of the data within it.
3037  *
3038  * First, state data for the box is processed in do_tpx_state_first.
3039  * This is followed by processing the topology in do_tpx_mtop.
3040  * Coordinate and velocity vectors are handled next in do_tpx_state_second.
3041  * The last file information processed is the collection of simulation parameters in do_tpx_ir.
3042  * When reading, a final processing step is undertaken at the end.
3043  *
3044  * \param[in] serializer Abstract serializer used to read/write data.
3045  * \param[in] tpx The file header data.
3046  * \param[in,out] ir Datastructures with simulation parameters.
3047  * \param[in,out] state Global state data.
3048  * \param[in,out] x Individual coordinates for processing, deprecated.
3049  * \param[in,out] v Individual velocities for processing, deprecated.
3050  * \param[in,out] mtop Global topology.
3051  */
3052 static int do_tpx_body(gmx::ISerializer* serializer,
3053                        TpxFileHeader*    tpx,
3054                        t_inputrec*       ir,
3055                        t_state*          state,
3056                        rvec*             x,
3057                        rvec*             v,
3058                        gmx_mtop_t*       mtop)
3059 {
3060     if (state)
3061     {
3062         do_tpx_state_first(serializer, tpx, state);
3063     }
3064     do_tpx_mtop(serializer, tpx, mtop);
3065     if (state)
3066     {
3067         do_tpx_state_second(serializer, tpx, state, x, v);
3068     }
3069     int ePBC = do_tpx_ir(serializer, tpx, ir);
3070     if (serializer->reading())
3071     {
3072         do_tpx_finalize(tpx, ir, state, mtop);
3073     }
3074     return ePBC;
3075 }
3076
3077 /*! \brief
3078  * Overload for do_tpx_body that defaults to state vectors being nullptr.
3079  *
3080  * \param[in] serializer Abstract serializer used to read/write data.
3081  * \param[in] tpx The file header data.
3082  * \param[in,out] ir Datastructures with simulation parameters.
3083  * \param[in,out] mtop Global topology.
3084  */
3085 static int do_tpx_body(gmx::ISerializer* serializer, TpxFileHeader* tpx, t_inputrec* ir, gmx_mtop_t* mtop)
3086 {
3087     return do_tpx_body(serializer, tpx, ir, nullptr, nullptr, nullptr, mtop);
3088 }
3089
3090 static t_fileio* open_tpx(const char* fn, const char* mode)
3091 {
3092     return gmx_fio_open(fn, mode);
3093 }
3094
3095 static void close_tpx(t_fileio* fio)
3096 {
3097     gmx_fio_close(fio);
3098 }
3099
3100 /*! \brief
3101  * Fill information into the header only from state before writing.
3102  *
3103  * Populating the header needs to be independent from writing the information
3104  * to file to allow things like writing the raw byte stream.
3105  *
3106  * \param[in] state The current simulation state. Can't write without it.
3107  * \param[in] ir Parameter and system information.
3108  * \param[in] mtop Global topology.
3109  * \returns Fully populated header.
3110  */
3111 static TpxFileHeader populateTpxHeader(const t_state& state, const t_inputrec* ir, const gmx_mtop_t* mtop)
3112 {
3113     TpxFileHeader header;
3114     header.natoms         = state.natoms;
3115     header.ngtc           = state.ngtc;
3116     header.fep_state      = state.fep_state;
3117     header.lambda         = state.lambda[efptFEP];
3118     header.bIr            = ir != nullptr;
3119     header.bTop           = mtop != nullptr;
3120     header.bX             = (state.flags & (1 << estX)) != 0;
3121     header.bV             = (state.flags & (1 << estV)) != 0;
3122     header.bF             = false;
3123     header.bBox           = true;
3124     header.fileVersion    = tpx_version;
3125     header.fileGeneration = tpx_generation;
3126     header.isDouble       = (sizeof(real) == sizeof(double));
3127
3128     return header;
3129 }
3130
3131 /*! \brief
3132  * Process the body of a TPR file as char buffer.
3133  *
3134  * Reads/writes the information in \p buffer from/to the \p serializer
3135  * provided to the function. Does not interact with the actual
3136  * TPR datastructures but with an in memory representation of the
3137  * data, so that this data can be efficiently read or written from/to
3138  * an original source.
3139  *
3140  * \param[in] serializer The abstract serializer used for reading or writing
3141  *                       the information in \p buffer.
3142  * \param[in,out] buffer Information from TPR file as char buffer.
3143  */
3144 static void doTpxBodyBuffer(gmx::ISerializer* serializer, gmx::ArrayRef<char> buffer)
3145 {
3146     serializer->doCharArray(buffer.data(), buffer.size());
3147 }
3148
3149 /*! \brief
3150  * Populates simulation datastructures.
3151  *
3152  * Here the information from the serialization interface \p serializer
3153  * is used to first populate the datastructures containing the simulation
3154  * information. Depending on the version found in the header \p tpx,
3155  * this is done using the new reading of the data as one block from disk,
3156  * followed by complete deserialization of the information read from there.
3157  * Otherwise, the datastructures are populated as before one by one from disk.
3158  * The second version is the default for the legacy tools that read the
3159  * coordinates and velocities separate from the state.
3160  *
3161  * After reading in the data, a separate buffer is populated from them
3162  * containing only \p ir and \p mtop that can be communicated directly
3163  * to nodes needing the information to set up a simulation.
3164  *
3165  * \param[in] tpx The file header.
3166  * \param[in] serializer The Serialization interface used to read the TPR.
3167  * \param[out] ir Input rec to populate.
3168  * \param[out] state State vectors to populate.
3169  * \param[out] x Coordinates to populate if needed.
3170  * \param[out] v Velocities to populate if needed.
3171  * \param[out] mtop Global topology to populate.
3172  *
3173  * \returns Partial de-serialized TPR used for communication to nodes.
3174  */
3175 static PartialDeserializedTprFile readTpxBody(TpxFileHeader*    tpx,
3176                                               gmx::ISerializer* serializer,
3177                                               t_inputrec*       ir,
3178                                               t_state*          state,
3179                                               rvec*             x,
3180                                               rvec*             v,
3181                                               gmx_mtop_t*       mtop)
3182 {
3183     PartialDeserializedTprFile partialDeserializedTpr;
3184     if (tpx->fileVersion >= tpxv_AddSizeField && tpx->fileGeneration >= 27)
3185     {
3186         partialDeserializedTpr.body.resize(tpx->sizeOfTprBody);
3187         partialDeserializedTpr.header = *tpx;
3188         doTpxBodyBuffer(serializer, partialDeserializedTpr.body);
3189
3190         partialDeserializedTpr.ePBC =
3191                 completeTprDeserialization(&partialDeserializedTpr, ir, state, x, v, mtop);
3192     }
3193     else
3194     {
3195         partialDeserializedTpr.ePBC = do_tpx_body(serializer, tpx, ir, state, x, v, mtop);
3196     }
3197     // Update header to system info for communication to nodes.
3198     // As we only need to communicate the inputrec and mtop to other nodes,
3199     // we prepare a new char buffer with the information we have already read
3200     // in on master.
3201     partialDeserializedTpr.header = populateTpxHeader(*state, ir, mtop);
3202     gmx::InMemorySerializer tprBodySerializer;
3203     do_tpx_body(&tprBodySerializer, &partialDeserializedTpr.header, ir, mtop);
3204     partialDeserializedTpr.body = tprBodySerializer.finishAndGetBuffer();
3205
3206     return partialDeserializedTpr;
3207 }
3208
3209 /************************************************************
3210  *
3211  *  The following routines are the exported ones
3212  *
3213  ************************************************************/
3214
3215 TpxFileHeader readTpxHeader(const char* fileName, bool canReadTopologyOnly)
3216 {
3217     t_fileio* fio;
3218
3219     fio = open_tpx(fileName, "r");
3220     gmx::FileIOXdrSerializer serializer(fio);
3221
3222     TpxFileHeader tpx;
3223     do_tpxheader(&serializer, &tpx, fileName, fio, canReadTopologyOnly);
3224     close_tpx(fio);
3225     return tpx;
3226 }
3227
3228 void write_tpx_state(const char* fn, const t_inputrec* ir, const t_state* state, const gmx_mtop_t* mtop)
3229 {
3230     /* To write a state, we first need to write the state information to a buffer before
3231      * we append the raw bytes to the file. For this, the header information needs to be
3232      * populated before we write the main body because it has some information that is
3233      * otherwise not available.
3234      */
3235
3236     t_fileio* fio;
3237
3238     TpxFileHeader tpx = populateTpxHeader(*state, ir, mtop);
3239
3240     gmx::InMemorySerializer tprBodySerializer;
3241
3242     do_tpx_body(&tprBodySerializer, &tpx, const_cast<t_inputrec*>(ir), const_cast<t_state*>(state),
3243                 nullptr, nullptr, const_cast<gmx_mtop_t*>(mtop));
3244
3245     std::vector<char> tprBody = tprBodySerializer.finishAndGetBuffer();
3246     tpx.sizeOfTprBody         = tprBody.size();
3247
3248     fio = open_tpx(fn, "w");
3249     gmx::FileIOXdrSerializer serializer(fio);
3250     do_tpxheader(&serializer, &tpx, fn, fio, ir == nullptr);
3251     doTpxBodyBuffer(&serializer, tprBody);
3252
3253     close_tpx(fio);
3254 }
3255
3256 int completeTprDeserialization(PartialDeserializedTprFile* partialDeserializedTpr,
3257                                t_inputrec*                 ir,
3258                                t_state*                    state,
3259                                rvec*                       x,
3260                                rvec*                       v,
3261                                gmx_mtop_t*                 mtop)
3262 {
3263     gmx::InMemoryDeserializer tprBodyDeserializer(partialDeserializedTpr->body,
3264                                                   partialDeserializedTpr->header.isDouble);
3265     return do_tpx_body(&tprBodyDeserializer, &partialDeserializedTpr->header, ir, state, x, v, mtop);
3266 }
3267
3268 int completeTprDeserialization(PartialDeserializedTprFile* partialDeserializedTpr,
3269                                t_inputrec*                 ir,
3270                                gmx_mtop_t*                 mtop)
3271 {
3272     return completeTprDeserialization(partialDeserializedTpr, ir, nullptr, nullptr, nullptr, mtop);
3273 }
3274
3275 PartialDeserializedTprFile read_tpx_state(const char* fn, t_inputrec* ir, t_state* state, gmx_mtop_t* mtop)
3276 {
3277     t_fileio* fio;
3278     fio = open_tpx(fn, "r");
3279     gmx::FileIOXdrSerializer   serializer(fio);
3280     PartialDeserializedTprFile partialDeserializedTpr;
3281     do_tpxheader(&serializer, &partialDeserializedTpr.header, fn, fio, ir == nullptr);
3282     partialDeserializedTpr =
3283             readTpxBody(&partialDeserializedTpr.header, &serializer, ir, state, nullptr, nullptr, mtop);
3284     close_tpx(fio);
3285     return partialDeserializedTpr;
3286 }
3287
3288 int read_tpx(const char* fn, t_inputrec* ir, matrix box, int* natoms, rvec* x, rvec* v, gmx_mtop_t* mtop)
3289 {
3290     t_fileio* fio;
3291     t_state   state;
3292
3293     TpxFileHeader tpx;
3294     fio = open_tpx(fn, "r");
3295     gmx::FileIOXdrSerializer serializer(fio);
3296     do_tpxheader(&serializer, &tpx, fn, fio, ir == nullptr);
3297     PartialDeserializedTprFile partialDeserializedTpr =
3298             readTpxBody(&tpx, &serializer, ir, &state, x, v, mtop);
3299     close_tpx(fio);
3300     if (mtop != nullptr && natoms != nullptr)
3301     {
3302         *natoms = mtop->natoms;
3303     }
3304     if (box)
3305     {
3306         copy_mat(state.box, box);
3307     }
3308     return partialDeserializedTpr.ePBC;
3309 }
3310
3311 int read_tpx_top(const char* fn, t_inputrec* ir, matrix box, int* natoms, rvec* x, rvec* v, t_topology* top)
3312 {
3313     gmx_mtop_t mtop;
3314     int        ePBC;
3315
3316     ePBC = read_tpx(fn, ir, box, natoms, x, v, &mtop);
3317
3318     *top = gmx_mtop_t_to_t_topology(&mtop, true);
3319
3320     return ePBC;
3321 }
3322
3323 gmx_bool fn2bTPX(const char* file)
3324 {
3325     return (efTPR == fn2ftp(file));
3326 }
3327
3328 void pr_tpxheader(FILE* fp, int indent, const char* title, const TpxFileHeader* sh)
3329 {
3330     if (available(fp, sh, indent, title))
3331     {
3332         indent = pr_title(fp, indent, title);
3333         pr_indent(fp, indent);
3334         fprintf(fp, "bIr    = %spresent\n", sh->bIr ? "" : "not ");
3335         pr_indent(fp, indent);
3336         fprintf(fp, "bBox   = %spresent\n", sh->bBox ? "" : "not ");
3337         pr_indent(fp, indent);
3338         fprintf(fp, "bTop   = %spresent\n", sh->bTop ? "" : "not ");
3339         pr_indent(fp, indent);
3340         fprintf(fp, "bX     = %spresent\n", sh->bX ? "" : "not ");
3341         pr_indent(fp, indent);
3342         fprintf(fp, "bV     = %spresent\n", sh->bV ? "" : "not ");
3343         pr_indent(fp, indent);
3344         fprintf(fp, "bF     = %spresent\n", sh->bF ? "" : "not ");
3345
3346         pr_indent(fp, indent);
3347         fprintf(fp, "natoms = %d\n", sh->natoms);
3348         pr_indent(fp, indent);
3349         fprintf(fp, "lambda = %e\n", sh->lambda);
3350         pr_indent(fp, indent);
3351         fprintf(fp, "buffer size = %" PRId64 "\n", sh->sizeOfTprBody);
3352     }
3353 }