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