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