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