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