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