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