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