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