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