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