Merge branch release-5-0 into release-5-1
[alexxy/gromacs.git] / src / gromacs / fileio / tpxio.c
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, 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 <stdlib.h>
44 #include <string.h>
45
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/fileio/filenm.h"
48 #include "gromacs/fileio/gmxfio.h"
49 #include "gromacs/legacyheaders/copyrite.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gromacs/legacyheaders/txtdump.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/topology/atomprop.h"
55 #include "gromacs/topology/block.h"
56 #include "gromacs/topology/mtop_util.h"
57 #include "gromacs/topology/symtab.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/utility/cstringutil.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
63
64 #define TPX_TAG_RELEASE  "release"
65
66 /*! \brief Tag string for the file format written to run input files
67  * written by this version of the code.
68  *
69  * Change this if you want to change the run input format in a feature
70  * branch. This ensures that there will not be different run input
71  * formats around which cannot be distinguished, while not causing
72  * problems rebasing the feature branch onto upstream changes. When
73  * merging with mainstream GROMACS, set this tag string back to
74  * TPX_TAG_RELEASE, and instead add an element to tpxv and set
75  * tpx_version to that.
76  */
77 static const char *tpx_tag = TPX_TAG_RELEASE;
78
79 /*! \brief Enum of values that describe the contents of a tpr file
80  * whose format matches a version number
81  *
82  * The enum helps the code be more self-documenting and ensure merges
83  * do not silently resolve when two patches make the same bump. When
84  * adding new functionality, add a new element to the end of this
85  * enumeration, change the definition of tpx_version, and write code
86  * below that does the right thing according to the value of
87  * file_version. */
88 enum tpxv {
89     tpxv_ComputationalElectrophysiology = 96,                /**< support for ion/water position swaps (computational electrophysiology) */
90     tpxv_Use64BitRandomSeed,                                 /**< change ld_seed from int to gmx_int64_t */
91     tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
92     tpxv_InteractiveMolecularDynamics,                       /**< interactive molecular dynamics (IMD) */
93     tpxv_RemoveObsoleteParameters1,                          /**< remove optimize_fft, dihre_fc, nstcheckpoint */
94     tpxv_PullCoordTypeGeom,                                  /**< add pull type and geometry per group and flat-bottom */
95     tpxv_PullGeomDirRel,                                     /**< add pull geometry direction-relative */
96     tpxv_IntermolecularBondeds                               /**< permit inter-molecular bonded interactions in the topology */
97 };
98
99 /*! \brief Version number of the file format written to run input
100  * files by this version of the code.
101  *
102  * The tpx_version number should be increased whenever the file format
103  * in the main development branch changes, generally to the highest
104  * value present in tpxv. Backward compatibility for reading old run
105  * input files is maintained by checking this version number against
106  * that of the file and then using the correct code path.
107  *
108  * When developing a feature branch that needs to change the run input
109  * file format, change tpx_tag instead. */
110 static const int tpx_version = tpxv_IntermolecularBondeds;
111
112
113 /* This number should only be increased when you edit the TOPOLOGY section
114  * or the HEADER of the tpx format.
115  * This way we can maintain forward compatibility too for all analysis tools
116  * and/or external programs that only need to know the atom/residue names,
117  * charges, and bond connectivity.
118  *
119  * It first appeared in tpx version 26, when I also moved the inputrecord
120  * to the end of the tpx file, so we can just skip it if we only
121  * want the topology.
122  *
123  * In particular, it must be increased when adding new elements to
124  * ftupd, so that old code can read new .tpr files.
125  */
126 static const int tpx_generation = 26;
127
128 /* This number should be the most recent backwards incompatible version
129  * I.e., if this number is 9, we cannot read tpx version 9 with this code.
130  */
131 static const int tpx_incompatible_version = 9;
132
133
134
135 /* Struct used to maintain tpx compatibility when function types are added */
136 typedef struct {
137     int fvnr;  /* file version number in which the function type first appeared */
138     int ftype; /* function type */
139 } t_ftupd;
140
141 /*
142  * The entries should be ordered in:
143  * 1. ascending file version number
144  * 2. ascending function type number
145  */
146 /*static const t_ftupd ftupd[] = {
147    { 20, F_CUBICBONDS        },
148    { 20, F_CONNBONDS         },
149    { 20, F_HARMONIC          },
150    { 20, F_EQM,              },
151    { 22, F_DISRESVIOL        },
152    { 22, F_ORIRES            },
153    { 22, F_ORIRESDEV         },
154    { 26, F_FOURDIHS          },
155    { 26, F_PIDIHS            },
156    { 26, F_DIHRES            },
157    { 26, F_DIHRESVIOL        },
158    { 30, F_CROSS_BOND_BONDS  },
159    { 30, F_CROSS_BOND_ANGLES },
160    { 30, F_UREY_BRADLEY      },
161    { 30, F_POLARIZATION      },
162    { 54, F_DHDL_CON          },
163    };*/
164 /*
165  * The entries should be ordered in:
166  * 1. ascending function type number
167  * 2. ascending file version number
168  */
169 /* question; what is the purpose of the commented code above? */
170 static const t_ftupd ftupd[] = {
171     { 20, F_CUBICBONDS        },
172     { 20, F_CONNBONDS         },
173     { 20, F_HARMONIC          },
174     { 34, F_FENEBONDS         },
175     { 43, F_TABBONDS          },
176     { 43, F_TABBONDSNC        },
177     { 70, F_RESTRBONDS        },
178     { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
179     { 76, F_LINEAR_ANGLES     },
180     { 30, F_CROSS_BOND_BONDS  },
181     { 30, F_CROSS_BOND_ANGLES },
182     { 30, F_UREY_BRADLEY      },
183     { 34, F_QUARTIC_ANGLES    },
184     { 43, F_TABANGLES         },
185     { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
186     { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
187     { 26, F_FOURDIHS          },
188     { 26, F_PIDIHS            },
189     { 43, F_TABDIHS           },
190     { 65, F_CMAP              },
191     { 60, F_GB12              },
192     { 61, F_GB13              },
193     { 61, F_GB14              },
194     { 72, F_GBPOL             },
195     { 72, F_NPSOLVATION       },
196     { 41, F_LJC14_Q           },
197     { 41, F_LJC_PAIRS_NB      },
198     { 32, F_BHAM_LR           },
199     { 32, F_RF_EXCL           },
200     { 32, F_COUL_RECIP        },
201     { 93, F_LJ_RECIP          },
202     { 46, F_DPD               },
203     { 30, F_POLARIZATION      },
204     { 36, F_THOLE_POL         },
205     { 90, F_FBPOSRES          },
206     { 22, F_DISRESVIOL        },
207     { 22, F_ORIRES            },
208     { 22, F_ORIRESDEV         },
209     { 26, F_DIHRES            },
210     { 26, F_DIHRESVIOL        },
211     { 49, F_VSITE4FDN         },
212     { 50, F_VSITEN            },
213     { 46, F_COM_PULL          },
214     { 20, F_EQM               },
215     { 46, F_ECONSERVED        },
216     { 69, F_VTEMP_NOLONGERUSED},
217     { 66, F_PDISPCORR         },
218     { 54, F_DVDL_CONSTR       },
219     { 76, F_ANHARM_POL        },
220     { 79, F_DVDL_COUL         },
221     { 79, F_DVDL_VDW,         },
222     { 79, F_DVDL_BONDED,      },
223     { 79, F_DVDL_RESTRAINT    },
224     { 79, F_DVDL_TEMPERATURE  },
225 };
226 #define NFTUPD asize(ftupd)
227
228 /* Needed for backward compatibility */
229 #define MAXNODES 256
230
231 /**************************************************************
232  *
233  * Now the higer level routines that do io of the structures and arrays
234  *
235  **************************************************************/
236 static void do_pullgrp_tpx_pre95(t_fileio     *fio,
237                                  t_pull_group *pgrp,
238                                  t_pull_coord *pcrd,
239                                  gmx_bool      bRead,
240                                  int           file_version)
241 {
242     int  i;
243     rvec tmp;
244
245     gmx_fio_do_int(fio, pgrp->nat);
246     if (bRead)
247     {
248         snew(pgrp->ind, pgrp->nat);
249     }
250     gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
251     gmx_fio_do_int(fio, pgrp->nweight);
252     if (bRead)
253     {
254         snew(pgrp->weight, pgrp->nweight);
255     }
256     gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
257     gmx_fio_do_int(fio, pgrp->pbcatom);
258     gmx_fio_do_rvec(fio, pcrd->vec);
259     clear_rvec(pcrd->origin);
260     gmx_fio_do_rvec(fio, tmp);
261     pcrd->init = tmp[0];
262     gmx_fio_do_real(fio, pcrd->rate);
263     gmx_fio_do_real(fio, pcrd->k);
264     if (file_version >= 56)
265     {
266         gmx_fio_do_real(fio, pcrd->kB);
267     }
268     else
269     {
270         pcrd->kB = pcrd->k;
271     }
272 }
273
274 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
275 {
276     int      i;
277
278     gmx_fio_do_int(fio, pgrp->nat);
279     if (bRead)
280     {
281         snew(pgrp->ind, pgrp->nat);
282     }
283     gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
284     gmx_fio_do_int(fio, pgrp->nweight);
285     if (bRead)
286     {
287         snew(pgrp->weight, pgrp->nweight);
288     }
289     gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
290     gmx_fio_do_int(fio, pgrp->pbcatom);
291 }
292
293 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd, int file_version,
294                           int ePullOld, int eGeomOld, ivec dimOld)
295 {
296     int      i;
297
298     gmx_fio_do_int(fio, pcrd->group[0]);
299     gmx_fio_do_int(fio, pcrd->group[1]);
300     if (file_version >= tpxv_PullCoordTypeGeom)
301     {
302         gmx_fio_do_int(fio,  pcrd->eType);
303         gmx_fio_do_int(fio,  pcrd->eGeom);
304         if (pcrd->eGeom == epullgDIRRELATIVE)
305         {
306             gmx_fio_do_int(fio, pcrd->group[2]);
307             gmx_fio_do_int(fio, pcrd->group[3]);
308         }
309         gmx_fio_do_ivec(fio, pcrd->dim);
310     }
311     else
312     {
313         pcrd->eType = ePullOld;
314         pcrd->eGeom = eGeomOld;
315         copy_ivec(dimOld, pcrd->dim);
316     }
317     gmx_fio_do_rvec(fio, pcrd->origin);
318     gmx_fio_do_rvec(fio, pcrd->vec);
319     if (file_version >= tpxv_PullCoordTypeGeom)
320     {
321         gmx_fio_do_gmx_bool(fio, pcrd->bStart);
322     }
323     else
324     {
325         /* This parameter is only printed, but not actually used by mdrun */
326         pcrd->bStart = FALSE;
327     }
328     gmx_fio_do_real(fio, pcrd->init);
329     gmx_fio_do_real(fio, pcrd->rate);
330     gmx_fio_do_real(fio, pcrd->k);
331     gmx_fio_do_real(fio, pcrd->kB);
332 }
333
334 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
335 {
336     /* i is used in the ndo_double macro*/
337     int      i;
338     real     fv;
339     real     rdum;
340     int      n_lambda = fepvals->n_lambda;
341
342     /* reset the lambda calculation window */
343     fepvals->lambda_start_n = 0;
344     fepvals->lambda_stop_n  = n_lambda;
345     if (file_version >= 79)
346     {
347         if (n_lambda > 0)
348         {
349             if (bRead)
350             {
351                 snew(expand->init_lambda_weights, n_lambda);
352             }
353             gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
354             gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
355         }
356
357         gmx_fio_do_int(fio, expand->nstexpanded);
358         gmx_fio_do_int(fio, expand->elmcmove);
359         gmx_fio_do_int(fio, expand->elamstats);
360         gmx_fio_do_int(fio, expand->lmc_repeats);
361         gmx_fio_do_int(fio, expand->gibbsdeltalam);
362         gmx_fio_do_int(fio, expand->lmc_forced_nstart);
363         gmx_fio_do_int(fio, expand->lmc_seed);
364         gmx_fio_do_real(fio, expand->mc_temp);
365         gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
366         gmx_fio_do_int(fio, expand->nstTij);
367         gmx_fio_do_int(fio, expand->minvarmin);
368         gmx_fio_do_int(fio, expand->c_range);
369         gmx_fio_do_real(fio, expand->wl_scale);
370         gmx_fio_do_real(fio, expand->wl_ratio);
371         gmx_fio_do_real(fio, expand->init_wl_delta);
372         gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
373         gmx_fio_do_int(fio, expand->elmceq);
374         gmx_fio_do_int(fio, expand->equil_steps);
375         gmx_fio_do_int(fio, expand->equil_samples);
376         gmx_fio_do_int(fio, expand->equil_n_at_lam);
377         gmx_fio_do_real(fio, expand->equil_wl_delta);
378         gmx_fio_do_real(fio, expand->equil_ratio);
379     }
380 }
381
382 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
383                            int file_version)
384 {
385     if (file_version >= 79)
386     {
387         gmx_fio_do_int(fio, simtemp->eSimTempScale);
388         gmx_fio_do_real(fio, simtemp->simtemp_high);
389         gmx_fio_do_real(fio, simtemp->simtemp_low);
390         if (n_lambda > 0)
391         {
392             if (bRead)
393             {
394                 snew(simtemp->temperatures, n_lambda);
395             }
396             gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
397         }
398     }
399 }
400
401 static void do_imd(t_fileio *fio, t_IMD *imd, gmx_bool bRead)
402 {
403     gmx_fio_do_int(fio, imd->nat);
404     if (bRead)
405     {
406         snew(imd->ind, imd->nat);
407     }
408     gmx_fio_ndo_int(fio, imd->ind, imd->nat);
409 }
410
411 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
412 {
413     /* i is defined in the ndo_double macro; use g to iterate. */
414     int      i, g;
415     real     fv;
416     real     rdum;
417
418     /* free energy values */
419
420     if (file_version >= 79)
421     {
422         gmx_fio_do_int(fio, fepvals->init_fep_state);
423         gmx_fio_do_double(fio, fepvals->init_lambda);
424         gmx_fio_do_double(fio, fepvals->delta_lambda);
425     }
426     else if (file_version >= 59)
427     {
428         gmx_fio_do_double(fio, fepvals->init_lambda);
429         gmx_fio_do_double(fio, fepvals->delta_lambda);
430     }
431     else
432     {
433         gmx_fio_do_real(fio, rdum);
434         fepvals->init_lambda = rdum;
435         gmx_fio_do_real(fio, rdum);
436         fepvals->delta_lambda = rdum;
437     }
438     if (file_version >= 79)
439     {
440         gmx_fio_do_int(fio, fepvals->n_lambda);
441         if (bRead)
442         {
443             snew(fepvals->all_lambda, efptNR);
444         }
445         for (g = 0; g < efptNR; g++)
446         {
447             if (fepvals->n_lambda > 0)
448             {
449                 if (bRead)
450                 {
451                     snew(fepvals->all_lambda[g], fepvals->n_lambda);
452                 }
453                 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
454                 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
455             }
456             else if (fepvals->init_lambda >= 0)
457             {
458                 fepvals->separate_dvdl[efptFEP] = TRUE;
459             }
460         }
461     }
462     else if (file_version >= 64)
463     {
464         gmx_fio_do_int(fio, fepvals->n_lambda);
465         if (bRead)
466         {
467             int g;
468
469             snew(fepvals->all_lambda, efptNR);
470             /* still allocate the all_lambda array's contents. */
471             for (g = 0; g < efptNR; g++)
472             {
473                 if (fepvals->n_lambda > 0)
474                 {
475                     snew(fepvals->all_lambda[g], fepvals->n_lambda);
476                 }
477             }
478         }
479         gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
480                            fepvals->n_lambda);
481         if (fepvals->init_lambda >= 0)
482         {
483             int g, h;
484
485             fepvals->separate_dvdl[efptFEP] = TRUE;
486
487             if (bRead)
488             {
489                 /* copy the contents of the efptFEP lambda component to all
490                    the other components */
491                 for (g = 0; g < efptNR; g++)
492                 {
493                     for (h = 0; h < fepvals->n_lambda; h++)
494                     {
495                         if (g != efptFEP)
496                         {
497                             fepvals->all_lambda[g][h] =
498                                 fepvals->all_lambda[efptFEP][h];
499                         }
500                     }
501                 }
502             }
503         }
504     }
505     else
506     {
507         fepvals->n_lambda     = 0;
508         fepvals->all_lambda   = NULL;
509         if (fepvals->init_lambda >= 0)
510         {
511             fepvals->separate_dvdl[efptFEP] = TRUE;
512         }
513     }
514     if (file_version >= 13)
515     {
516         gmx_fio_do_real(fio, fepvals->sc_alpha);
517     }
518     else
519     {
520         fepvals->sc_alpha = 0;
521     }
522     if (file_version >= 38)
523     {
524         gmx_fio_do_int(fio, fepvals->sc_power);
525     }
526     else
527     {
528         fepvals->sc_power = 2;
529     }
530     if (file_version >= 79)
531     {
532         gmx_fio_do_real(fio, fepvals->sc_r_power);
533     }
534     else
535     {
536         fepvals->sc_r_power = 6.0;
537     }
538     if (file_version >= 15)
539     {
540         gmx_fio_do_real(fio, fepvals->sc_sigma);
541     }
542     else
543     {
544         fepvals->sc_sigma = 0.3;
545     }
546     if (bRead)
547     {
548         if (file_version >= 71)
549         {
550             fepvals->sc_sigma_min = fepvals->sc_sigma;
551         }
552         else
553         {
554             fepvals->sc_sigma_min = 0;
555         }
556     }
557     if (file_version >= 79)
558     {
559         gmx_fio_do_int(fio, fepvals->bScCoul);
560     }
561     else
562     {
563         fepvals->bScCoul = TRUE;
564     }
565     if (file_version >= 64)
566     {
567         gmx_fio_do_int(fio, fepvals->nstdhdl);
568     }
569     else
570     {
571         fepvals->nstdhdl = 1;
572     }
573
574     if (file_version >= 73)
575     {
576         gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
577         gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
578     }
579     else
580     {
581         fepvals->separate_dhdl_file = esepdhdlfileYES;
582         fepvals->dhdl_derivatives   = edhdlderivativesYES;
583     }
584     if (file_version >= 71)
585     {
586         gmx_fio_do_int(fio, fepvals->dh_hist_size);
587         gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
588     }
589     else
590     {
591         fepvals->dh_hist_size    = 0;
592         fepvals->dh_hist_spacing = 0.1;
593     }
594     if (file_version >= 79)
595     {
596         gmx_fio_do_int(fio, fepvals->edHdLPrintEnergy);
597     }
598     else
599     {
600         fepvals->edHdLPrintEnergy = edHdLPrintEnergyNO;
601     }
602
603     /* handle lambda_neighbors */
604     if ((file_version >= 83 && file_version < 90) || file_version >= 92)
605     {
606         gmx_fio_do_int(fio, fepvals->lambda_neighbors);
607         if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
608              (fepvals->init_lambda < 0) )
609         {
610             fepvals->lambda_start_n = (fepvals->init_fep_state -
611                                        fepvals->lambda_neighbors);
612             fepvals->lambda_stop_n = (fepvals->init_fep_state +
613                                       fepvals->lambda_neighbors + 1);
614             if (fepvals->lambda_start_n < 0)
615             {
616                 fepvals->lambda_start_n = 0;;
617             }
618             if (fepvals->lambda_stop_n >= fepvals->n_lambda)
619             {
620                 fepvals->lambda_stop_n = fepvals->n_lambda;
621             }
622         }
623         else
624         {
625             fepvals->lambda_start_n = 0;
626             fepvals->lambda_stop_n  = fepvals->n_lambda;
627         }
628     }
629     else
630     {
631         fepvals->lambda_start_n = 0;
632         fepvals->lambda_stop_n  = fepvals->n_lambda;
633     }
634 }
635
636 static void do_pull(t_fileio *fio, pull_params_t *pull, gmx_bool bRead,
637                     int file_version, int ePullOld)
638 {
639     int  eGeomOld;
640     ivec dimOld;
641     int  g;
642
643     if (file_version >= 95)
644     {
645         gmx_fio_do_int(fio, pull->ngroup);
646     }
647     gmx_fio_do_int(fio, pull->ncoord);
648     if (file_version < 95)
649     {
650         pull->ngroup = pull->ncoord + 1;
651     }
652     if (file_version < tpxv_PullCoordTypeGeom)
653     {
654         real dum;
655
656         gmx_fio_do_int(fio, eGeomOld);
657         gmx_fio_do_ivec(fio, dimOld);
658         /* The inner cylinder radius, now removed */
659         gmx_fio_do_real(fio, dum);
660     }
661     gmx_fio_do_real(fio, pull->cylinder_r);
662     gmx_fio_do_real(fio, pull->constr_tol);
663     if (file_version >= 95)
664     {
665         gmx_fio_do_int(fio, pull->bPrintCOM1);
666         /* With file_version < 95 this value is set below */
667     }
668     if (file_version >= tpxv_PullCoordTypeGeom)
669     {
670         gmx_fio_do_int(fio, pull->bPrintCOM2);
671         gmx_fio_do_int(fio, pull->bPrintRefValue);
672         gmx_fio_do_int(fio, pull->bPrintComp);
673     }
674     else
675     {
676         pull->bPrintCOM2     = FALSE;
677         pull->bPrintRefValue = FALSE;
678         pull->bPrintComp     = TRUE;
679     }
680     gmx_fio_do_int(fio, pull->nstxout);
681     gmx_fio_do_int(fio, pull->nstfout);
682     if (bRead)
683     {
684         snew(pull->group, pull->ngroup);
685         snew(pull->coord, pull->ncoord);
686     }
687     if (file_version < 95)
688     {
689         /* epullgPOS for position pulling, before epullgDIRPBC was removed */
690         if (eGeomOld == epullgDIRPBC)
691         {
692             gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
693         }
694         if (eGeomOld > epullgDIRPBC)
695         {
696             eGeomOld -= 1;
697         }
698
699         for (g = 0; g < pull->ngroup; g++)
700         {
701             /* We read and ignore a pull coordinate for group 0 */
702             do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[max(g-1, 0)],
703                                  bRead, file_version);
704             if (g > 0)
705             {
706                 pull->coord[g-1].group[0] = 0;
707                 pull->coord[g-1].group[1] = g;
708             }
709         }
710
711         pull->bPrintCOM1 = (pull->group[0].nat > 0);
712     }
713     else
714     {
715         for (g = 0; g < pull->ngroup; g++)
716         {
717             do_pull_group(fio, &pull->group[g], bRead);
718         }
719         for (g = 0; g < pull->ncoord; g++)
720         {
721             do_pull_coord(fio, &pull->coord[g],
722                           file_version, ePullOld, eGeomOld, dimOld);
723         }
724     }
725 }
726
727
728 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
729 {
730     int      i;
731
732     gmx_fio_do_int(fio, rotg->eType);
733     gmx_fio_do_int(fio, rotg->bMassW);
734     gmx_fio_do_int(fio, rotg->nat);
735     if (bRead)
736     {
737         snew(rotg->ind, rotg->nat);
738     }
739     gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
740     if (bRead)
741     {
742         snew(rotg->x_ref, rotg->nat);
743     }
744     gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
745     gmx_fio_do_rvec(fio, rotg->vec);
746     gmx_fio_do_rvec(fio, rotg->pivot);
747     gmx_fio_do_real(fio, rotg->rate);
748     gmx_fio_do_real(fio, rotg->k);
749     gmx_fio_do_real(fio, rotg->slab_dist);
750     gmx_fio_do_real(fio, rotg->min_gaussian);
751     gmx_fio_do_real(fio, rotg->eps);
752     gmx_fio_do_int(fio, rotg->eFittype);
753     gmx_fio_do_int(fio, rotg->PotAngle_nstep);
754     gmx_fio_do_real(fio, rotg->PotAngle_step);
755 }
756
757 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
758 {
759     int g;
760
761     gmx_fio_do_int(fio, rot->ngrp);
762     gmx_fio_do_int(fio, rot->nstrout);
763     gmx_fio_do_int(fio, rot->nstsout);
764     if (bRead)
765     {
766         snew(rot->grp, rot->ngrp);
767     }
768     for (g = 0; g < rot->ngrp; g++)
769     {
770         do_rotgrp(fio, &rot->grp[g], bRead);
771     }
772 }
773
774
775 static void do_swapcoords(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead)
776 {
777     int i, j;
778
779
780     gmx_fio_do_int(fio, swap->nat);
781     gmx_fio_do_int(fio, swap->nat_sol);
782     for (j = 0; j < 2; j++)
783     {
784         gmx_fio_do_int(fio, swap->nat_split[j]);
785         gmx_fio_do_int(fio, swap->massw_split[j]);
786     }
787     gmx_fio_do_int(fio, swap->nstswap);
788     gmx_fio_do_int(fio, swap->nAverage);
789     gmx_fio_do_real(fio, swap->threshold);
790     gmx_fio_do_real(fio, swap->cyl0r);
791     gmx_fio_do_real(fio, swap->cyl0u);
792     gmx_fio_do_real(fio, swap->cyl0l);
793     gmx_fio_do_real(fio, swap->cyl1r);
794     gmx_fio_do_real(fio, swap->cyl1u);
795     gmx_fio_do_real(fio, swap->cyl1l);
796
797     if (bRead)
798     {
799         snew(swap->ind, swap->nat);
800         snew(swap->ind_sol, swap->nat_sol);
801         for (j = 0; j < 2; j++)
802         {
803             snew(swap->ind_split[j], swap->nat_split[j]);
804         }
805     }
806
807     gmx_fio_ndo_int(fio, swap->ind, swap->nat);
808     gmx_fio_ndo_int(fio, swap->ind_sol, swap->nat_sol);
809     for (j = 0; j < 2; j++)
810     {
811         gmx_fio_ndo_int(fio, swap->ind_split[j], swap->nat_split[j]);
812     }
813
814     for (j = 0; j < eCompNR; j++)
815     {
816         gmx_fio_do_int(fio, swap->nanions[j]);
817         gmx_fio_do_int(fio, swap->ncations[j]);
818     }
819
820 }
821
822
823 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
824                         int file_version, real *fudgeQQ)
825 {
826     int      i, j, k, *tmp, idum = 0;
827     real     rdum, bd_temp;
828     rvec     vdum;
829     gmx_bool bSimAnn, bdum = 0;
830     real     zerotemptime, finish_t, init_temp, finish_temp;
831
832     if (file_version != tpx_version)
833     {
834         /* Give a warning about features that are not accessible */
835         fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
836                 file_version, tpx_version);
837     }
838
839     if (bRead)
840     {
841         init_inputrec(ir);
842     }
843
844     if (file_version == 0)
845     {
846         return;
847     }
848
849     /* Basic inputrec stuff */
850     gmx_fio_do_int(fio, ir->eI);
851     if (file_version >= 62)
852     {
853         gmx_fio_do_int64(fio, ir->nsteps);
854     }
855     else
856     {
857         gmx_fio_do_int(fio, idum);
858         ir->nsteps = idum;
859     }
860
861     if (file_version > 25)
862     {
863         if (file_version >= 62)
864         {
865             gmx_fio_do_int64(fio, ir->init_step);
866         }
867         else
868         {
869             gmx_fio_do_int(fio, idum);
870             ir->init_step = idum;
871         }
872     }
873     else
874     {
875         ir->init_step = 0;
876     }
877
878     if (file_version >= 58)
879     {
880         gmx_fio_do_int(fio, ir->simulation_part);
881     }
882     else
883     {
884         ir->simulation_part = 1;
885     }
886
887     if (file_version >= 67)
888     {
889         gmx_fio_do_int(fio, ir->nstcalcenergy);
890     }
891     else
892     {
893         ir->nstcalcenergy = 1;
894     }
895     if (file_version < 53)
896     {
897         /* The pbc info has been moved out of do_inputrec,
898          * since we always want it, also without reading the inputrec.
899          */
900         gmx_fio_do_int(fio, ir->ePBC);
901         if ((file_version <= 15) && (ir->ePBC == 2))
902         {
903             ir->ePBC = epbcNONE;
904         }
905         if (file_version >= 45)
906         {
907             gmx_fio_do_int(fio, ir->bPeriodicMols);
908         }
909         else
910         {
911             if (ir->ePBC == 2)
912             {
913                 ir->ePBC          = epbcXYZ;
914                 ir->bPeriodicMols = TRUE;
915             }
916             else
917             {
918                 ir->bPeriodicMols = FALSE;
919             }
920         }
921     }
922     if (file_version >= 81)
923     {
924         gmx_fio_do_int(fio, ir->cutoff_scheme);
925         if (file_version < 94)
926         {
927             ir->cutoff_scheme = 1 - ir->cutoff_scheme;
928         }
929     }
930     else
931     {
932         ir->cutoff_scheme = ecutsGROUP;
933     }
934     gmx_fio_do_int(fio, ir->ns_type);
935     gmx_fio_do_int(fio, ir->nstlist);
936     gmx_fio_do_int(fio, idum); /* used to be ndelta; not used anymore */
937     if (file_version < 41)
938     {
939         gmx_fio_do_int(fio, idum);
940         gmx_fio_do_int(fio, idum);
941     }
942     if (file_version >= 45)
943     {
944         gmx_fio_do_real(fio, ir->rtpi);
945     }
946     else
947     {
948         ir->rtpi = 0.05;
949     }
950     gmx_fio_do_int(fio, ir->nstcomm);
951     if (file_version > 34)
952     {
953         gmx_fio_do_int(fio, ir->comm_mode);
954     }
955     else if (ir->nstcomm < 0)
956     {
957         ir->comm_mode = ecmANGULAR;
958     }
959     else
960     {
961         ir->comm_mode = ecmLINEAR;
962     }
963     ir->nstcomm = abs(ir->nstcomm);
964
965     /* ignore nstcheckpoint */
966     if (file_version > 25 && file_version < tpxv_RemoveObsoleteParameters1)
967     {
968         gmx_fio_do_int(fio, idum);
969     }
970
971     gmx_fio_do_int(fio, ir->nstcgsteep);
972
973     if (file_version >= 30)
974     {
975         gmx_fio_do_int(fio, ir->nbfgscorr);
976     }
977     else if (bRead)
978     {
979         ir->nbfgscorr = 10;
980     }
981
982     gmx_fio_do_int(fio, ir->nstlog);
983     gmx_fio_do_int(fio, ir->nstxout);
984     gmx_fio_do_int(fio, ir->nstvout);
985     gmx_fio_do_int(fio, ir->nstfout);
986     gmx_fio_do_int(fio, ir->nstenergy);
987     gmx_fio_do_int(fio, ir->nstxout_compressed);
988     if (file_version >= 59)
989     {
990         gmx_fio_do_double(fio, ir->init_t);
991         gmx_fio_do_double(fio, ir->delta_t);
992     }
993     else
994     {
995         gmx_fio_do_real(fio, rdum);
996         ir->init_t = rdum;
997         gmx_fio_do_real(fio, rdum);
998         ir->delta_t = rdum;
999     }
1000     gmx_fio_do_real(fio, ir->x_compression_precision);
1001     if (file_version < 19)
1002     {
1003         gmx_fio_do_int(fio, idum);
1004         gmx_fio_do_int(fio, idum);
1005     }
1006     if (file_version < 18)
1007     {
1008         gmx_fio_do_int(fio, idum);
1009     }
1010     if (file_version >= 81)
1011     {
1012         gmx_fio_do_real(fio, ir->verletbuf_tol);
1013     }
1014     else
1015     {
1016         ir->verletbuf_tol = 0;
1017     }
1018     gmx_fio_do_real(fio, ir->rlist);
1019     if (file_version >= 67)
1020     {
1021         gmx_fio_do_real(fio, ir->rlistlong);
1022     }
1023     if (file_version >= 82 && file_version != 90)
1024     {
1025         gmx_fio_do_int(fio, ir->nstcalclr);
1026     }
1027     else
1028     {
1029         /* Calculate at NS steps */
1030         ir->nstcalclr = ir->nstlist;
1031     }
1032     gmx_fio_do_int(fio, ir->coulombtype);
1033     if (file_version < 32 && ir->coulombtype == eelRF)
1034     {
1035         ir->coulombtype = eelRF_NEC;
1036     }
1037     if (file_version >= 81)
1038     {
1039         gmx_fio_do_int(fio, ir->coulomb_modifier);
1040     }
1041     else
1042     {
1043         ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1044     }
1045     gmx_fio_do_real(fio, ir->rcoulomb_switch);
1046     gmx_fio_do_real(fio, ir->rcoulomb);
1047     gmx_fio_do_int(fio, ir->vdwtype);
1048     if (file_version >= 81)
1049     {
1050         gmx_fio_do_int(fio, ir->vdw_modifier);
1051     }
1052     else
1053     {
1054         ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1055     }
1056     gmx_fio_do_real(fio, ir->rvdw_switch);
1057     gmx_fio_do_real(fio, ir->rvdw);
1058     if (file_version < 67)
1059     {
1060         ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1061     }
1062     gmx_fio_do_int(fio, ir->eDispCorr);
1063     gmx_fio_do_real(fio, ir->epsilon_r);
1064     if (file_version >= 37)
1065     {
1066         gmx_fio_do_real(fio, ir->epsilon_rf);
1067     }
1068     else
1069     {
1070         if (EEL_RF(ir->coulombtype))
1071         {
1072             ir->epsilon_rf = ir->epsilon_r;
1073             ir->epsilon_r  = 1.0;
1074         }
1075         else
1076         {
1077             ir->epsilon_rf = 1.0;
1078         }
1079     }
1080     if (file_version >= 29)
1081     {
1082         gmx_fio_do_real(fio, ir->tabext);
1083     }
1084     else
1085     {
1086         ir->tabext = 1.0;
1087     }
1088
1089     if (file_version > 25)
1090     {
1091         gmx_fio_do_int(fio, ir->gb_algorithm);
1092         gmx_fio_do_int(fio, ir->nstgbradii);
1093         gmx_fio_do_real(fio, ir->rgbradii);
1094         gmx_fio_do_real(fio, ir->gb_saltconc);
1095         gmx_fio_do_int(fio, ir->implicit_solvent);
1096     }
1097     else
1098     {
1099         ir->gb_algorithm     = egbSTILL;
1100         ir->nstgbradii       = 1;
1101         ir->rgbradii         = 1.0;
1102         ir->gb_saltconc      = 0;
1103         ir->implicit_solvent = eisNO;
1104     }
1105     if (file_version >= 55)
1106     {
1107         gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1108         gmx_fio_do_real(fio, ir->gb_obc_alpha);
1109         gmx_fio_do_real(fio, ir->gb_obc_beta);
1110         gmx_fio_do_real(fio, ir->gb_obc_gamma);
1111         if (file_version >= 60)
1112         {
1113             gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1114             gmx_fio_do_int(fio, ir->sa_algorithm);
1115         }
1116         else
1117         {
1118             ir->gb_dielectric_offset = 0.009;
1119             ir->sa_algorithm         = esaAPPROX;
1120         }
1121         gmx_fio_do_real(fio, ir->sa_surface_tension);
1122
1123         /* Override sa_surface_tension if it is not changed in the mpd-file */
1124         if (ir->sa_surface_tension < 0)
1125         {
1126             if (ir->gb_algorithm == egbSTILL)
1127             {
1128                 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1129             }
1130             else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1131             {
1132                 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1133             }
1134         }
1135
1136     }
1137     else
1138     {
1139         /* Better use sensible values than insane (0.0) ones... */
1140         ir->gb_epsilon_solvent = 80;
1141         ir->gb_obc_alpha       = 1.0;
1142         ir->gb_obc_beta        = 0.8;
1143         ir->gb_obc_gamma       = 4.85;
1144         ir->sa_surface_tension = 2.092;
1145     }
1146
1147
1148     if (file_version >= 81)
1149     {
1150         gmx_fio_do_real(fio, ir->fourier_spacing);
1151     }
1152     else
1153     {
1154         ir->fourier_spacing = 0.0;
1155     }
1156     gmx_fio_do_int(fio, ir->nkx);
1157     gmx_fio_do_int(fio, ir->nky);
1158     gmx_fio_do_int(fio, ir->nkz);
1159     gmx_fio_do_int(fio, ir->pme_order);
1160     gmx_fio_do_real(fio, ir->ewald_rtol);
1161
1162     if (file_version >= 93)
1163     {
1164         gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1165     }
1166     else
1167     {
1168         ir->ewald_rtol_lj = ir->ewald_rtol;
1169     }
1170
1171     if (file_version >= 24)
1172     {
1173         gmx_fio_do_int(fio, ir->ewald_geometry);
1174     }
1175     else
1176     {
1177         ir->ewald_geometry = eewg3D;
1178     }
1179
1180     if (file_version <= 17)
1181     {
1182         ir->epsilon_surface = 0;
1183         if (file_version == 17)
1184         {
1185             gmx_fio_do_int(fio, idum);
1186         }
1187     }
1188     else
1189     {
1190         gmx_fio_do_real(fio, ir->epsilon_surface);
1191     }
1192
1193     /* ignore bOptFFT */
1194     if (file_version < tpxv_RemoveObsoleteParameters1)
1195     {
1196         gmx_fio_do_gmx_bool(fio, bdum);
1197     }
1198
1199     if (file_version >= 93)
1200     {
1201         gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1202     }
1203     gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1204     gmx_fio_do_int(fio, ir->etc);
1205     /* before version 18, ir->etc was a gmx_bool (ir->btc),
1206      * but the values 0 and 1 still mean no and
1207      * berendsen temperature coupling, respectively.
1208      */
1209     if (file_version >= 79)
1210     {
1211         gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1212     }
1213     if (file_version >= 71)
1214     {
1215         gmx_fio_do_int(fio, ir->nsttcouple);
1216     }
1217     else
1218     {
1219         ir->nsttcouple = ir->nstcalcenergy;
1220     }
1221     if (file_version <= 15)
1222     {
1223         gmx_fio_do_int(fio, idum);
1224     }
1225     if (file_version <= 17)
1226     {
1227         gmx_fio_do_int(fio, ir->epct);
1228         if (file_version <= 15)
1229         {
1230             if (ir->epct == 5)
1231             {
1232                 ir->epct = epctSURFACETENSION;
1233             }
1234             gmx_fio_do_int(fio, idum);
1235         }
1236         ir->epct -= 1;
1237         /* we have removed the NO alternative at the beginning */
1238         if (ir->epct == -1)
1239         {
1240             ir->epc  = epcNO;
1241             ir->epct = epctISOTROPIC;
1242         }
1243         else
1244         {
1245             ir->epc = epcBERENDSEN;
1246         }
1247     }
1248     else
1249     {
1250         gmx_fio_do_int(fio, ir->epc);
1251         gmx_fio_do_int(fio, ir->epct);
1252     }
1253     if (file_version >= 71)
1254     {
1255         gmx_fio_do_int(fio, ir->nstpcouple);
1256     }
1257     else
1258     {
1259         ir->nstpcouple = ir->nstcalcenergy;
1260     }
1261     gmx_fio_do_real(fio, ir->tau_p);
1262     if (file_version <= 15)
1263     {
1264         gmx_fio_do_rvec(fio, vdum);
1265         clear_mat(ir->ref_p);
1266         for (i = 0; i < DIM; i++)
1267         {
1268             ir->ref_p[i][i] = vdum[i];
1269         }
1270     }
1271     else
1272     {
1273         gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1274         gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1275         gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1276     }
1277     if (file_version <= 15)
1278     {
1279         gmx_fio_do_rvec(fio, vdum);
1280         clear_mat(ir->compress);
1281         for (i = 0; i < DIM; i++)
1282         {
1283             ir->compress[i][i] = vdum[i];
1284         }
1285     }
1286     else
1287     {
1288         gmx_fio_do_rvec(fio, ir->compress[XX]);
1289         gmx_fio_do_rvec(fio, ir->compress[YY]);
1290         gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1291     }
1292     if (file_version >= 47)
1293     {
1294         gmx_fio_do_int(fio, ir->refcoord_scaling);
1295         gmx_fio_do_rvec(fio, ir->posres_com);
1296         gmx_fio_do_rvec(fio, ir->posres_comB);
1297     }
1298     else
1299     {
1300         ir->refcoord_scaling = erscNO;
1301         clear_rvec(ir->posres_com);
1302         clear_rvec(ir->posres_comB);
1303     }
1304     if ((file_version > 25) && (file_version < 79))
1305     {
1306         gmx_fio_do_int(fio, ir->andersen_seed);
1307     }
1308     else
1309     {
1310         ir->andersen_seed = 0;
1311     }
1312     if (file_version < 26)
1313     {
1314         gmx_fio_do_gmx_bool(fio, bSimAnn);
1315         gmx_fio_do_real(fio, zerotemptime);
1316     }
1317
1318     if (file_version < 37)
1319     {
1320         gmx_fio_do_real(fio, rdum);
1321     }
1322
1323     gmx_fio_do_real(fio, ir->shake_tol);
1324     if (file_version < 54)
1325     {
1326         gmx_fio_do_real(fio, *fudgeQQ);
1327     }
1328
1329     gmx_fio_do_int(fio, ir->efep);
1330     if (file_version <= 14 && ir->efep != efepNO)
1331     {
1332         ir->efep = efepYES;
1333     }
1334     do_fepvals(fio, ir->fepvals, bRead, file_version);
1335
1336     if (file_version >= 79)
1337     {
1338         gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1339         if (ir->bSimTemp)
1340         {
1341             ir->bSimTemp = TRUE;
1342         }
1343     }
1344     else
1345     {
1346         ir->bSimTemp = FALSE;
1347     }
1348     if (ir->bSimTemp)
1349     {
1350         do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1351     }
1352
1353     if (file_version >= 79)
1354     {
1355         gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1356         if (ir->bExpanded)
1357         {
1358             ir->bExpanded = TRUE;
1359         }
1360         else
1361         {
1362             ir->bExpanded = FALSE;
1363         }
1364     }
1365     if (ir->bExpanded)
1366     {
1367         do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1368     }
1369     if (file_version >= 57)
1370     {
1371         gmx_fio_do_int(fio, ir->eDisre);
1372     }
1373     gmx_fio_do_int(fio, ir->eDisreWeighting);
1374     if (file_version < 22)
1375     {
1376         if (ir->eDisreWeighting == 0)
1377         {
1378             ir->eDisreWeighting = edrwEqual;
1379         }
1380         else
1381         {
1382             ir->eDisreWeighting = edrwConservative;
1383         }
1384     }
1385     gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1386     gmx_fio_do_real(fio, ir->dr_fc);
1387     gmx_fio_do_real(fio, ir->dr_tau);
1388     gmx_fio_do_int(fio, ir->nstdisreout);
1389     if (file_version >= 22)
1390     {
1391         gmx_fio_do_real(fio, ir->orires_fc);
1392         gmx_fio_do_real(fio, ir->orires_tau);
1393         gmx_fio_do_int(fio, ir->nstorireout);
1394     }
1395     else
1396     {
1397         ir->orires_fc   = 0;
1398         ir->orires_tau  = 0;
1399         ir->nstorireout = 0;
1400     }
1401
1402     /* ignore dihre_fc */
1403     if (file_version >= 26 && file_version < 79)
1404     {
1405         gmx_fio_do_real(fio, rdum);
1406         if (file_version < 56)
1407         {
1408             gmx_fio_do_real(fio, rdum);
1409             gmx_fio_do_int(fio, idum);
1410         }
1411     }
1412
1413     gmx_fio_do_real(fio, ir->em_stepsize);
1414     gmx_fio_do_real(fio, ir->em_tol);
1415     if (file_version >= 22)
1416     {
1417         gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1418     }
1419     else if (bRead)
1420     {
1421         ir->bShakeSOR = TRUE;
1422     }
1423     if (file_version >= 11)
1424     {
1425         gmx_fio_do_int(fio, ir->niter);
1426     }
1427     else if (bRead)
1428     {
1429         ir->niter = 25;
1430         fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1431                 ir->niter);
1432     }
1433     if (file_version >= 21)
1434     {
1435         gmx_fio_do_real(fio, ir->fc_stepsize);
1436     }
1437     else
1438     {
1439         ir->fc_stepsize = 0;
1440     }
1441     gmx_fio_do_int(fio, ir->eConstrAlg);
1442     gmx_fio_do_int(fio, ir->nProjOrder);
1443     gmx_fio_do_real(fio, ir->LincsWarnAngle);
1444     if (file_version <= 14)
1445     {
1446         gmx_fio_do_int(fio, idum);
1447     }
1448     if (file_version >= 26)
1449     {
1450         gmx_fio_do_int(fio, ir->nLincsIter);
1451     }
1452     else if (bRead)
1453     {
1454         ir->nLincsIter = 1;
1455         fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1456                 ir->nLincsIter);
1457     }
1458     if (file_version < 33)
1459     {
1460         gmx_fio_do_real(fio, bd_temp);
1461     }
1462     gmx_fio_do_real(fio, ir->bd_fric);
1463     if (file_version >= tpxv_Use64BitRandomSeed)
1464     {
1465         gmx_fio_do_int64(fio, ir->ld_seed);
1466     }
1467     else
1468     {
1469         gmx_fio_do_int(fio, idum);
1470         ir->ld_seed = idum;
1471     }
1472     if (file_version >= 33)
1473     {
1474         for (i = 0; i < DIM; i++)
1475         {
1476             gmx_fio_do_rvec(fio, ir->deform[i]);
1477         }
1478     }
1479     else
1480     {
1481         for (i = 0; i < DIM; i++)
1482         {
1483             clear_rvec(ir->deform[i]);
1484         }
1485     }
1486     if (file_version >= 14)
1487     {
1488         gmx_fio_do_real(fio, ir->cos_accel);
1489     }
1490     else if (bRead)
1491     {
1492         ir->cos_accel = 0;
1493     }
1494     gmx_fio_do_int(fio, ir->userint1);
1495     gmx_fio_do_int(fio, ir->userint2);
1496     gmx_fio_do_int(fio, ir->userint3);
1497     gmx_fio_do_int(fio, ir->userint4);
1498     gmx_fio_do_real(fio, ir->userreal1);
1499     gmx_fio_do_real(fio, ir->userreal2);
1500     gmx_fio_do_real(fio, ir->userreal3);
1501     gmx_fio_do_real(fio, ir->userreal4);
1502
1503     /* AdResS stuff */
1504     if (file_version >= 77)
1505     {
1506         gmx_fio_do_gmx_bool(fio, ir->bAdress);
1507         if (ir->bAdress)
1508         {
1509             if (bRead)
1510             {
1511                 snew(ir->adress, 1);
1512             }
1513             gmx_fio_do_int(fio, ir->adress->type);
1514             gmx_fio_do_real(fio, ir->adress->const_wf);
1515             gmx_fio_do_real(fio, ir->adress->ex_width);
1516             gmx_fio_do_real(fio, ir->adress->hy_width);
1517             gmx_fio_do_int(fio, ir->adress->icor);
1518             gmx_fio_do_int(fio, ir->adress->site);
1519             gmx_fio_do_rvec(fio, ir->adress->refs);
1520             gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1521             gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1522             gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1523             gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1524
1525             if (bRead)
1526             {
1527                 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1528             }
1529             if (ir->adress->n_tf_grps > 0)
1530             {
1531                 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1532             }
1533             if (bRead)
1534             {
1535                 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1536             }
1537             if (ir->adress->n_energy_grps > 0)
1538             {
1539                 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1540             }
1541         }
1542     }
1543     else
1544     {
1545         ir->bAdress = FALSE;
1546     }
1547
1548     /* pull stuff */
1549     if (file_version >= 48)
1550     {
1551         int ePullOld = 0;
1552
1553         if (file_version >= tpxv_PullCoordTypeGeom)
1554         {
1555             gmx_fio_do_gmx_bool(fio, ir->bPull);
1556         }
1557         else
1558         {
1559             gmx_fio_do_int(fio, ePullOld);
1560             ir->bPull = (ePullOld > 0);
1561             /* We removed the first ePull=ePullNo for the enum */
1562             ePullOld -= 1;
1563         }
1564         if (ir->bPull)
1565         {
1566             if (bRead)
1567             {
1568                 snew(ir->pull, 1);
1569             }
1570             do_pull(fio, ir->pull, bRead, file_version, ePullOld);
1571         }
1572     }
1573     else
1574     {
1575         ir->bPull = FALSE;
1576     }
1577
1578     /* Enforced rotation */
1579     if (file_version >= 74)
1580     {
1581         gmx_fio_do_int(fio, ir->bRot);
1582         if (ir->bRot == TRUE)
1583         {
1584             if (bRead)
1585             {
1586                 snew(ir->rot, 1);
1587             }
1588             do_rot(fio, ir->rot, bRead);
1589         }
1590     }
1591     else
1592     {
1593         ir->bRot = FALSE;
1594     }
1595
1596     /* Interactive molecular dynamics */
1597     if (file_version >= tpxv_InteractiveMolecularDynamics)
1598     {
1599         gmx_fio_do_int(fio, ir->bIMD);
1600         if (TRUE == ir->bIMD)
1601         {
1602             if (bRead)
1603             {
1604                 snew(ir->imd, 1);
1605             }
1606             do_imd(fio, ir->imd, bRead);
1607         }
1608     }
1609     else
1610     {
1611         /* We don't support IMD sessions for old .tpr files */
1612         ir->bIMD = FALSE;
1613     }
1614
1615     /* grpopts stuff */
1616     gmx_fio_do_int(fio, ir->opts.ngtc);
1617     if (file_version >= 69)
1618     {
1619         gmx_fio_do_int(fio, ir->opts.nhchainlength);
1620     }
1621     else
1622     {
1623         ir->opts.nhchainlength = 1;
1624     }
1625     gmx_fio_do_int(fio, ir->opts.ngacc);
1626     gmx_fio_do_int(fio, ir->opts.ngfrz);
1627     gmx_fio_do_int(fio, ir->opts.ngener);
1628
1629     if (bRead)
1630     {
1631         snew(ir->opts.nrdf,   ir->opts.ngtc);
1632         snew(ir->opts.ref_t,  ir->opts.ngtc);
1633         snew(ir->opts.annealing, ir->opts.ngtc);
1634         snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1635         snew(ir->opts.anneal_time, ir->opts.ngtc);
1636         snew(ir->opts.anneal_temp, ir->opts.ngtc);
1637         snew(ir->opts.tau_t,  ir->opts.ngtc);
1638         snew(ir->opts.nFreeze, ir->opts.ngfrz);
1639         snew(ir->opts.acc,    ir->opts.ngacc);
1640         snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1641     }
1642     if (ir->opts.ngtc > 0)
1643     {
1644         if (bRead && file_version < 13)
1645         {
1646             snew(tmp, ir->opts.ngtc);
1647             gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1648             for (i = 0; i < ir->opts.ngtc; i++)
1649             {
1650                 ir->opts.nrdf[i] = tmp[i];
1651             }
1652             sfree(tmp);
1653         }
1654         else
1655         {
1656             gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1657         }
1658         gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1659         gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1660         if (file_version < 33 && ir->eI == eiBD)
1661         {
1662             for (i = 0; i < ir->opts.ngtc; i++)
1663             {
1664                 ir->opts.tau_t[i] = bd_temp;
1665             }
1666         }
1667     }
1668     if (ir->opts.ngfrz > 0)
1669     {
1670         gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1671     }
1672     if (ir->opts.ngacc > 0)
1673     {
1674         gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1675     }
1676     if (file_version >= 12)
1677     {
1678         gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1679                         ir->opts.ngener*ir->opts.ngener);
1680     }
1681
1682     if (bRead && file_version < 26)
1683     {
1684         for (i = 0; i < ir->opts.ngtc; i++)
1685         {
1686             if (bSimAnn)
1687             {
1688                 ir->opts.annealing[i]      = eannSINGLE;
1689                 ir->opts.anneal_npoints[i] = 2;
1690                 snew(ir->opts.anneal_time[i], 2);
1691                 snew(ir->opts.anneal_temp[i], 2);
1692                 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1693                 finish_t                   = ir->init_t + ir->nsteps * ir->delta_t;
1694                 init_temp                  = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1695                 finish_temp                = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1696                 ir->opts.anneal_time[i][0] = ir->init_t;
1697                 ir->opts.anneal_time[i][1] = finish_t;
1698                 ir->opts.anneal_temp[i][0] = init_temp;
1699                 ir->opts.anneal_temp[i][1] = finish_temp;
1700             }
1701             else
1702             {
1703                 ir->opts.annealing[i]      = eannNO;
1704                 ir->opts.anneal_npoints[i] = 0;
1705             }
1706         }
1707     }
1708     else
1709     {
1710         /* file version 26 or later */
1711         /* First read the lists with annealing and npoints for each group */
1712         gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1713         gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1714         for (j = 0; j < (ir->opts.ngtc); j++)
1715         {
1716             k = ir->opts.anneal_npoints[j];
1717             if (bRead)
1718             {
1719                 snew(ir->opts.anneal_time[j], k);
1720                 snew(ir->opts.anneal_temp[j], k);
1721             }
1722             gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1723             gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1724         }
1725     }
1726     /* Walls */
1727     if (file_version >= 45)
1728     {
1729         gmx_fio_do_int(fio, ir->nwall);
1730         gmx_fio_do_int(fio, ir->wall_type);
1731         if (file_version >= 50)
1732         {
1733             gmx_fio_do_real(fio, ir->wall_r_linpot);
1734         }
1735         else
1736         {
1737             ir->wall_r_linpot = -1;
1738         }
1739         gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1740         gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1741         gmx_fio_do_real(fio, ir->wall_density[0]);
1742         gmx_fio_do_real(fio, ir->wall_density[1]);
1743         gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1744     }
1745     else
1746     {
1747         ir->nwall            = 0;
1748         ir->wall_type        = 0;
1749         ir->wall_atomtype[0] = -1;
1750         ir->wall_atomtype[1] = -1;
1751         ir->wall_density[0]  = 0;
1752         ir->wall_density[1]  = 0;
1753         ir->wall_ewald_zfac  = 3;
1754     }
1755     /* Cosine stuff for electric fields */
1756     for (j = 0; (j < DIM); j++)
1757     {
1758         gmx_fio_do_int(fio, ir->ex[j].n);
1759         gmx_fio_do_int(fio, ir->et[j].n);
1760         if (bRead)
1761         {
1762             snew(ir->ex[j].a,  ir->ex[j].n);
1763             snew(ir->ex[j].phi, ir->ex[j].n);
1764             snew(ir->et[j].a,  ir->et[j].n);
1765             snew(ir->et[j].phi, ir->et[j].n);
1766         }
1767         gmx_fio_ndo_real(fio, ir->ex[j].a,  ir->ex[j].n);
1768         gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1769         gmx_fio_ndo_real(fio, ir->et[j].a,  ir->et[j].n);
1770         gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1771     }
1772
1773     /* Swap ions */
1774     if (file_version >= tpxv_ComputationalElectrophysiology)
1775     {
1776         gmx_fio_do_int(fio, ir->eSwapCoords);
1777         if (ir->eSwapCoords != eswapNO)
1778         {
1779             if (bRead)
1780             {
1781                 snew(ir->swap, 1);
1782             }
1783             do_swapcoords(fio, ir->swap, bRead);
1784         }
1785     }
1786
1787     /* QMMM stuff */
1788     if (file_version >= 39)
1789     {
1790         gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1791         gmx_fio_do_int(fio, ir->QMMMscheme);
1792         gmx_fio_do_real(fio, ir->scalefactor);
1793         gmx_fio_do_int(fio, ir->opts.ngQM);
1794         if (bRead)
1795         {
1796             snew(ir->opts.QMmethod,    ir->opts.ngQM);
1797             snew(ir->opts.QMbasis,     ir->opts.ngQM);
1798             snew(ir->opts.QMcharge,    ir->opts.ngQM);
1799             snew(ir->opts.QMmult,      ir->opts.ngQM);
1800             snew(ir->opts.bSH,         ir->opts.ngQM);
1801             snew(ir->opts.CASorbitals, ir->opts.ngQM);
1802             snew(ir->opts.CASelectrons, ir->opts.ngQM);
1803             snew(ir->opts.SAon,        ir->opts.ngQM);
1804             snew(ir->opts.SAoff,       ir->opts.ngQM);
1805             snew(ir->opts.SAsteps,     ir->opts.ngQM);
1806             snew(ir->opts.bOPT,        ir->opts.ngQM);
1807             snew(ir->opts.bTS,         ir->opts.ngQM);
1808         }
1809         if (ir->opts.ngQM > 0)
1810         {
1811             gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1812             gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1813             gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1814             gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1815             gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1816             gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1817             gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1818             gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1819             gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1820             gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1821             gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1822             gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1823         }
1824         /* end of QMMM stuff */
1825     }
1826 }
1827
1828
1829 static void do_harm(t_fileio *fio, t_iparams *iparams)
1830 {
1831     gmx_fio_do_real(fio, iparams->harmonic.rA);
1832     gmx_fio_do_real(fio, iparams->harmonic.krA);
1833     gmx_fio_do_real(fio, iparams->harmonic.rB);
1834     gmx_fio_do_real(fio, iparams->harmonic.krB);
1835 }
1836
1837 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1838                 gmx_bool bRead, int file_version)
1839 {
1840     int      idum;
1841     real     rdum;
1842
1843     if (!bRead)
1844     {
1845         gmx_fio_set_comment(fio, interaction_function[ftype].name);
1846     }
1847     switch (ftype)
1848     {
1849         case F_ANGLES:
1850         case F_G96ANGLES:
1851         case F_BONDS:
1852         case F_G96BONDS:
1853         case F_HARMONIC:
1854         case F_IDIHS:
1855             do_harm(fio, iparams);
1856             if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1857             {
1858                 /* Correct incorrect storage of parameters */
1859                 iparams->pdihs.phiB = iparams->pdihs.phiA;
1860                 iparams->pdihs.cpB  = iparams->pdihs.cpA;
1861             }
1862             break;
1863         case F_RESTRANGLES:
1864             gmx_fio_do_real(fio, iparams->harmonic.rA);
1865             gmx_fio_do_real(fio, iparams->harmonic.krA);
1866             break;
1867         case F_LINEAR_ANGLES:
1868             gmx_fio_do_real(fio, iparams->linangle.klinA);
1869             gmx_fio_do_real(fio, iparams->linangle.aA);
1870             gmx_fio_do_real(fio, iparams->linangle.klinB);
1871             gmx_fio_do_real(fio, iparams->linangle.aB);
1872             break;
1873         case F_FENEBONDS:
1874             gmx_fio_do_real(fio, iparams->fene.bm);
1875             gmx_fio_do_real(fio, iparams->fene.kb);
1876             break;
1877
1878         case F_RESTRBONDS:
1879             gmx_fio_do_real(fio, iparams->restraint.lowA);
1880             gmx_fio_do_real(fio, iparams->restraint.up1A);
1881             gmx_fio_do_real(fio, iparams->restraint.up2A);
1882             gmx_fio_do_real(fio, iparams->restraint.kA);
1883             gmx_fio_do_real(fio, iparams->restraint.lowB);
1884             gmx_fio_do_real(fio, iparams->restraint.up1B);
1885             gmx_fio_do_real(fio, iparams->restraint.up2B);
1886             gmx_fio_do_real(fio, iparams->restraint.kB);
1887             break;
1888         case F_TABBONDS:
1889         case F_TABBONDSNC:
1890         case F_TABANGLES:
1891         case F_TABDIHS:
1892             gmx_fio_do_real(fio, iparams->tab.kA);
1893             gmx_fio_do_int(fio, iparams->tab.table);
1894             gmx_fio_do_real(fio, iparams->tab.kB);
1895             break;
1896         case F_CROSS_BOND_BONDS:
1897             gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1898             gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1899             gmx_fio_do_real(fio, iparams->cross_bb.krr);
1900             break;
1901         case F_CROSS_BOND_ANGLES:
1902             gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1903             gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1904             gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1905             gmx_fio_do_real(fio, iparams->cross_ba.krt);
1906             break;
1907         case F_UREY_BRADLEY:
1908             gmx_fio_do_real(fio, iparams->u_b.thetaA);
1909             gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1910             gmx_fio_do_real(fio, iparams->u_b.r13A);
1911             gmx_fio_do_real(fio, iparams->u_b.kUBA);
1912             if (file_version >= 79)
1913             {
1914                 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1915                 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1916                 gmx_fio_do_real(fio, iparams->u_b.r13B);
1917                 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1918             }
1919             else
1920             {
1921                 iparams->u_b.thetaB  = iparams->u_b.thetaA;
1922                 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1923                 iparams->u_b.r13B    = iparams->u_b.r13A;
1924                 iparams->u_b.kUBB    = iparams->u_b.kUBA;
1925             }
1926             break;
1927         case F_QUARTIC_ANGLES:
1928             gmx_fio_do_real(fio, iparams->qangle.theta);
1929             gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1930             break;
1931         case F_BHAM:
1932             gmx_fio_do_real(fio, iparams->bham.a);
1933             gmx_fio_do_real(fio, iparams->bham.b);
1934             gmx_fio_do_real(fio, iparams->bham.c);
1935             break;
1936         case F_MORSE:
1937             gmx_fio_do_real(fio, iparams->morse.b0A);
1938             gmx_fio_do_real(fio, iparams->morse.cbA);
1939             gmx_fio_do_real(fio, iparams->morse.betaA);
1940             if (file_version >= 79)
1941             {
1942                 gmx_fio_do_real(fio, iparams->morse.b0B);
1943                 gmx_fio_do_real(fio, iparams->morse.cbB);
1944                 gmx_fio_do_real(fio, iparams->morse.betaB);
1945             }
1946             else
1947             {
1948                 iparams->morse.b0B   = iparams->morse.b0A;
1949                 iparams->morse.cbB   = iparams->morse.cbA;
1950                 iparams->morse.betaB = iparams->morse.betaA;
1951             }
1952             break;
1953         case F_CUBICBONDS:
1954             gmx_fio_do_real(fio, iparams->cubic.b0);
1955             gmx_fio_do_real(fio, iparams->cubic.kb);
1956             gmx_fio_do_real(fio, iparams->cubic.kcub);
1957             break;
1958         case F_CONNBONDS:
1959             break;
1960         case F_POLARIZATION:
1961             gmx_fio_do_real(fio, iparams->polarize.alpha);
1962             break;
1963         case F_ANHARM_POL:
1964             gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1965             gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1966             gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1967             break;
1968         case F_WATER_POL:
1969             if (file_version < 31)
1970             {
1971                 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1972             }
1973             gmx_fio_do_real(fio, iparams->wpol.al_x);
1974             gmx_fio_do_real(fio, iparams->wpol.al_y);
1975             gmx_fio_do_real(fio, iparams->wpol.al_z);
1976             gmx_fio_do_real(fio, iparams->wpol.rOH);
1977             gmx_fio_do_real(fio, iparams->wpol.rHH);
1978             gmx_fio_do_real(fio, iparams->wpol.rOD);
1979             break;
1980         case F_THOLE_POL:
1981             gmx_fio_do_real(fio, iparams->thole.a);
1982             gmx_fio_do_real(fio, iparams->thole.alpha1);
1983             gmx_fio_do_real(fio, iparams->thole.alpha2);
1984             gmx_fio_do_real(fio, iparams->thole.rfac);
1985             break;
1986         case F_LJ:
1987             gmx_fio_do_real(fio, iparams->lj.c6);
1988             gmx_fio_do_real(fio, iparams->lj.c12);
1989             break;
1990         case F_LJ14:
1991             gmx_fio_do_real(fio, iparams->lj14.c6A);
1992             gmx_fio_do_real(fio, iparams->lj14.c12A);
1993             gmx_fio_do_real(fio, iparams->lj14.c6B);
1994             gmx_fio_do_real(fio, iparams->lj14.c12B);
1995             break;
1996         case F_LJC14_Q:
1997             gmx_fio_do_real(fio, iparams->ljc14.fqq);
1998             gmx_fio_do_real(fio, iparams->ljc14.qi);
1999             gmx_fio_do_real(fio, iparams->ljc14.qj);
2000             gmx_fio_do_real(fio, iparams->ljc14.c6);
2001             gmx_fio_do_real(fio, iparams->ljc14.c12);
2002             break;
2003         case F_LJC_PAIRS_NB:
2004             gmx_fio_do_real(fio, iparams->ljcnb.qi);
2005             gmx_fio_do_real(fio, iparams->ljcnb.qj);
2006             gmx_fio_do_real(fio, iparams->ljcnb.c6);
2007             gmx_fio_do_real(fio, iparams->ljcnb.c12);
2008             break;
2009         case F_PDIHS:
2010         case F_PIDIHS:
2011         case F_ANGRES:
2012         case F_ANGRESZ:
2013             gmx_fio_do_real(fio, iparams->pdihs.phiA);
2014             gmx_fio_do_real(fio, iparams->pdihs.cpA);
2015             if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
2016             {
2017                 /* Read the incorrectly stored multiplicity */
2018                 gmx_fio_do_real(fio, iparams->harmonic.rB);
2019                 gmx_fio_do_real(fio, iparams->harmonic.krB);
2020                 iparams->pdihs.phiB = iparams->pdihs.phiA;
2021                 iparams->pdihs.cpB  = iparams->pdihs.cpA;
2022             }
2023             else
2024             {
2025                 gmx_fio_do_real(fio, iparams->pdihs.phiB);
2026                 gmx_fio_do_real(fio, iparams->pdihs.cpB);
2027                 gmx_fio_do_int(fio, iparams->pdihs.mult);
2028             }
2029             break;
2030         case F_RESTRDIHS:
2031             gmx_fio_do_real(fio, iparams->pdihs.phiA);
2032             gmx_fio_do_real(fio, iparams->pdihs.cpA);
2033             break;
2034         case F_DISRES:
2035             gmx_fio_do_int(fio, iparams->disres.label);
2036             gmx_fio_do_int(fio, iparams->disres.type);
2037             gmx_fio_do_real(fio, iparams->disres.low);
2038             gmx_fio_do_real(fio, iparams->disres.up1);
2039             gmx_fio_do_real(fio, iparams->disres.up2);
2040             gmx_fio_do_real(fio, iparams->disres.kfac);
2041             break;
2042         case F_ORIRES:
2043             gmx_fio_do_int(fio, iparams->orires.ex);
2044             gmx_fio_do_int(fio, iparams->orires.label);
2045             gmx_fio_do_int(fio, iparams->orires.power);
2046             gmx_fio_do_real(fio, iparams->orires.c);
2047             gmx_fio_do_real(fio, iparams->orires.obs);
2048             gmx_fio_do_real(fio, iparams->orires.kfac);
2049             break;
2050         case F_DIHRES:
2051             if (file_version < 82)
2052             {
2053                 gmx_fio_do_int(fio, idum);
2054                 gmx_fio_do_int(fio, idum);
2055             }
2056             gmx_fio_do_real(fio, iparams->dihres.phiA);
2057             gmx_fio_do_real(fio, iparams->dihres.dphiA);
2058             gmx_fio_do_real(fio, iparams->dihres.kfacA);
2059             if (file_version >= 82)
2060             {
2061                 gmx_fio_do_real(fio, iparams->dihres.phiB);
2062                 gmx_fio_do_real(fio, iparams->dihres.dphiB);
2063                 gmx_fio_do_real(fio, iparams->dihres.kfacB);
2064             }
2065             else
2066             {
2067                 iparams->dihres.phiB  = iparams->dihres.phiA;
2068                 iparams->dihres.dphiB = iparams->dihres.dphiA;
2069                 iparams->dihres.kfacB = iparams->dihres.kfacA;
2070             }
2071             break;
2072         case F_POSRES:
2073             gmx_fio_do_rvec(fio, iparams->posres.pos0A);
2074             gmx_fio_do_rvec(fio, iparams->posres.fcA);
2075             if (bRead && file_version < 27)
2076             {
2077                 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
2078                 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
2079             }
2080             else
2081             {
2082                 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
2083                 gmx_fio_do_rvec(fio, iparams->posres.fcB);
2084             }
2085             break;
2086         case F_FBPOSRES:
2087             gmx_fio_do_int(fio, iparams->fbposres.geom);
2088             gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
2089             gmx_fio_do_real(fio, iparams->fbposres.r);
2090             gmx_fio_do_real(fio, iparams->fbposres.k);
2091             break;
2092         case F_CBTDIHS:
2093             gmx_fio_ndo_real(fio, iparams->cbtdihs.cbtcA, NR_CBTDIHS);
2094             break;
2095         case F_RBDIHS:
2096             gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2097             if (file_version >= 25)
2098             {
2099                 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2100             }
2101             break;
2102         case F_FOURDIHS:
2103             /* Fourier dihedrals are internally represented
2104              * as Ryckaert-Bellemans since those are faster to compute.
2105              */
2106             gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2107             gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2108             break;
2109         case F_CONSTR:
2110         case F_CONSTRNC:
2111             gmx_fio_do_real(fio, iparams->constr.dA);
2112             gmx_fio_do_real(fio, iparams->constr.dB);
2113             break;
2114         case F_SETTLE:
2115             gmx_fio_do_real(fio, iparams->settle.doh);
2116             gmx_fio_do_real(fio, iparams->settle.dhh);
2117             break;
2118         case F_VSITE2:
2119             gmx_fio_do_real(fio, iparams->vsite.a);
2120             break;
2121         case F_VSITE3:
2122         case F_VSITE3FD:
2123         case F_VSITE3FAD:
2124             gmx_fio_do_real(fio, iparams->vsite.a);
2125             gmx_fio_do_real(fio, iparams->vsite.b);
2126             break;
2127         case F_VSITE3OUT:
2128         case F_VSITE4FD:
2129         case F_VSITE4FDN:
2130             gmx_fio_do_real(fio, iparams->vsite.a);
2131             gmx_fio_do_real(fio, iparams->vsite.b);
2132             gmx_fio_do_real(fio, iparams->vsite.c);
2133             break;
2134         case F_VSITEN:
2135             gmx_fio_do_int(fio, iparams->vsiten.n);
2136             gmx_fio_do_real(fio, iparams->vsiten.a);
2137             break;
2138         case F_GB12:
2139         case F_GB13:
2140         case F_GB14:
2141             /* We got rid of some parameters in version 68 */
2142             if (bRead && file_version < 68)
2143             {
2144                 gmx_fio_do_real(fio, rdum);
2145                 gmx_fio_do_real(fio, rdum);
2146                 gmx_fio_do_real(fio, rdum);
2147                 gmx_fio_do_real(fio, rdum);
2148             }
2149             gmx_fio_do_real(fio, iparams->gb.sar);
2150             gmx_fio_do_real(fio, iparams->gb.st);
2151             gmx_fio_do_real(fio, iparams->gb.pi);
2152             gmx_fio_do_real(fio, iparams->gb.gbr);
2153             gmx_fio_do_real(fio, iparams->gb.bmlt);
2154             break;
2155         case F_CMAP:
2156             gmx_fio_do_int(fio, iparams->cmap.cmapA);
2157             gmx_fio_do_int(fio, iparams->cmap.cmapB);
2158             break;
2159         default:
2160             gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
2161                       ftype, interaction_function[ftype].name, __FILE__, __LINE__);
2162     }
2163     if (!bRead)
2164     {
2165         gmx_fio_unset_comment(fio);
2166     }
2167 }
2168
2169 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
2170                      int ftype)
2171 {
2172     int      i, k, idum;
2173
2174     if (!bRead)
2175     {
2176         gmx_fio_set_comment(fio, interaction_function[ftype].name);
2177     }
2178     if (file_version < 44)
2179     {
2180         for (i = 0; i < MAXNODES; i++)
2181         {
2182             gmx_fio_do_int(fio, idum);
2183         }
2184     }
2185     gmx_fio_do_int(fio, ilist->nr);
2186     if (bRead)
2187     {
2188         snew(ilist->iatoms, ilist->nr);
2189     }
2190     gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2191     if (!bRead)
2192     {
2193         gmx_fio_unset_comment(fio);
2194     }
2195 }
2196
2197 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2198                         gmx_bool bRead, int file_version)
2199 {
2200     int          idum, i, j;
2201     unsigned int k;
2202
2203     gmx_fio_do_int(fio, ffparams->atnr);
2204     if (file_version < 57)
2205     {
2206         gmx_fio_do_int(fio, idum);
2207     }
2208     gmx_fio_do_int(fio, ffparams->ntypes);
2209     if (bRead && debug)
2210     {
2211         fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2212                 ffparams->atnr, ffparams->ntypes);
2213     }
2214     if (bRead)
2215     {
2216         snew(ffparams->functype, ffparams->ntypes);
2217         snew(ffparams->iparams, ffparams->ntypes);
2218     }
2219     /* Read/write all the function types */
2220     gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2221     if (bRead && debug)
2222     {
2223         pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2224     }
2225
2226     if (file_version >= 66)
2227     {
2228         gmx_fio_do_double(fio, ffparams->reppow);
2229     }
2230     else
2231     {
2232         ffparams->reppow = 12.0;
2233     }
2234
2235     if (file_version >= 57)
2236     {
2237         gmx_fio_do_real(fio, ffparams->fudgeQQ);
2238     }
2239
2240     /* Check whether all these function types are supported by the code.
2241      * In practice the code is backwards compatible, which means that the
2242      * numbering may have to be altered from old numbering to new numbering
2243      */
2244     for (i = 0; (i < ffparams->ntypes); i++)
2245     {
2246         if (bRead)
2247         {
2248             /* Loop over file versions */
2249             for (k = 0; (k < NFTUPD); k++)
2250             {
2251                 /* Compare the read file_version to the update table */
2252                 if ((file_version < ftupd[k].fvnr) &&
2253                     (ffparams->functype[i] >= ftupd[k].ftype))
2254                 {
2255                     ffparams->functype[i] += 1;
2256                     if (debug)
2257                     {
2258                         fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2259                                 i, ffparams->functype[i],
2260                                 interaction_function[ftupd[k].ftype].longname);
2261                         fflush(debug);
2262                     }
2263                 }
2264             }
2265         }
2266
2267         do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2268                    file_version);
2269         if (bRead && debug)
2270         {
2271             pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2272         }
2273     }
2274 }
2275
2276 static void add_settle_atoms(t_ilist *ilist)
2277 {
2278     int i;
2279
2280     /* Settle used to only store the first atom: add the other two */
2281     srenew(ilist->iatoms, 2*ilist->nr);
2282     for (i = ilist->nr/2-1; i >= 0; i--)
2283     {
2284         ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2285         ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2286         ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2287         ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2288     }
2289     ilist->nr = 2*ilist->nr;
2290 }
2291
2292 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2293                       int file_version)
2294 {
2295     int          i, j, renum[F_NRE];
2296     gmx_bool     bClear;
2297     unsigned int k;
2298
2299     for (j = 0; (j < F_NRE); j++)
2300     {
2301         bClear = FALSE;
2302         if (bRead)
2303         {
2304             for (k = 0; k < NFTUPD; k++)
2305             {
2306                 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2307                 {
2308                     bClear = TRUE;
2309                 }
2310             }
2311         }
2312         if (bClear)
2313         {
2314             ilist[j].nr     = 0;
2315             ilist[j].iatoms = NULL;
2316         }
2317         else
2318         {
2319             do_ilist(fio, &ilist[j], bRead, file_version, j);
2320             if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2321             {
2322                 add_settle_atoms(&ilist[j]);
2323             }
2324         }
2325         /*
2326            if (bRead && gmx_debug_at)
2327            pr_ilist(debug,0,interaction_function[j].longname,
2328                functype,&ilist[j],TRUE);
2329          */
2330     }
2331 }
2332
2333 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2334                     gmx_bool bRead, int file_version)
2335 {
2336     do_ffparams(fio, ffparams, bRead, file_version);
2337
2338     if (file_version >= 54)
2339     {
2340         gmx_fio_do_real(fio, ffparams->fudgeQQ);
2341     }
2342
2343     do_ilists(fio, molt->ilist, bRead, file_version);
2344 }
2345
2346 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2347 {
2348     int      i, idum, dum_nra, *dum_a;
2349
2350     if (file_version < 44)
2351     {
2352         for (i = 0; i < MAXNODES; i++)
2353         {
2354             gmx_fio_do_int(fio, idum);
2355         }
2356     }
2357     gmx_fio_do_int(fio, block->nr);
2358     if (file_version < 51)
2359     {
2360         gmx_fio_do_int(fio, dum_nra);
2361     }
2362     if (bRead)
2363     {
2364         if ((block->nalloc_index > 0) && (NULL != block->index))
2365         {
2366             sfree(block->index);
2367         }
2368         block->nalloc_index = block->nr+1;
2369         snew(block->index, block->nalloc_index);
2370     }
2371     gmx_fio_ndo_int(fio, block->index, block->nr+1);
2372
2373     if (file_version < 51 && dum_nra > 0)
2374     {
2375         snew(dum_a, dum_nra);
2376         gmx_fio_ndo_int(fio, dum_a, dum_nra);
2377         sfree(dum_a);
2378     }
2379 }
2380
2381 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2382                       int file_version)
2383 {
2384     int      i, idum;
2385
2386     if (file_version < 44)
2387     {
2388         for (i = 0; i < MAXNODES; i++)
2389         {
2390             gmx_fio_do_int(fio, idum);
2391         }
2392     }
2393     gmx_fio_do_int(fio, block->nr);
2394     gmx_fio_do_int(fio, block->nra);
2395     if (bRead)
2396     {
2397         block->nalloc_index = block->nr+1;
2398         snew(block->index, block->nalloc_index);
2399         block->nalloc_a = block->nra;
2400         snew(block->a, block->nalloc_a);
2401     }
2402     gmx_fio_ndo_int(fio, block->index, block->nr+1);
2403     gmx_fio_ndo_int(fio, block->a, block->nra);
2404 }
2405
2406 /* This is a primitive routine to make it possible to translate atomic numbers
2407  * to element names when reading TPR files, without making the Gromacs library
2408  * directory a dependency on mdrun (which is the case if we need elements.dat).
2409  */
2410 static const char *
2411 atomicnumber_to_element(int atomicnumber)
2412 {
2413     const char * p;
2414
2415     /* This does not have to be complete, so we only include elements likely
2416      * to occur in PDB files.
2417      */
2418     switch (atomicnumber)
2419     {
2420         case 1:  p = "H";  break;
2421         case 5:  p = "B";  break;
2422         case 6:  p = "C";  break;
2423         case 7:  p = "N";  break;
2424         case 8:  p = "O";  break;
2425         case 9:  p = "F";  break;
2426         case 11: p = "Na"; break;
2427         case 12: p = "Mg"; break;
2428         case 15: p = "P";  break;
2429         case 16: p = "S";  break;
2430         case 17: p = "Cl"; break;
2431         case 18: p = "Ar"; break;
2432         case 19: p = "K";  break;
2433         case 20: p = "Ca"; break;
2434         case 25: p = "Mn"; break;
2435         case 26: p = "Fe"; break;
2436         case 28: p = "Ni"; break;
2437         case 29: p = "Cu"; break;
2438         case 30: p = "Zn"; break;
2439         case 35: p = "Br"; break;
2440         case 47: p = "Ag"; break;
2441         default: p = "";   break;
2442     }
2443     return p;
2444 }
2445
2446
2447 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2448                     int file_version, gmx_groups_t *groups, int atnr)
2449 {
2450     int    i, myngrp;
2451     char * p_elem;
2452
2453     gmx_fio_do_real(fio, atom->m);
2454     gmx_fio_do_real(fio, atom->q);
2455     gmx_fio_do_real(fio, atom->mB);
2456     gmx_fio_do_real(fio, atom->qB);
2457     gmx_fio_do_ushort(fio, atom->type);
2458     gmx_fio_do_ushort(fio, atom->typeB);
2459     gmx_fio_do_int(fio, atom->ptype);
2460     gmx_fio_do_int(fio, atom->resind);
2461     if (file_version >= 52)
2462     {
2463         gmx_fio_do_int(fio, atom->atomnumber);
2464         if (bRead)
2465         {
2466             /* Set element string from atomic number if present.
2467              * This routine returns an empty string if the name is not found.
2468              */
2469             strncpy(atom->elem, atomicnumber_to_element(atom->atomnumber), 4);
2470             /* avoid warnings about potentially unterminated string */
2471             atom->elem[3] = '\0';
2472         }
2473     }
2474     else if (bRead)
2475     {
2476         atom->atomnumber = NOTSET;
2477     }
2478     if (file_version < 23)
2479     {
2480         myngrp = 8;
2481     }
2482     else if (file_version < 39)
2483     {
2484         myngrp = 9;
2485     }
2486     else
2487     {
2488         myngrp = ngrp;
2489     }
2490
2491     if (file_version < 57)
2492     {
2493         unsigned char uchar[egcNR];
2494         gmx_fio_ndo_uchar(fio, uchar, myngrp);
2495         for (i = myngrp; (i < ngrp); i++)
2496         {
2497             uchar[i] = 0;
2498         }
2499         /* Copy the old data format to the groups struct */
2500         for (i = 0; i < ngrp; i++)
2501         {
2502             groups->grpnr[i][atnr] = uchar[i];
2503         }
2504     }
2505 }
2506
2507 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2508                     int file_version)
2509 {
2510     int      i, j, myngrp;
2511
2512     if (file_version < 23)
2513     {
2514         myngrp = 8;
2515     }
2516     else if (file_version < 39)
2517     {
2518         myngrp = 9;
2519     }
2520     else
2521     {
2522         myngrp = ngrp;
2523     }
2524
2525     for (j = 0; (j < ngrp); j++)
2526     {
2527         if (j < myngrp)
2528         {
2529             gmx_fio_do_int(fio, grps[j].nr);
2530             if (bRead)
2531             {
2532                 snew(grps[j].nm_ind, grps[j].nr);
2533             }
2534             gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2535         }
2536         else
2537         {
2538             grps[j].nr = 1;
2539             snew(grps[j].nm_ind, grps[j].nr);
2540         }
2541     }
2542 }
2543
2544 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2545 {
2546     int ls;
2547
2548     if (bRead)
2549     {
2550         gmx_fio_do_int(fio, ls);
2551         *nm = get_symtab_handle(symtab, ls);
2552     }
2553     else
2554     {
2555         ls = lookup_symtab(symtab, *nm);
2556         gmx_fio_do_int(fio, ls);
2557     }
2558 }
2559
2560 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2561                       t_symtab *symtab)
2562 {
2563     int  j;
2564
2565     for (j = 0; (j < nstr); j++)
2566     {
2567         do_symstr(fio, &(nm[j]), bRead, symtab);
2568     }
2569 }
2570
2571 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2572                        t_symtab *symtab, int file_version)
2573 {
2574     int  j;
2575
2576     for (j = 0; (j < n); j++)
2577     {
2578         do_symstr(fio, &(ri[j].name), bRead, symtab);
2579         if (file_version >= 63)
2580         {
2581             gmx_fio_do_int(fio, ri[j].nr);
2582             gmx_fio_do_uchar(fio, ri[j].ic);
2583         }
2584         else
2585         {
2586             ri[j].nr = j + 1;
2587             ri[j].ic = ' ';
2588         }
2589     }
2590 }
2591
2592 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2593                      int file_version,
2594                      gmx_groups_t *groups)
2595 {
2596     int i;
2597
2598     gmx_fio_do_int(fio, atoms->nr);
2599     gmx_fio_do_int(fio, atoms->nres);
2600     if (file_version < 57)
2601     {
2602         gmx_fio_do_int(fio, groups->ngrpname);
2603         for (i = 0; i < egcNR; i++)
2604         {
2605             groups->ngrpnr[i] = atoms->nr;
2606             snew(groups->grpnr[i], groups->ngrpnr[i]);
2607         }
2608     }
2609     if (bRead)
2610     {
2611         snew(atoms->atom, atoms->nr);
2612         snew(atoms->atomname, atoms->nr);
2613         snew(atoms->atomtype, atoms->nr);
2614         snew(atoms->atomtypeB, atoms->nr);
2615         snew(atoms->resinfo, atoms->nres);
2616         if (file_version < 57)
2617         {
2618             snew(groups->grpname, groups->ngrpname);
2619         }
2620         atoms->pdbinfo = NULL;
2621     }
2622     for (i = 0; (i < atoms->nr); i++)
2623     {
2624         do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2625     }
2626     do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2627     if (bRead && (file_version <= 20))
2628     {
2629         for (i = 0; i < atoms->nr; i++)
2630         {
2631             atoms->atomtype[i]  = put_symtab(symtab, "?");
2632             atoms->atomtypeB[i] = put_symtab(symtab, "?");
2633         }
2634     }
2635     else
2636     {
2637         do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2638         do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2639     }
2640     do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2641
2642     if (file_version < 57)
2643     {
2644         do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2645
2646         do_grps(fio, egcNR, groups->grps, bRead, file_version);
2647     }
2648 }
2649
2650 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2651                       gmx_bool bRead, t_symtab *symtab,
2652                       int file_version)
2653 {
2654     int      g, n, i;
2655
2656     do_grps(fio, egcNR, groups->grps, bRead, file_version);
2657     gmx_fio_do_int(fio, groups->ngrpname);
2658     if (bRead)
2659     {
2660         snew(groups->grpname, groups->ngrpname);
2661     }
2662     do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2663     for (g = 0; g < egcNR; g++)
2664     {
2665         gmx_fio_do_int(fio, groups->ngrpnr[g]);
2666         if (groups->ngrpnr[g] == 0)
2667         {
2668             if (bRead)
2669             {
2670                 groups->grpnr[g] = NULL;
2671             }
2672         }
2673         else
2674         {
2675             if (bRead)
2676             {
2677                 snew(groups->grpnr[g], groups->ngrpnr[g]);
2678             }
2679             gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2680         }
2681     }
2682 }
2683
2684 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2685                          int file_version)
2686 {
2687     int      i, j;
2688
2689     if (file_version > 25)
2690     {
2691         gmx_fio_do_int(fio, atomtypes->nr);
2692         j = atomtypes->nr;
2693         if (bRead)
2694         {
2695             snew(atomtypes->radius, j);
2696             snew(atomtypes->vol, j);
2697             snew(atomtypes->surftens, j);
2698             snew(atomtypes->atomnumber, j);
2699             snew(atomtypes->gb_radius, j);
2700             snew(atomtypes->S_hct, j);
2701         }
2702         gmx_fio_ndo_real(fio, atomtypes->radius, j);
2703         gmx_fio_ndo_real(fio, atomtypes->vol, j);
2704         gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2705         if (file_version >= 40)
2706         {
2707             gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2708         }
2709         if (file_version >= 60)
2710         {
2711             gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2712             gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2713         }
2714     }
2715     else
2716     {
2717         /* File versions prior to 26 cannot do GBSA,
2718          * so they dont use this structure
2719          */
2720         atomtypes->nr         = 0;
2721         atomtypes->radius     = NULL;
2722         atomtypes->vol        = NULL;
2723         atomtypes->surftens   = NULL;
2724         atomtypes->atomnumber = NULL;
2725         atomtypes->gb_radius  = NULL;
2726         atomtypes->S_hct      = NULL;
2727     }
2728 }
2729
2730 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2731 {
2732     int       i, nr;
2733     t_symbuf *symbuf;
2734     char      buf[STRLEN];
2735
2736     gmx_fio_do_int(fio, symtab->nr);
2737     nr     = symtab->nr;
2738     if (bRead)
2739     {
2740         snew(symtab->symbuf, 1);
2741         symbuf          = symtab->symbuf;
2742         symbuf->bufsize = nr;
2743         snew(symbuf->buf, nr);
2744         for (i = 0; (i < nr); i++)
2745         {
2746             gmx_fio_do_string(fio, buf);
2747             symbuf->buf[i] = gmx_strdup(buf);
2748         }
2749     }
2750     else
2751     {
2752         symbuf = symtab->symbuf;
2753         while (symbuf != NULL)
2754         {
2755             for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2756             {
2757                 gmx_fio_do_string(fio, symbuf->buf[i]);
2758             }
2759             nr    -= i;
2760             symbuf = symbuf->next;
2761         }
2762         if (nr != 0)
2763         {
2764             gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2765         }
2766     }
2767 }
2768
2769 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2770 {
2771     int i, j, ngrid, gs, nelem;
2772
2773     gmx_fio_do_int(fio, cmap_grid->ngrid);
2774     gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2775
2776     ngrid = cmap_grid->ngrid;
2777     gs    = cmap_grid->grid_spacing;
2778     nelem = gs * gs;
2779
2780     if (bRead)
2781     {
2782         snew(cmap_grid->cmapdata, ngrid);
2783
2784         for (i = 0; i < cmap_grid->ngrid; i++)
2785         {
2786             snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2787         }
2788     }
2789
2790     for (i = 0; i < cmap_grid->ngrid; i++)
2791     {
2792         for (j = 0; j < nelem; j++)
2793         {
2794             gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2795             gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2796             gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2797             gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2798         }
2799     }
2800 }
2801
2802
2803 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2804 {
2805     int  m, a, a0, a1, r;
2806     char c, chainid;
2807     int  chainnum;
2808
2809     /* We always assign a new chain number, but save the chain id characters
2810      * for larger molecules.
2811      */
2812 #define CHAIN_MIN_ATOMS 15
2813
2814     chainnum = 0;
2815     chainid  = 'A';
2816     for (m = 0; m < mols->nr; m++)
2817     {
2818         a0 = mols->index[m];
2819         a1 = mols->index[m+1];
2820         if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2821         {
2822             c = chainid;
2823             chainid++;
2824         }
2825         else
2826         {
2827             c = ' ';
2828         }
2829         for (a = a0; a < a1; a++)
2830         {
2831             atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2832             atoms->resinfo[atoms->atom[a].resind].chainid  = c;
2833         }
2834         chainnum++;
2835     }
2836
2837     /* Blank out the chain id if there was only one chain */
2838     if (chainid == 'B')
2839     {
2840         for (r = 0; r < atoms->nres; r++)
2841         {
2842             atoms->resinfo[r].chainid = ' ';
2843         }
2844     }
2845 }
2846
2847 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2848                        t_symtab *symtab, int file_version,
2849                        gmx_groups_t *groups)
2850 {
2851     int i;
2852
2853     if (file_version >= 57)
2854     {
2855         do_symstr(fio, &(molt->name), bRead, symtab);
2856     }
2857
2858     do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2859
2860     if (bRead && gmx_debug_at)
2861     {
2862         pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2863     }
2864
2865     if (file_version >= 57)
2866     {
2867         do_ilists(fio, molt->ilist, bRead, file_version);
2868
2869         do_block(fio, &molt->cgs, bRead, file_version);
2870         if (bRead && gmx_debug_at)
2871         {
2872             pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2873         }
2874     }
2875
2876     /* This used to be in the atoms struct */
2877     do_blocka(fio, &molt->excls, bRead, file_version);
2878 }
2879
2880 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2881 {
2882     int i;
2883
2884     gmx_fio_do_int(fio, molb->type);
2885     gmx_fio_do_int(fio, molb->nmol);
2886     gmx_fio_do_int(fio, molb->natoms_mol);
2887     /* Position restraint coordinates */
2888     gmx_fio_do_int(fio, molb->nposres_xA);
2889     if (molb->nposres_xA > 0)
2890     {
2891         if (bRead)
2892         {
2893             snew(molb->posres_xA, molb->nposres_xA);
2894         }
2895         gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2896     }
2897     gmx_fio_do_int(fio, molb->nposres_xB);
2898     if (molb->nposres_xB > 0)
2899     {
2900         if (bRead)
2901         {
2902             snew(molb->posres_xB, molb->nposres_xB);
2903         }
2904         gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2905     }
2906
2907 }
2908
2909 static t_block mtop_mols(gmx_mtop_t *mtop)
2910 {
2911     int     mb, m, a, mol;
2912     t_block mols;
2913
2914     mols.nr = 0;
2915     for (mb = 0; mb < mtop->nmolblock; mb++)
2916     {
2917         mols.nr += mtop->molblock[mb].nmol;
2918     }
2919     mols.nalloc_index = mols.nr + 1;
2920     snew(mols.index, mols.nalloc_index);
2921
2922     a             = 0;
2923     m             = 0;
2924     mols.index[m] = a;
2925     for (mb = 0; mb < mtop->nmolblock; mb++)
2926     {
2927         for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2928         {
2929             a += mtop->molblock[mb].natoms_mol;
2930             m++;
2931             mols.index[m] = a;
2932         }
2933     }
2934
2935     return mols;
2936 }
2937
2938 static void add_posres_molblock(gmx_mtop_t *mtop)
2939 {
2940     t_ilist        *il, *ilfb;
2941     int             am, i, mol, a;
2942     gmx_bool        bFE;
2943     gmx_molblock_t *molb;
2944     t_iparams      *ip;
2945
2946     /* posres reference positions are stored in ip->posres (if present) and
2947        in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2948        posres.pos0A are identical to fbposres.pos0. */
2949     il   = &mtop->moltype[0].ilist[F_POSRES];
2950     ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2951     if (il->nr == 0 && ilfb->nr == 0)
2952     {
2953         return;
2954     }
2955     am  = 0;
2956     bFE = FALSE;
2957     for (i = 0; i < il->nr; i += 2)
2958     {
2959         ip = &mtop->ffparams.iparams[il->iatoms[i]];
2960         am = max(am, il->iatoms[i+1]);
2961         if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2962             ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2963             ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2964         {
2965             bFE = TRUE;
2966         }
2967     }
2968     /* This loop is required if we have only flat-bottomed posres:
2969        - set am
2970        - bFE == FALSE (no B-state for flat-bottomed posres) */
2971     if (il->nr == 0)
2972     {
2973         for (i = 0; i < ilfb->nr; i += 2)
2974         {
2975             ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2976             am = max(am, ilfb->iatoms[i+1]);
2977         }
2978     }
2979     /* Make the posres coordinate block end at a molecule end */
2980     mol = 0;
2981     while (am >= mtop->mols.index[mol+1])
2982     {
2983         mol++;
2984     }
2985     molb             = &mtop->molblock[0];
2986     molb->nposres_xA = mtop->mols.index[mol+1];
2987     snew(molb->posres_xA, molb->nposres_xA);
2988     if (bFE)
2989     {
2990         molb->nposres_xB = molb->nposres_xA;
2991         snew(molb->posres_xB, molb->nposres_xB);
2992     }
2993     else
2994     {
2995         molb->nposres_xB = 0;
2996     }
2997     for (i = 0; i < il->nr; i += 2)
2998     {
2999         ip                     = &mtop->ffparams.iparams[il->iatoms[i]];
3000         a                      = il->iatoms[i+1];
3001         molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
3002         molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
3003         molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
3004         if (bFE)
3005         {
3006             molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
3007             molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
3008             molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
3009         }
3010     }
3011     if (il->nr == 0)
3012     {
3013         /* If only flat-bottomed posres are present, take reference pos from them.
3014            Here: bFE == FALSE      */
3015         for (i = 0; i < ilfb->nr; i += 2)
3016         {
3017             ip                     = &mtop->ffparams.iparams[ilfb->iatoms[i]];
3018             a                      = ilfb->iatoms[i+1];
3019             molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
3020             molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
3021             molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
3022         }
3023     }
3024 }
3025
3026 static void set_disres_npair(gmx_mtop_t *mtop)
3027 {
3028     t_iparams            *ip;
3029     gmx_mtop_ilistloop_t  iloop;
3030     t_ilist              *ilist, *il;
3031     int                   nmol, i, npair;
3032     t_iatom              *a;
3033
3034     ip = mtop->ffparams.iparams;
3035
3036     iloop     = gmx_mtop_ilistloop_init(mtop);
3037     while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
3038     {
3039         il = &ilist[F_DISRES];
3040
3041         if (il->nr > 0)
3042         {
3043             a     = il->iatoms;
3044             npair = 0;
3045             for (i = 0; i < il->nr; i += 3)
3046             {
3047                 npair++;
3048                 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
3049                 {
3050                     ip[a[i]].disres.npair = npair;
3051                     npair                 = 0;
3052                 }
3053             }
3054         }
3055     }
3056 }
3057
3058 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
3059                     int file_version)
3060 {
3061     int            mt, mb, i;
3062     t_blocka       dumb;
3063
3064     if (bRead)
3065     {
3066         init_mtop(mtop);
3067     }
3068     do_symtab(fio, &(mtop->symtab), bRead);
3069     if (bRead && debug)
3070     {
3071         pr_symtab(debug, 0, "symtab", &mtop->symtab);
3072     }
3073
3074     do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
3075
3076     if (file_version >= 57)
3077     {
3078         do_ffparams(fio, &mtop->ffparams, bRead, file_version);
3079
3080         gmx_fio_do_int(fio, mtop->nmoltype);
3081     }
3082     else
3083     {
3084         mtop->nmoltype = 1;
3085     }
3086     if (bRead)
3087     {
3088         snew(mtop->moltype, mtop->nmoltype);
3089         if (file_version < 57)
3090         {
3091             mtop->moltype[0].name = mtop->name;
3092         }
3093     }
3094     for (mt = 0; mt < mtop->nmoltype; mt++)
3095     {
3096         do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
3097                    &mtop->groups);
3098     }
3099
3100     if (file_version >= 57)
3101     {
3102         gmx_fio_do_int(fio, mtop->nmolblock);
3103     }
3104     else
3105     {
3106         mtop->nmolblock = 1;
3107     }
3108     if (bRead)
3109     {
3110         snew(mtop->molblock, mtop->nmolblock);
3111     }
3112     if (file_version >= 57)
3113     {
3114         for (mb = 0; mb < mtop->nmolblock; mb++)
3115         {
3116             do_molblock(fio, &mtop->molblock[mb], bRead);
3117         }
3118         gmx_fio_do_int(fio, mtop->natoms);
3119     }
3120     else
3121     {
3122         mtop->molblock[0].type       = 0;
3123         mtop->molblock[0].nmol       = 1;
3124         mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
3125         mtop->molblock[0].nposres_xA = 0;
3126         mtop->molblock[0].nposres_xB = 0;
3127     }
3128
3129     if (file_version >= tpxv_IntermolecularBondeds)
3130     {
3131         gmx_fio_do_gmx_bool(fio, mtop->bIntermolecularInteractions);
3132         if (mtop->bIntermolecularInteractions)
3133         {
3134             if (bRead)
3135             {
3136                 snew(mtop->intermolecular_ilist, F_NRE);
3137             }
3138             do_ilists(fio, mtop->intermolecular_ilist, bRead, file_version);
3139         }
3140     }
3141     else
3142     {
3143         mtop->bIntermolecularInteractions = FALSE;
3144     }
3145
3146     do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
3147     if (bRead && debug)
3148     {
3149         pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
3150     }
3151
3152     if (file_version < 57)
3153     {
3154         /* Debug statements are inside do_idef */
3155         do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
3156         mtop->natoms = mtop->moltype[0].atoms.nr;
3157     }
3158
3159     if (file_version >= 65)
3160     {
3161         do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
3162     }
3163     else
3164     {
3165         mtop->ffparams.cmap_grid.ngrid        = 0;
3166         mtop->ffparams.cmap_grid.grid_spacing = 0;
3167         mtop->ffparams.cmap_grid.cmapdata     = NULL;
3168     }
3169
3170     if (file_version >= 57)
3171     {
3172         do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
3173     }
3174
3175     if (file_version < 57)
3176     {
3177         do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
3178         if (bRead && gmx_debug_at)
3179         {
3180             pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
3181         }
3182         do_block(fio, &mtop->mols, bRead, file_version);
3183         /* Add the posres coordinates to the molblock */
3184         add_posres_molblock(mtop);
3185     }
3186     if (bRead)
3187     {
3188         if (file_version >= 57)
3189         {
3190             done_block(&mtop->mols);
3191             mtop->mols = mtop_mols(mtop);
3192         }
3193         if (gmx_debug_at)
3194         {
3195             pr_block(debug, 0, "mols", &mtop->mols, TRUE);
3196         }
3197     }
3198
3199     if (file_version < 51)
3200     {
3201         /* Here used to be the shake blocks */
3202         do_blocka(fio, &dumb, bRead, file_version);
3203         if (dumb.nr > 0)
3204         {
3205             sfree(dumb.index);
3206         }
3207         if (dumb.nra > 0)
3208         {
3209             sfree(dumb.a);
3210         }
3211     }
3212
3213     if (bRead)
3214     {
3215         close_symtab(&(mtop->symtab));
3216     }
3217 }
3218
3219 /* If TopOnlyOK is TRUE then we can read even future versions
3220  * of tpx files, provided the file_generation hasn't changed.
3221  * If it is FALSE, we need the inputrecord too, and bail out
3222  * if the file is newer than the program.
3223  *
3224  * The version and generation if the topology (see top of this file)
3225  * are returned in the two last arguments.
3226  *
3227  * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
3228  */
3229 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
3230                          gmx_bool TopOnlyOK, int *file_version,
3231                          int *file_generation)
3232 {
3233     char      buf[STRLEN];
3234     char      file_tag[STRLEN];
3235     gmx_bool  bDouble;
3236     int       precision;
3237     int       fver, fgen;
3238     int       idum = 0;
3239     real      rdum = 0;
3240
3241     gmx_fio_checktype(fio);
3242     gmx_fio_setdebug(fio, bDebugMode());
3243
3244     /* XDR binary topology file */
3245     precision = sizeof(real);
3246     if (bRead)
3247     {
3248         gmx_fio_do_string(fio, buf);
3249         if (strncmp(buf, "VERSION", 7))
3250         {
3251             gmx_fatal(FARGS, "Can not read file %s,\n"
3252                       "             this file is from a GROMACS version which is older than 2.0\n"
3253                       "             Make a new one with grompp or use a gro or pdb file, if possible",
3254                       gmx_fio_getname(fio));
3255         }
3256         gmx_fio_do_int(fio, precision);
3257         bDouble = (precision == sizeof(double));
3258         if ((precision != sizeof(float)) && !bDouble)
3259         {
3260             gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3261                       "instead of %d or %d",
3262                       gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3263         }
3264         gmx_fio_setprecision(fio, bDouble);
3265         fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3266                 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3267     }
3268     else
3269     {
3270         gmx_fio_write_string(fio, GromacsVersion());
3271         bDouble = (precision == sizeof(double));
3272         gmx_fio_setprecision(fio, bDouble);
3273         gmx_fio_do_int(fio, precision);
3274         fver = tpx_version;
3275         sprintf(file_tag, "%s", tpx_tag);
3276         fgen = tpx_generation;
3277     }
3278
3279     /* Check versions! */
3280     gmx_fio_do_int(fio, fver);
3281
3282     /* This is for backward compatibility with development versions 77-79
3283      * where the tag was, mistakenly, placed before the generation,
3284      * which would cause a segv instead of a proper error message
3285      * when reading the topology only from tpx with <77 code.
3286      */
3287     if (fver >= 77 && fver <= 79)
3288     {
3289         gmx_fio_do_string(fio, file_tag);
3290     }
3291
3292     if (fver >= 26)
3293     {
3294         gmx_fio_do_int(fio, fgen);
3295     }
3296     else
3297     {
3298         fgen = 0;
3299     }
3300
3301     if (fver >= 81)
3302     {
3303         gmx_fio_do_string(fio, file_tag);
3304     }
3305     if (bRead)
3306     {
3307         if (fver < 77)
3308         {
3309             /* Versions before 77 don't have the tag, set it to release */
3310             sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3311         }
3312
3313         if (strcmp(file_tag, tpx_tag) != 0)
3314         {
3315             fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3316                     file_tag, tpx_tag);
3317
3318             /* We only support reading tpx files with the same tag as the code
3319              * or tpx files with the release tag and with lower version number.
3320              */
3321             if (strcmp(file_tag, TPX_TAG_RELEASE) != 0 && fver < tpx_version)
3322             {
3323                 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3324                           gmx_fio_getname(fio), fver, file_tag,
3325                           tpx_version, tpx_tag);
3326             }
3327         }
3328     }
3329
3330     if (file_version != NULL)
3331     {
3332         *file_version = fver;
3333     }
3334     if (file_generation != NULL)
3335     {
3336         *file_generation = fgen;
3337     }
3338
3339
3340     if ((fver <= tpx_incompatible_version) ||
3341         ((fver > tpx_version) && !TopOnlyOK) ||
3342         (fgen > tpx_generation) ||
3343         tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3344     {
3345         gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3346                   gmx_fio_getname(fio), fver, tpx_version);
3347     }
3348
3349     gmx_fio_do_int(fio, tpx->natoms);
3350     if (fver >= 28)
3351     {
3352         gmx_fio_do_int(fio, tpx->ngtc);
3353     }
3354     else
3355     {
3356         tpx->ngtc = 0;
3357     }
3358     if (fver < 62)
3359     {
3360         gmx_fio_do_int(fio, idum);
3361         gmx_fio_do_real(fio, rdum);
3362     }
3363     /*a better decision will eventually (5.0 or later) need to be made
3364        on how to treat the alchemical state of the system, which can now
3365        vary through a simulation, and cannot be completely described
3366        though a single lambda variable, or even a single state
3367        index. Eventually, should probably be a vector. MRS*/
3368     if (fver >= 79)
3369     {
3370         gmx_fio_do_int(fio, tpx->fep_state);
3371     }
3372     gmx_fio_do_real(fio, tpx->lambda);
3373     gmx_fio_do_int(fio, tpx->bIr);
3374     gmx_fio_do_int(fio, tpx->bTop);
3375     gmx_fio_do_int(fio, tpx->bX);
3376     gmx_fio_do_int(fio, tpx->bV);
3377     gmx_fio_do_int(fio, tpx->bF);
3378     gmx_fio_do_int(fio, tpx->bBox);
3379
3380     if ((fgen > tpx_generation))
3381     {
3382         /* This can only happen if TopOnlyOK=TRUE */
3383         tpx->bIr = FALSE;
3384     }
3385 }
3386
3387 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3388                   t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3389                   gmx_bool bXVallocated)
3390 {
3391     t_tpxheader     tpx;
3392     t_inputrec      dum_ir;
3393     gmx_mtop_t      dum_top;
3394     gmx_bool        TopOnlyOK;
3395     int             file_version, file_generation;
3396     int             i;
3397     rvec           *xptr, *vptr;
3398     int             ePBC;
3399     gmx_bool        bPeriodicMols;
3400
3401     if (!bRead)
3402     {
3403         tpx.natoms    = state->natoms;
3404         tpx.ngtc      = state->ngtc; /* need to add nnhpres here? */
3405         tpx.fep_state = state->fep_state;
3406         tpx.lambda    = state->lambda[efptFEP];
3407         tpx.bIr       = (ir       != NULL);
3408         tpx.bTop      = (mtop     != NULL);
3409         tpx.bX        = (state->x != NULL);
3410         tpx.bV        = (state->v != NULL);
3411         tpx.bF        = (f        != NULL);
3412         tpx.bBox      = TRUE;
3413     }
3414
3415     TopOnlyOK = (ir == NULL);
3416
3417     do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3418
3419     if (bRead)
3420     {
3421         state->flags  = 0;
3422         /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3423         /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3424         if (bXVallocated)
3425         {
3426             xptr = state->x;
3427             vptr = state->v;
3428             init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3429             state->natoms = tpx.natoms;
3430             state->nalloc = tpx.natoms;
3431             state->x      = xptr;
3432             state->v      = vptr;
3433         }
3434         else
3435         {
3436             init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3437         }
3438     }
3439
3440 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3441
3442     do_test(fio, tpx.bBox, state->box);
3443     if (tpx.bBox)
3444     {
3445         gmx_fio_ndo_rvec(fio, state->box, DIM);
3446         if (file_version >= 51)
3447         {
3448             gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3449         }
3450         else
3451         {
3452             /* We initialize box_rel after reading the inputrec */
3453             clear_mat(state->box_rel);
3454         }
3455         if (file_version >= 28)
3456         {
3457             gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3458             if (file_version < 56)
3459             {
3460                 matrix mdum;
3461                 gmx_fio_ndo_rvec(fio, mdum, DIM);
3462             }
3463         }
3464     }
3465
3466     if (state->ngtc > 0 && file_version >= 28)
3467     {
3468         real *dumv;
3469         /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3470         /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3471         /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3472         snew(dumv, state->ngtc);
3473         if (file_version < 69)
3474         {
3475             gmx_fio_ndo_real(fio, dumv, state->ngtc);
3476         }
3477         /* These used to be the Berendsen tcoupl_lambda's */
3478         gmx_fio_ndo_real(fio, dumv, state->ngtc);
3479         sfree(dumv);
3480     }
3481
3482     /* Prior to tpx version 26, the inputrec was here.
3483      * I moved it to enable partial forward-compatibility
3484      * for analysis/viewer programs.
3485      */
3486     if (file_version < 26)
3487     {
3488         do_test(fio, tpx.bIr, ir);
3489         if (tpx.bIr)
3490         {
3491             if (ir)
3492             {
3493                 do_inputrec(fio, ir, bRead, file_version,
3494                             mtop ? &mtop->ffparams.fudgeQQ : NULL);
3495                 if (bRead && debug)
3496                 {
3497                     pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3498                 }
3499             }
3500             else
3501             {
3502                 do_inputrec(fio, &dum_ir, bRead, file_version,
3503                             mtop ? &mtop->ffparams.fudgeQQ : NULL);
3504                 if (bRead && debug)
3505                 {
3506                     pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3507                 }
3508                 done_inputrec(&dum_ir);
3509             }
3510
3511         }
3512     }
3513
3514     do_test(fio, tpx.bTop, mtop);
3515     if (tpx.bTop)
3516     {
3517         if (mtop)
3518         {
3519             do_mtop(fio, mtop, bRead, file_version);
3520         }
3521         else
3522         {
3523             do_mtop(fio, &dum_top, bRead, file_version);
3524             done_mtop(&dum_top, TRUE);
3525         }
3526     }
3527     do_test(fio, tpx.bX, state->x);
3528     if (tpx.bX)
3529     {
3530         if (bRead)
3531         {
3532             state->flags |= (1<<estX);
3533         }
3534         gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3535     }
3536
3537     do_test(fio, tpx.bV, state->v);
3538     if (tpx.bV)
3539     {
3540         if (bRead)
3541         {
3542             state->flags |= (1<<estV);
3543         }
3544         gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3545     }
3546
3547     do_test(fio, tpx.bF, f);
3548     if (tpx.bF)
3549     {
3550         gmx_fio_ndo_rvec(fio, f, state->natoms);
3551     }
3552
3553     /* Starting with tpx version 26, we have the inputrec
3554      * at the end of the file, so we can ignore it
3555      * if the file is never than the software (but still the
3556      * same generation - see comments at the top of this file.
3557      *
3558      *
3559      */
3560     ePBC          = -1;
3561     bPeriodicMols = FALSE;
3562     if (file_version >= 26)
3563     {
3564         do_test(fio, tpx.bIr, ir);
3565         if (tpx.bIr)
3566         {
3567             if (file_version >= 53)
3568             {
3569                 /* Removed the pbc info from do_inputrec, since we always want it */
3570                 if (!bRead)
3571                 {
3572                     ePBC          = ir->ePBC;
3573                     bPeriodicMols = ir->bPeriodicMols;
3574                 }
3575                 gmx_fio_do_int(fio, ePBC);
3576                 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3577             }
3578             if (file_generation <= tpx_generation && ir)
3579             {
3580                 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3581                 if (bRead && debug)
3582                 {
3583                     pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3584                 }
3585                 if (file_version < 51)
3586                 {
3587                     set_box_rel(ir, state);
3588                 }
3589                 if (file_version < 53)
3590                 {
3591                     ePBC          = ir->ePBC;
3592                     bPeriodicMols = ir->bPeriodicMols;
3593                 }
3594             }
3595             if (bRead && ir && file_version >= 53)
3596             {
3597                 /* We need to do this after do_inputrec, since that initializes ir */
3598                 ir->ePBC          = ePBC;
3599                 ir->bPeriodicMols = bPeriodicMols;
3600             }
3601         }
3602     }
3603
3604     if (bRead)
3605     {
3606         if (tpx.bIr && ir)
3607         {
3608             if (state->ngtc == 0)
3609             {
3610                 /* Reading old version without tcoupl state data: set it */
3611                 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3612             }
3613             if (tpx.bTop && mtop)
3614             {
3615                 if (file_version < 57)
3616                 {
3617                     if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3618                     {
3619                         ir->eDisre = edrSimple;
3620                     }
3621                     else
3622                     {
3623                         ir->eDisre = edrNone;
3624                     }
3625                 }
3626                 set_disres_npair(mtop);
3627             }
3628         }
3629
3630         if (tpx.bTop && mtop)
3631         {
3632             gmx_mtop_finalize(mtop);
3633         }
3634
3635         if (file_version >= 57)
3636         {
3637             char *env;
3638             int   ienv;
3639             env = getenv("GMX_NOCHARGEGROUPS");
3640             if (env != NULL)
3641             {
3642                 sscanf(env, "%d", &ienv);
3643                 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3644                         ienv);
3645                 if (ienv > 0)
3646                 {
3647                     fprintf(stderr,
3648                             "Will make single atomic charge groups in non-solvent%s\n",
3649                             ienv > 1 ? " and solvent" : "");
3650                     gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3651                 }
3652                 fprintf(stderr, "\n");
3653             }
3654         }
3655     }
3656
3657     return ePBC;
3658 }
3659
3660 /************************************************************
3661  *
3662  *  The following routines are the exported ones
3663  *
3664  ************************************************************/
3665
3666 t_fileio *open_tpx(const char *fn, const char *mode)
3667 {
3668     return gmx_fio_open(fn, mode);
3669 }
3670
3671 void close_tpx(t_fileio *fio)
3672 {
3673     gmx_fio_close(fio);
3674 }
3675
3676 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3677                     int *file_version, int *file_generation)
3678 {
3679     t_fileio *fio;
3680
3681     fio = open_tpx(fn, "r");
3682     do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3683     close_tpx(fio);
3684 }
3685
3686 void write_tpx_state(const char *fn,
3687                      t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3688 {
3689     t_fileio *fio;
3690
3691     fio = open_tpx(fn, "w");
3692     do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3693     close_tpx(fio);
3694 }
3695
3696 void read_tpx_state(const char *fn,
3697                     t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3698 {
3699     t_fileio *fio;
3700
3701     fio = open_tpx(fn, "r");
3702     do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3703     close_tpx(fio);
3704 }
3705
3706 int read_tpx(const char *fn,
3707              t_inputrec *ir, matrix box, int *natoms,
3708              rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3709 {
3710     t_fileio *fio;
3711     t_state   state;
3712     int       ePBC;
3713
3714     state.x = x;
3715     state.v = v;
3716     fio     = open_tpx(fn, "r");
3717     ePBC    = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3718     close_tpx(fio);
3719     *natoms = state.natoms;
3720     if (box)
3721     {
3722         copy_mat(state.box, box);
3723     }
3724     state.x = NULL;
3725     state.v = NULL;
3726     done_state(&state);
3727
3728     return ePBC;
3729 }
3730
3731 int read_tpx_top(const char *fn,
3732                  t_inputrec *ir, matrix box, int *natoms,
3733                  rvec *x, rvec *v, rvec *f, t_topology *top)
3734 {
3735     gmx_mtop_t  mtop;
3736     t_topology *ltop;
3737     int         ePBC;
3738
3739     ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3740
3741     *top = gmx_mtop_t_to_t_topology(&mtop);
3742
3743     return ePBC;
3744 }
3745
3746 gmx_bool fn2bTPX(const char *file)
3747 {
3748     return (efTPR == fn2ftp(file));
3749 }
3750
3751 static void done_gmx_groups_t(gmx_groups_t *g)
3752 {
3753     int i;
3754
3755     for (i = 0; (i < egcNR); i++)
3756     {
3757         if (NULL != g->grps[i].nm_ind)
3758         {
3759             sfree(g->grps[i].nm_ind);
3760             g->grps[i].nm_ind = NULL;
3761         }
3762         if (NULL != g->grpnr[i])
3763         {
3764             sfree(g->grpnr[i]);
3765             g->grpnr[i] = NULL;
3766         }
3767     }
3768     /* The contents of this array is in symtab, don't free it here */
3769     sfree(g->grpname);
3770 }
3771
3772 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3773                        rvec **x, rvec **v, matrix box, gmx_bool bMass)
3774 {
3775     t_tpxheader      header;
3776     int              natoms, i, version, generation;
3777     gmx_bool         bTop, bXNULL = FALSE;
3778     gmx_mtop_t      *mtop;
3779     t_topology      *topconv;
3780     gmx_atomprop_t   aps;
3781
3782     bTop  = fn2bTPX(infile);
3783     *ePBC = -1;
3784     if (bTop)
3785     {
3786         read_tpxheader(infile, &header, TRUE, &version, &generation);
3787         if (x)
3788         {
3789             snew(*x, header.natoms);
3790         }
3791         if (v)
3792         {
3793             snew(*v, header.natoms);
3794         }
3795         snew(mtop, 1);
3796         *ePBC = read_tpx(infile, NULL, box, &natoms,
3797                          (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3798         *top = gmx_mtop_t_to_t_topology(mtop);
3799         /* In this case we need to throw away the group data too */
3800         done_gmx_groups_t(&mtop->groups);
3801         sfree(mtop);
3802         strcpy(title, *top->name);
3803         tpx_make_chain_identifiers(&top->atoms, &top->mols);
3804     }
3805     else
3806     {
3807         get_stx_coordnum(infile, &natoms);
3808         init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3809         if (x == NULL)
3810         {
3811             snew(x, 1);
3812             bXNULL = TRUE;
3813         }
3814         snew(*x, natoms);
3815         if (v)
3816         {
3817             snew(*v, natoms);
3818         }
3819         read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3820         if (bXNULL)
3821         {
3822             sfree(*x);
3823             sfree(x);
3824         }
3825         if (bMass)
3826         {
3827             aps = gmx_atomprop_init();
3828             for (i = 0; (i < natoms); i++)
3829             {
3830                 if (!gmx_atomprop_query(aps, epropMass,
3831                                         *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3832                                         *top->atoms.atomname[i],
3833                                         &(top->atoms.atom[i].m)))
3834                 {
3835                     if (debug)
3836                     {
3837                         fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3838                                 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3839                                 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3840                                 *top->atoms.atomname[i]);
3841                     }
3842                 }
3843             }
3844             gmx_atomprop_destroy(aps);
3845         }
3846         top->idef.ntypes = -1;
3847     }
3848
3849     return bTop;
3850 }