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