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