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