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