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