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