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