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