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