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