Fixes #1194 error reading 4.5.5 tpr files with dihres
[alexxy/gromacs.git] / src / gmxlib / tpxio.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team,
6  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012,2013, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 /* This file is completely threadsafe - keep it that way! */
43 #ifdef GMX_THREAD_MPI
44 #include <thread_mpi.h>
45 #endif
46
47
48 #include <ctype.h>
49 #include "sysstuff.h"
50 #include "smalloc.h"
51 #include "string2.h"
52 #include "gmx_fatal.h"
53 #include "macros.h"
54 #include "names.h"
55 #include "symtab.h"
56 #include "futil.h"
57 #include "filenm.h"
58 #include "gmxfio.h"
59 #include "topsort.h"
60 #include "tpxio.h"
61 #include "txtdump.h"
62 #include "confio.h"
63 #include "atomprop.h"
64 #include "copyrite.h"
65 #include "vec.h"
66 #include "mtop_util.h"
67
68 #define TPX_TAG_RELEASE  "release"
69
70 /* This is the tag string which is stored in the tpx file.
71  * Change this if you want to change the tpx format in a feature branch.
72  * This ensures that there will not be different tpx formats around which
73  * can not be distinguished.
74  */
75 static const char *tpx_tag = TPX_TAG_RELEASE;
76
77 /* This number should be increased whenever the file format changes! */
78 static const int tpx_version = 83;
79
80 /* This number should only be increased when you edit the TOPOLOGY section
81  * or the HEADER of the tpx format.
82  * This way we can maintain forward compatibility too for all analysis tools
83  * and/or external programs that only need to know the atom/residue names,
84  * charges, and bond connectivity.
85  *
86  * It first appeared in tpx version 26, when I also moved the inputrecord
87  * to the end of the tpx file, so we can just skip it if we only
88  * want the topology.
89  */
90 static const int tpx_generation = 24;
91
92 /* This number should be the most recent backwards incompatible version
93  * I.e., if this number is 9, we cannot read tpx version 9 with this code.
94  */
95 static const int tpx_incompatible_version = 9;
96
97
98
99 /* Struct used to maintain tpx compatibility when function types are added */
100 typedef struct {
101     int fvnr;  /* file version number in which the function type first appeared */
102     int ftype; /* function type */
103 } t_ftupd;
104
105 /*
106  * The entries should be ordered in:
107  * 1. ascending file version number
108  * 2. ascending function type number
109  */
110 /*static const t_ftupd ftupd[] = {
111    { 20, F_CUBICBONDS        },
112    { 20, F_CONNBONDS         },
113    { 20, F_HARMONIC          },
114    { 20, F_EQM,              },
115    { 22, F_DISRESVIOL        },
116    { 22, F_ORIRES            },
117    { 22, F_ORIRESDEV         },
118    { 26, F_FOURDIHS          },
119    { 26, F_PIDIHS            },
120    { 26, F_DIHRES            },
121    { 26, F_DIHRESVIOL        },
122    { 30, F_CROSS_BOND_BONDS  },
123    { 30, F_CROSS_BOND_ANGLES },
124    { 30, F_UREY_BRADLEY      },
125    { 30, F_POLARIZATION      },
126    { 54, F_DHDL_CON          },
127    };*/
128 /*
129  * The entries should be ordered in:
130  * 1. ascending function type number
131  * 2. ascending file version number
132  */
133 /* question; what is the purpose of the commented code above? */
134 static const t_ftupd ftupd[] = {
135     { 20, F_CUBICBONDS        },
136     { 20, F_CONNBONDS         },
137     { 20, F_HARMONIC          },
138     { 34, F_FENEBONDS         },
139     { 43, F_TABBONDS          },
140     { 43, F_TABBONDSNC        },
141     { 70, F_RESTRBONDS        },
142     { 76, F_LINEAR_ANGLES     },
143     { 30, F_CROSS_BOND_BONDS  },
144     { 30, F_CROSS_BOND_ANGLES },
145     { 30, F_UREY_BRADLEY      },
146     { 34, F_QUARTIC_ANGLES    },
147     { 43, F_TABANGLES         },
148     { 26, F_FOURDIHS          },
149     { 26, F_PIDIHS            },
150     { 43, F_TABDIHS           },
151     { 65, F_CMAP              },
152     { 60, F_GB12              },
153     { 61, F_GB13              },
154     { 61, F_GB14              },
155     { 72, F_GBPOL             },
156     { 72, F_NPSOLVATION       },
157     { 41, F_LJC14_Q           },
158     { 41, F_LJC_PAIRS_NB      },
159     { 32, F_BHAM_LR           },
160     { 32, F_RF_EXCL           },
161     { 32, F_COUL_RECIP        },
162     { 46, F_DPD               },
163     { 30, F_POLARIZATION      },
164     { 36, F_THOLE_POL         },
165     { 22, F_DISRESVIOL        },
166     { 22, F_ORIRES            },
167     { 22, F_ORIRESDEV         },
168     { 26, F_DIHRES            },
169     { 26, F_DIHRESVIOL        },
170     { 49, F_VSITE4FDN         },
171     { 50, F_VSITEN            },
172     { 46, F_COM_PULL          },
173     { 20, F_EQM               },
174     { 46, F_ECONSERVED        },
175     { 69, F_VTEMP_NOLONGERUSED},
176     { 66, F_PDISPCORR         },
177     { 54, F_DVDL_CONSTR       },
178     { 76, F_ANHARM_POL        },
179     { 79, F_DVDL_COUL         },
180     { 79, F_DVDL_VDW,         },
181     { 79, F_DVDL_BONDED,      },
182     { 79, F_DVDL_RESTRAINT    },
183     { 79, F_DVDL_TEMPERATURE  },
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)
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 >= 80)
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 >= 80)
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)
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 >= 80)
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 < 82)
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 >= 82)
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_RBDIHS:
1826             bDum = gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1827             if (file_version >= 25)
1828             {
1829                 bDum = gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1830             }
1831             break;
1832         case F_FOURDIHS:
1833             /* Fourier dihedrals are internally represented
1834              * as Ryckaert-Bellemans since those are faster to compute.
1835              */
1836             bDum = gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1837             bDum = gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1838             break;
1839         case F_CONSTR:
1840         case F_CONSTRNC:
1841             gmx_fio_do_real(fio, iparams->constr.dA);
1842             gmx_fio_do_real(fio, iparams->constr.dB);
1843             break;
1844         case F_SETTLE:
1845             gmx_fio_do_real(fio, iparams->settle.doh);
1846             gmx_fio_do_real(fio, iparams->settle.dhh);
1847             break;
1848         case F_VSITE2:
1849             gmx_fio_do_real(fio, iparams->vsite.a);
1850             break;
1851         case F_VSITE3:
1852         case F_VSITE3FD:
1853         case F_VSITE3FAD:
1854             gmx_fio_do_real(fio, iparams->vsite.a);
1855             gmx_fio_do_real(fio, iparams->vsite.b);
1856             break;
1857         case F_VSITE3OUT:
1858         case F_VSITE4FD:
1859         case F_VSITE4FDN:
1860             gmx_fio_do_real(fio, iparams->vsite.a);
1861             gmx_fio_do_real(fio, iparams->vsite.b);
1862             gmx_fio_do_real(fio, iparams->vsite.c);
1863             break;
1864         case F_VSITEN:
1865             gmx_fio_do_int(fio, iparams->vsiten.n);
1866             gmx_fio_do_real(fio, iparams->vsiten.a);
1867             break;
1868         case F_GB12:
1869         case F_GB13:
1870         case F_GB14:
1871             /* We got rid of some parameters in version 68 */
1872             if (bRead && file_version < 68)
1873             {
1874                 gmx_fio_do_real(fio, rdum);
1875                 gmx_fio_do_real(fio, rdum);
1876                 gmx_fio_do_real(fio, rdum);
1877                 gmx_fio_do_real(fio, rdum);
1878             }
1879             gmx_fio_do_real(fio, iparams->gb.sar);
1880             gmx_fio_do_real(fio, iparams->gb.st);
1881             gmx_fio_do_real(fio, iparams->gb.pi);
1882             gmx_fio_do_real(fio, iparams->gb.gbr);
1883             gmx_fio_do_real(fio, iparams->gb.bmlt);
1884             break;
1885         case F_CMAP:
1886             gmx_fio_do_int(fio, iparams->cmap.cmapA);
1887             gmx_fio_do_int(fio, iparams->cmap.cmapB);
1888             break;
1889         default:
1890             gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1891                       ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1892     }
1893     if (!bRead)
1894     {
1895         gmx_fio_unset_comment(fio);
1896     }
1897 }
1898
1899 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
1900                      int ftype)
1901 {
1902     int      i, k, idum;
1903     gmx_bool bDum = TRUE;
1904
1905     if (!bRead)
1906     {
1907         gmx_fio_set_comment(fio, interaction_function[ftype].name);
1908     }
1909     if (file_version < 44)
1910     {
1911         for (i = 0; i < MAXNODES; i++)
1912         {
1913             gmx_fio_do_int(fio, idum);
1914         }
1915     }
1916     gmx_fio_do_int(fio, ilist->nr);
1917     if (bRead)
1918     {
1919         snew(ilist->iatoms, ilist->nr);
1920     }
1921     bDum = gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
1922     if (!bRead)
1923     {
1924         gmx_fio_unset_comment(fio);
1925     }
1926 }
1927
1928 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1929                         gmx_bool bRead, int file_version)
1930 {
1931     int          idum, i, j;
1932     gmx_bool     bDum = TRUE;
1933     unsigned int k;
1934
1935     gmx_fio_do_int(fio, ffparams->atnr);
1936     if (file_version < 57)
1937     {
1938         gmx_fio_do_int(fio, idum);
1939     }
1940     gmx_fio_do_int(fio, ffparams->ntypes);
1941     if (bRead && debug)
1942     {
1943         fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
1944                 ffparams->atnr, ffparams->ntypes);
1945     }
1946     if (bRead)
1947     {
1948         snew(ffparams->functype, ffparams->ntypes);
1949         snew(ffparams->iparams, ffparams->ntypes);
1950     }
1951     /* Read/write all the function types */
1952     bDum = gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
1953     if (bRead && debug)
1954     {
1955         pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
1956     }
1957
1958     if (file_version >= 66)
1959     {
1960         gmx_fio_do_double(fio, ffparams->reppow);
1961     }
1962     else
1963     {
1964         ffparams->reppow = 12.0;
1965     }
1966
1967     if (file_version >= 57)
1968     {
1969         gmx_fio_do_real(fio, ffparams->fudgeQQ);
1970     }
1971
1972     /* Check whether all these function types are supported by the code.
1973      * In practice the code is backwards compatible, which means that the
1974      * numbering may have to be altered from old numbering to new numbering
1975      */
1976     for (i = 0; (i < ffparams->ntypes); i++)
1977     {
1978         if (bRead)
1979         {
1980             /* Loop over file versions */
1981             for (k = 0; (k < NFTUPD); k++)
1982             {
1983                 /* Compare the read file_version to the update table */
1984                 if ((file_version < ftupd[k].fvnr) &&
1985                     (ffparams->functype[i] >= ftupd[k].ftype))
1986                 {
1987                     ffparams->functype[i] += 1;
1988                     if (debug)
1989                     {
1990                         fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
1991                                 i, ffparams->functype[i],
1992                                 interaction_function[ftupd[k].ftype].longname);
1993                         fflush(debug);
1994                     }
1995                 }
1996             }
1997         }
1998
1999         do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2000                    file_version);
2001         if (bRead && debug)
2002         {
2003             pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2004         }
2005     }
2006 }
2007
2008 static void add_settle_atoms(t_ilist *ilist)
2009 {
2010     int i;
2011
2012     /* Settle used to only store the first atom: add the other two */
2013     srenew(ilist->iatoms, 2*ilist->nr);
2014     for (i = ilist->nr/2-1; i >= 0; i--)
2015     {
2016         ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2017         ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2018         ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2019         ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2020     }
2021     ilist->nr = 2*ilist->nr;
2022 }
2023
2024 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2025                       int file_version)
2026 {
2027     int          i, j, renum[F_NRE];
2028     gmx_bool     bDum = TRUE, bClear;
2029     unsigned int k;
2030
2031     for (j = 0; (j < F_NRE); j++)
2032     {
2033         bClear = FALSE;
2034         if (bRead)
2035         {
2036             for (k = 0; k < NFTUPD; k++)
2037             {
2038                 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2039                 {
2040                     bClear = TRUE;
2041                 }
2042             }
2043         }
2044         if (bClear)
2045         {
2046             ilist[j].nr     = 0;
2047             ilist[j].iatoms = NULL;
2048         }
2049         else
2050         {
2051             do_ilist(fio, &ilist[j], bRead, file_version, j);
2052             if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2053             {
2054                 add_settle_atoms(&ilist[j]);
2055             }
2056         }
2057         /*
2058            if (bRead && gmx_debug_at)
2059            pr_ilist(debug,0,interaction_function[j].longname,
2060                functype,&ilist[j],TRUE);
2061          */
2062     }
2063 }
2064
2065 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2066                     gmx_bool bRead, int file_version)
2067 {
2068     do_ffparams(fio, ffparams, bRead, file_version);
2069
2070     if (file_version >= 54)
2071     {
2072         gmx_fio_do_real(fio, ffparams->fudgeQQ);
2073     }
2074
2075     do_ilists(fio, molt->ilist, bRead, file_version);
2076 }
2077
2078 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2079 {
2080     int      i, idum, dum_nra, *dum_a;
2081     gmx_bool bDum = TRUE;
2082
2083     if (file_version < 44)
2084     {
2085         for (i = 0; i < MAXNODES; i++)
2086         {
2087             gmx_fio_do_int(fio, idum);
2088         }
2089     }
2090     gmx_fio_do_int(fio, block->nr);
2091     if (file_version < 51)
2092     {
2093         gmx_fio_do_int(fio, dum_nra);
2094     }
2095     if (bRead)
2096     {
2097         block->nalloc_index = block->nr+1;
2098         snew(block->index, block->nalloc_index);
2099     }
2100     bDum = gmx_fio_ndo_int(fio, block->index, block->nr+1);
2101
2102     if (file_version < 51 && dum_nra > 0)
2103     {
2104         snew(dum_a, dum_nra);
2105         bDum = gmx_fio_ndo_int(fio, dum_a, dum_nra);
2106         sfree(dum_a);
2107     }
2108 }
2109
2110 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2111                       int file_version)
2112 {
2113     int      i, idum;
2114     gmx_bool bDum = TRUE;
2115
2116     if (file_version < 44)
2117     {
2118         for (i = 0; i < MAXNODES; i++)
2119         {
2120             gmx_fio_do_int(fio, idum);
2121         }
2122     }
2123     gmx_fio_do_int(fio, block->nr);
2124     gmx_fio_do_int(fio, block->nra);
2125     if (bRead)
2126     {
2127         block->nalloc_index = block->nr+1;
2128         snew(block->index, block->nalloc_index);
2129         block->nalloc_a = block->nra;
2130         snew(block->a, block->nalloc_a);
2131     }
2132     bDum = gmx_fio_ndo_int(fio, block->index, block->nr+1);
2133     bDum = gmx_fio_ndo_int(fio, block->a, block->nra);
2134 }
2135
2136 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2137                     int file_version, gmx_groups_t *groups, int atnr)
2138 {
2139     int i, myngrp;
2140
2141     gmx_fio_do_real(fio, atom->m);
2142     gmx_fio_do_real(fio, atom->q);
2143     gmx_fio_do_real(fio, atom->mB);
2144     gmx_fio_do_real(fio, atom->qB);
2145     gmx_fio_do_ushort(fio, atom->type);
2146     gmx_fio_do_ushort(fio, atom->typeB);
2147     gmx_fio_do_int(fio, atom->ptype);
2148     gmx_fio_do_int(fio, atom->resind);
2149     if (file_version >= 52)
2150     {
2151         gmx_fio_do_int(fio, atom->atomnumber);
2152     }
2153     else if (bRead)
2154     {
2155         atom->atomnumber = NOTSET;
2156     }
2157     if (file_version < 23)
2158     {
2159         myngrp = 8;
2160     }
2161     else if (file_version < 39)
2162     {
2163         myngrp = 9;
2164     }
2165     else
2166     {
2167         myngrp = ngrp;
2168     }
2169
2170     if (file_version < 57)
2171     {
2172         unsigned char uchar[egcNR];
2173         gmx_fio_ndo_uchar(fio, uchar, myngrp);
2174         for (i = myngrp; (i < ngrp); i++)
2175         {
2176             uchar[i] = 0;
2177         }
2178         /* Copy the old data format to the groups struct */
2179         for (i = 0; i < ngrp; i++)
2180         {
2181             groups->grpnr[i][atnr] = uchar[i];
2182         }
2183     }
2184 }
2185
2186 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2187                     int file_version)
2188 {
2189     int      i, j, myngrp;
2190     gmx_bool bDum = TRUE;
2191
2192     if (file_version < 23)
2193     {
2194         myngrp = 8;
2195     }
2196     else if (file_version < 39)
2197     {
2198         myngrp = 9;
2199     }
2200     else
2201     {
2202         myngrp = ngrp;
2203     }
2204
2205     for (j = 0; (j < ngrp); j++)
2206     {
2207         if (j < myngrp)
2208         {
2209             gmx_fio_do_int(fio, grps[j].nr);
2210             if (bRead)
2211             {
2212                 snew(grps[j].nm_ind, grps[j].nr);
2213             }
2214             bDum = gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2215         }
2216         else
2217         {
2218             grps[j].nr = 1;
2219             snew(grps[j].nm_ind, grps[j].nr);
2220         }
2221     }
2222 }
2223
2224 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2225 {
2226     int ls;
2227
2228     if (bRead)
2229     {
2230         gmx_fio_do_int(fio, ls);
2231         *nm = get_symtab_handle(symtab, ls);
2232     }
2233     else
2234     {
2235         ls = lookup_symtab(symtab, *nm);
2236         gmx_fio_do_int(fio, ls);
2237     }
2238 }
2239
2240 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2241                       t_symtab *symtab)
2242 {
2243     int  j;
2244
2245     for (j = 0; (j < nstr); j++)
2246     {
2247         do_symstr(fio, &(nm[j]), bRead, symtab);
2248     }
2249 }
2250
2251 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2252                        t_symtab *symtab, int file_version)
2253 {
2254     int  j;
2255
2256     for (j = 0; (j < n); j++)
2257     {
2258         do_symstr(fio, &(ri[j].name), bRead, symtab);
2259         if (file_version >= 63)
2260         {
2261             gmx_fio_do_int(fio, ri[j].nr);
2262             gmx_fio_do_uchar(fio, ri[j].ic);
2263         }
2264         else
2265         {
2266             ri[j].nr = j + 1;
2267             ri[j].ic = ' ';
2268         }
2269     }
2270 }
2271
2272 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2273                      int file_version,
2274                      gmx_groups_t *groups)
2275 {
2276     int i;
2277
2278     gmx_fio_do_int(fio, atoms->nr);
2279     gmx_fio_do_int(fio, atoms->nres);
2280     if (file_version < 57)
2281     {
2282         gmx_fio_do_int(fio, groups->ngrpname);
2283         for (i = 0; i < egcNR; i++)
2284         {
2285             groups->ngrpnr[i] = atoms->nr;
2286             snew(groups->grpnr[i], groups->ngrpnr[i]);
2287         }
2288     }
2289     if (bRead)
2290     {
2291         snew(atoms->atom, atoms->nr);
2292         snew(atoms->atomname, atoms->nr);
2293         snew(atoms->atomtype, atoms->nr);
2294         snew(atoms->atomtypeB, atoms->nr);
2295         snew(atoms->resinfo, atoms->nres);
2296         if (file_version < 57)
2297         {
2298             snew(groups->grpname, groups->ngrpname);
2299         }
2300         atoms->pdbinfo = NULL;
2301     }
2302     for (i = 0; (i < atoms->nr); i++)
2303     {
2304         do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2305     }
2306     do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2307     if (bRead && (file_version <= 20))
2308     {
2309         for (i = 0; i < atoms->nr; i++)
2310         {
2311             atoms->atomtype[i]  = put_symtab(symtab, "?");
2312             atoms->atomtypeB[i] = put_symtab(symtab, "?");
2313         }
2314     }
2315     else
2316     {
2317         do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2318         do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2319     }
2320     do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2321
2322     if (file_version < 57)
2323     {
2324         do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2325
2326         do_grps(fio, egcNR, groups->grps, bRead, file_version);
2327     }
2328 }
2329
2330 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2331                       gmx_bool bRead, t_symtab *symtab,
2332                       int file_version)
2333 {
2334     int      g, n, i;
2335     gmx_bool bDum = TRUE;
2336
2337     do_grps(fio, egcNR, groups->grps, bRead, file_version);
2338     gmx_fio_do_int(fio, groups->ngrpname);
2339     if (bRead)
2340     {
2341         snew(groups->grpname, groups->ngrpname);
2342     }
2343     do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2344     for (g = 0; g < egcNR; g++)
2345     {
2346         gmx_fio_do_int(fio, groups->ngrpnr[g]);
2347         if (groups->ngrpnr[g] == 0)
2348         {
2349             if (bRead)
2350             {
2351                 groups->grpnr[g] = NULL;
2352             }
2353         }
2354         else
2355         {
2356             if (bRead)
2357             {
2358                 snew(groups->grpnr[g], groups->ngrpnr[g]);
2359             }
2360             bDum = gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2361         }
2362     }
2363 }
2364
2365 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2366                          t_symtab *symtab, int file_version)
2367 {
2368     int      i, j;
2369     gmx_bool bDum = TRUE;
2370
2371     if (file_version > 25)
2372     {
2373         gmx_fio_do_int(fio, atomtypes->nr);
2374         j = atomtypes->nr;
2375         if (bRead)
2376         {
2377             snew(atomtypes->radius, j);
2378             snew(atomtypes->vol, j);
2379             snew(atomtypes->surftens, j);
2380             snew(atomtypes->atomnumber, j);
2381             snew(atomtypes->gb_radius, j);
2382             snew(atomtypes->S_hct, j);
2383         }
2384         bDum = gmx_fio_ndo_real(fio, atomtypes->radius, j);
2385         bDum = gmx_fio_ndo_real(fio, atomtypes->vol, j);
2386         bDum = gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2387         if (file_version >= 40)
2388         {
2389             bDum = gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2390         }
2391         if (file_version >= 60)
2392         {
2393             bDum = gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2394             bDum = gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2395         }
2396     }
2397     else
2398     {
2399         /* File versions prior to 26 cannot do GBSA,
2400          * so they dont use this structure
2401          */
2402         atomtypes->nr         = 0;
2403         atomtypes->radius     = NULL;
2404         atomtypes->vol        = NULL;
2405         atomtypes->surftens   = NULL;
2406         atomtypes->atomnumber = NULL;
2407         atomtypes->gb_radius  = NULL;
2408         atomtypes->S_hct      = NULL;
2409     }
2410 }
2411
2412 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2413 {
2414     int       i, nr;
2415     t_symbuf *symbuf;
2416     char      buf[STRLEN];
2417
2418     gmx_fio_do_int(fio, symtab->nr);
2419     nr     = symtab->nr;
2420     if (bRead)
2421     {
2422         snew(symtab->symbuf, 1);
2423         symbuf          = symtab->symbuf;
2424         symbuf->bufsize = nr;
2425         snew(symbuf->buf, nr);
2426         for (i = 0; (i < nr); i++)
2427         {
2428             gmx_fio_do_string(fio, buf);
2429             symbuf->buf[i] = strdup(buf);
2430         }
2431     }
2432     else
2433     {
2434         symbuf = symtab->symbuf;
2435         while (symbuf != NULL)
2436         {
2437             for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2438             {
2439                 gmx_fio_do_string(fio, symbuf->buf[i]);
2440             }
2441             nr    -= i;
2442             symbuf = symbuf->next;
2443         }
2444         if (nr != 0)
2445         {
2446             gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2447         }
2448     }
2449 }
2450
2451 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2452 {
2453     int i, j, ngrid, gs, nelem;
2454
2455     gmx_fio_do_int(fio, cmap_grid->ngrid);
2456     gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2457
2458     ngrid = cmap_grid->ngrid;
2459     gs    = cmap_grid->grid_spacing;
2460     nelem = gs * gs;
2461
2462     if (bRead)
2463     {
2464         snew(cmap_grid->cmapdata, ngrid);
2465
2466         for (i = 0; i < cmap_grid->ngrid; i++)
2467         {
2468             snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2469         }
2470     }
2471
2472     for (i = 0; i < cmap_grid->ngrid; i++)
2473     {
2474         for (j = 0; j < nelem; j++)
2475         {
2476             gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2477             gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2478             gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2479             gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2480         }
2481     }
2482 }
2483
2484
2485 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2486 {
2487     int  m, a, a0, a1, r;
2488     char c, chainid;
2489     int  chainnum;
2490
2491     /* We always assign a new chain number, but save the chain id characters
2492      * for larger molecules.
2493      */
2494 #define CHAIN_MIN_ATOMS 15
2495
2496     chainnum = 0;
2497     chainid  = 'A';
2498     for (m = 0; m < mols->nr; m++)
2499     {
2500         a0 = mols->index[m];
2501         a1 = mols->index[m+1];
2502         if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2503         {
2504             c = chainid;
2505             chainid++;
2506         }
2507         else
2508         {
2509             c = ' ';
2510         }
2511         for (a = a0; a < a1; a++)
2512         {
2513             atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2514             atoms->resinfo[atoms->atom[a].resind].chainid  = c;
2515         }
2516         chainnum++;
2517     }
2518
2519     /* Blank out the chain id if there was only one chain */
2520     if (chainid == 'B')
2521     {
2522         for (r = 0; r < atoms->nres; r++)
2523         {
2524             atoms->resinfo[r].chainid = ' ';
2525         }
2526     }
2527 }
2528
2529 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2530                        t_symtab *symtab, int file_version,
2531                        gmx_groups_t *groups)
2532 {
2533     int i;
2534
2535     if (file_version >= 57)
2536     {
2537         do_symstr(fio, &(molt->name), bRead, symtab);
2538     }
2539
2540     do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2541
2542     if (bRead && gmx_debug_at)
2543     {
2544         pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2545     }
2546
2547     if (file_version >= 57)
2548     {
2549         do_ilists(fio, molt->ilist, bRead, file_version);
2550
2551         do_block(fio, &molt->cgs, bRead, file_version);
2552         if (bRead && gmx_debug_at)
2553         {
2554             pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2555         }
2556     }
2557
2558     /* This used to be in the atoms struct */
2559     do_blocka(fio, &molt->excls, bRead, file_version);
2560 }
2561
2562 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead,
2563                         int file_version)
2564 {
2565     int i;
2566
2567     gmx_fio_do_int(fio, molb->type);
2568     gmx_fio_do_int(fio, molb->nmol);
2569     gmx_fio_do_int(fio, molb->natoms_mol);
2570     /* Position restraint coordinates */
2571     gmx_fio_do_int(fio, molb->nposres_xA);
2572     if (molb->nposres_xA > 0)
2573     {
2574         if (bRead)
2575         {
2576             snew(molb->posres_xA, molb->nposres_xA);
2577         }
2578         gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2579     }
2580     gmx_fio_do_int(fio, molb->nposres_xB);
2581     if (molb->nposres_xB > 0)
2582     {
2583         if (bRead)
2584         {
2585             snew(molb->posres_xB, molb->nposres_xB);
2586         }
2587         gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2588     }
2589
2590 }
2591
2592 static t_block mtop_mols(gmx_mtop_t *mtop)
2593 {
2594     int     mb, m, a, mol;
2595     t_block mols;
2596
2597     mols.nr = 0;
2598     for (mb = 0; mb < mtop->nmolblock; mb++)
2599     {
2600         mols.nr += mtop->molblock[mb].nmol;
2601     }
2602     mols.nalloc_index = mols.nr + 1;
2603     snew(mols.index, mols.nalloc_index);
2604
2605     a             = 0;
2606     m             = 0;
2607     mols.index[m] = a;
2608     for (mb = 0; mb < mtop->nmolblock; mb++)
2609     {
2610         for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2611         {
2612             a += mtop->molblock[mb].natoms_mol;
2613             m++;
2614             mols.index[m] = a;
2615         }
2616     }
2617
2618     return mols;
2619 }
2620
2621 static void add_posres_molblock(gmx_mtop_t *mtop)
2622 {
2623     t_ilist        *il;
2624     int             am, i, mol, a;
2625     gmx_bool        bFE;
2626     gmx_molblock_t *molb;
2627     t_iparams      *ip;
2628
2629     il = &mtop->moltype[0].ilist[F_POSRES];
2630     if (il->nr == 0)
2631     {
2632         return;
2633     }
2634     am  = 0;
2635     bFE = FALSE;
2636     for (i = 0; i < il->nr; i += 2)
2637     {
2638         ip = &mtop->ffparams.iparams[il->iatoms[i]];
2639         am = max(am, il->iatoms[i+1]);
2640         if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2641             ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2642             ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2643         {
2644             bFE = TRUE;
2645         }
2646     }
2647     /* Make the posres coordinate block end at a molecule end */
2648     mol = 0;
2649     while (am >= mtop->mols.index[mol+1])
2650     {
2651         mol++;
2652     }
2653     molb             = &mtop->molblock[0];
2654     molb->nposres_xA = mtop->mols.index[mol+1];
2655     snew(molb->posres_xA, molb->nposres_xA);
2656     if (bFE)
2657     {
2658         molb->nposres_xB = molb->nposres_xA;
2659         snew(molb->posres_xB, molb->nposres_xB);
2660     }
2661     else
2662     {
2663         molb->nposres_xB = 0;
2664     }
2665     for (i = 0; i < il->nr; i += 2)
2666     {
2667         ip                     = &mtop->ffparams.iparams[il->iatoms[i]];
2668         a                      = il->iatoms[i+1];
2669         molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2670         molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2671         molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2672         if (bFE)
2673         {
2674             molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2675             molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2676             molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2677         }
2678     }
2679 }
2680
2681 static void set_disres_npair(gmx_mtop_t *mtop)
2682 {
2683     int        mt, i, npair;
2684     t_iparams *ip;
2685     t_ilist   *il;
2686     t_iatom   *a;
2687
2688     ip = mtop->ffparams.iparams;
2689
2690     for (mt = 0; mt < mtop->nmoltype; mt++)
2691     {
2692         il = &mtop->moltype[mt].ilist[F_DISRES];
2693         if (il->nr > 0)
2694         {
2695             a     = il->iatoms;
2696             npair = 0;
2697             for (i = 0; i < il->nr; i += 3)
2698             {
2699                 npair++;
2700                 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2701                 {
2702                     ip[a[i]].disres.npair = npair;
2703                     npair                 = 0;
2704                 }
2705             }
2706         }
2707     }
2708 }
2709
2710 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2711                     int file_version)
2712 {
2713     int      mt, mb, i;
2714     t_blocka dumb;
2715
2716     if (bRead)
2717     {
2718         init_mtop(mtop);
2719     }
2720     do_symtab(fio, &(mtop->symtab), bRead);
2721     if (bRead && debug)
2722     {
2723         pr_symtab(debug, 0, "symtab", &mtop->symtab);
2724     }
2725
2726     do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2727
2728     if (file_version >= 57)
2729     {
2730         do_ffparams(fio, &mtop->ffparams, bRead, file_version);
2731
2732         gmx_fio_do_int(fio, mtop->nmoltype);
2733     }
2734     else
2735     {
2736         mtop->nmoltype = 1;
2737     }
2738     if (bRead)
2739     {
2740         snew(mtop->moltype, mtop->nmoltype);
2741         if (file_version < 57)
2742         {
2743             mtop->moltype[0].name = mtop->name;
2744         }
2745     }
2746     for (mt = 0; mt < mtop->nmoltype; mt++)
2747     {
2748         do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
2749                    &mtop->groups);
2750     }
2751
2752     if (file_version >= 57)
2753     {
2754         gmx_fio_do_int(fio, mtop->nmolblock);
2755     }
2756     else
2757     {
2758         mtop->nmolblock = 1;
2759     }
2760     if (bRead)
2761     {
2762         snew(mtop->molblock, mtop->nmolblock);
2763     }
2764     if (file_version >= 57)
2765     {
2766         for (mb = 0; mb < mtop->nmolblock; mb++)
2767         {
2768             do_molblock(fio, &mtop->molblock[mb], bRead, file_version);
2769         }
2770         gmx_fio_do_int(fio, mtop->natoms);
2771     }
2772     else
2773     {
2774         mtop->molblock[0].type       = 0;
2775         mtop->molblock[0].nmol       = 1;
2776         mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
2777         mtop->molblock[0].nposres_xA = 0;
2778         mtop->molblock[0].nposres_xB = 0;
2779     }
2780
2781     do_atomtypes (fio, &(mtop->atomtypes), bRead, &(mtop->symtab), file_version);
2782     if (bRead && debug)
2783     {
2784         pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
2785     }
2786
2787     if (file_version < 57)
2788     {
2789         /* Debug statements are inside do_idef */
2790         do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
2791         mtop->natoms = mtop->moltype[0].atoms.nr;
2792     }
2793
2794     if (file_version >= 65)
2795     {
2796         do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
2797     }
2798     else
2799     {
2800         mtop->ffparams.cmap_grid.ngrid        = 0;
2801         mtop->ffparams.cmap_grid.grid_spacing = 0;
2802         mtop->ffparams.cmap_grid.cmapdata     = NULL;
2803     }
2804
2805     if (file_version >= 57)
2806     {
2807         do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
2808     }
2809
2810     if (file_version < 57)
2811     {
2812         do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
2813         if (bRead && gmx_debug_at)
2814         {
2815             pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
2816         }
2817         do_block(fio, &mtop->mols, bRead, file_version);
2818         /* Add the posres coordinates to the molblock */
2819         add_posres_molblock(mtop);
2820     }
2821     if (bRead)
2822     {
2823         if (file_version >= 57)
2824         {
2825             mtop->mols = mtop_mols(mtop);
2826         }
2827         if (gmx_debug_at)
2828         {
2829             pr_block(debug, 0, "mols", &mtop->mols, TRUE);
2830         }
2831     }
2832
2833     if (file_version < 51)
2834     {
2835         /* Here used to be the shake blocks */
2836         do_blocka(fio, &dumb, bRead, file_version);
2837         if (dumb.nr > 0)
2838         {
2839             sfree(dumb.index);
2840         }
2841         if (dumb.nra > 0)
2842         {
2843             sfree(dumb.a);
2844         }
2845     }
2846
2847     if (bRead)
2848     {
2849         close_symtab(&(mtop->symtab));
2850     }
2851 }
2852
2853 /* If TopOnlyOK is TRUE then we can read even future versions
2854  * of tpx files, provided the file_generation hasn't changed.
2855  * If it is FALSE, we need the inputrecord too, and bail out
2856  * if the file is newer than the program.
2857  *
2858  * The version and generation if the topology (see top of this file)
2859  * are returned in the two last arguments.
2860  *
2861  * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2862  */
2863 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
2864                          gmx_bool TopOnlyOK, int *file_version,
2865                          int *file_generation)
2866 {
2867     char      buf[STRLEN];
2868     char      file_tag[STRLEN];
2869     gmx_bool  bDouble;
2870     int       precision;
2871     int       fver, fgen;
2872     int       idum = 0;
2873     real      rdum = 0;
2874
2875     gmx_fio_checktype(fio);
2876     gmx_fio_setdebug(fio, bDebugMode());
2877
2878     /* NEW! XDR tpb file */
2879     precision = sizeof(real);
2880     if (bRead)
2881     {
2882         gmx_fio_do_string(fio, buf);
2883         if (strncmp(buf, "VERSION", 7))
2884         {
2885             gmx_fatal(FARGS, "Can not read file %s,\n"
2886                       "             this file is from a Gromacs version which is older than 2.0\n"
2887                       "             Make a new one with grompp or use a gro or pdb file, if possible",
2888                       gmx_fio_getname(fio));
2889         }
2890         gmx_fio_do_int(fio, precision);
2891         bDouble = (precision == sizeof(double));
2892         if ((precision != sizeof(float)) && !bDouble)
2893         {
2894             gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
2895                       "instead of %d or %d",
2896                       gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
2897         }
2898         gmx_fio_setprecision(fio, bDouble);
2899         fprintf(stderr, "Reading file %s, %s (%s precision)\n",
2900                 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
2901     }
2902     else
2903     {
2904         gmx_fio_write_string(fio, GromacsVersion());
2905         bDouble = (precision == sizeof(double));
2906         gmx_fio_setprecision(fio, bDouble);
2907         gmx_fio_do_int(fio, precision);
2908         fver = tpx_version;
2909         sprintf(file_tag, "%s", tpx_tag);
2910         fgen = tpx_generation;
2911     }
2912
2913     /* Check versions! */
2914     gmx_fio_do_int(fio, fver);
2915
2916     /* This is for backward compatibility with development versions 77-79
2917      * where the tag was, mistakenly, placed before the generation,
2918      * which would cause a segv instead of a proper error message
2919      * when reading the topology only from tpx with <77 code.
2920      */
2921     if (fver >= 77 && fver <= 79)
2922     {
2923         gmx_fio_do_string(fio, file_tag);
2924     }
2925
2926     if (fver >= 26)
2927     {
2928         gmx_fio_do_int(fio, fgen);
2929     }
2930     else
2931     {
2932         fgen = 0;
2933     }
2934
2935     if (fver >= 80)
2936     {
2937         gmx_fio_do_string(fio, file_tag);
2938     }
2939     if (bRead)
2940     {
2941         if (fver < 77)
2942         {
2943             /* Versions before 77 don't have the tag, set it to release */
2944             sprintf(file_tag, "%s", TPX_TAG_RELEASE);
2945         }
2946
2947         if (strcmp(file_tag, tpx_tag) != 0)
2948         {
2949             fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
2950                     file_tag, tpx_tag);
2951
2952             /* We only support reading tpx files with the same tag as the code
2953              * or tpx files with the release tag and with lower version number.
2954              */
2955             if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
2956             {
2957                 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2958                           gmx_fio_getname(fio), fver, file_tag,
2959                           tpx_version, tpx_tag);
2960             }
2961         }
2962     }
2963
2964     if (file_version != NULL)
2965     {
2966         *file_version = fver;
2967     }
2968     if (file_generation != NULL)
2969     {
2970         *file_generation = fgen;
2971     }
2972
2973
2974     if ((fver <= tpx_incompatible_version) ||
2975         ((fver > tpx_version) && !TopOnlyOK) ||
2976         (fgen > tpx_generation))
2977     {
2978         gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
2979                   gmx_fio_getname(fio), fver, tpx_version);
2980     }
2981
2982     do_section(fio, eitemHEADER, bRead);
2983     gmx_fio_do_int(fio, tpx->natoms);
2984     if (fver >= 28)
2985     {
2986         gmx_fio_do_int(fio, tpx->ngtc);
2987     }
2988     else
2989     {
2990         tpx->ngtc = 0;
2991     }
2992     if (fver < 62)
2993     {
2994         gmx_fio_do_int(fio, idum);
2995         gmx_fio_do_real(fio, rdum);
2996     }
2997     /*a better decision will eventually (5.0 or later) need to be made
2998        on how to treat the alchemical state of the system, which can now
2999        vary through a simulation, and cannot be completely described
3000        though a single lambda variable, or even a single state
3001        index. Eventually, should probably be a vector. MRS*/
3002     if (fver >= 79)
3003     {
3004         gmx_fio_do_int(fio, tpx->fep_state);
3005     }
3006     gmx_fio_do_real(fio, tpx->lambda);
3007     gmx_fio_do_int(fio, tpx->bIr);
3008     gmx_fio_do_int(fio, tpx->bTop);
3009     gmx_fio_do_int(fio, tpx->bX);
3010     gmx_fio_do_int(fio, tpx->bV);
3011     gmx_fio_do_int(fio, tpx->bF);
3012     gmx_fio_do_int(fio, tpx->bBox);
3013
3014     if ((fgen > tpx_generation))
3015     {
3016         /* This can only happen if TopOnlyOK=TRUE */
3017         tpx->bIr = FALSE;
3018     }
3019 }
3020
3021 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3022                   t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3023                   gmx_bool bXVallocated)
3024 {
3025     t_tpxheader     tpx;
3026     t_inputrec      dum_ir;
3027     gmx_mtop_t      dum_top;
3028     gmx_bool        TopOnlyOK, bDum = TRUE;
3029     int             file_version, file_generation;
3030     int             i;
3031     rvec           *xptr, *vptr;
3032     int             ePBC;
3033     gmx_bool        bPeriodicMols;
3034
3035     if (!bRead)
3036     {
3037         tpx.natoms    = state->natoms;
3038         tpx.ngtc      = state->ngtc; /* need to add nnhpres here? */
3039         tpx.fep_state = state->fep_state;
3040         tpx.lambda    = state->lambda[efptFEP];
3041         tpx.bIr       = (ir       != NULL);
3042         tpx.bTop      = (mtop     != NULL);
3043         tpx.bX        = (state->x != NULL);
3044         tpx.bV        = (state->v != NULL);
3045         tpx.bF        = (f        != NULL);
3046         tpx.bBox      = TRUE;
3047     }
3048
3049     TopOnlyOK = (ir == NULL);
3050
3051     do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3052
3053     if (bRead)
3054     {
3055         state->flags  = 0;
3056         /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3057         /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3058         if (bXVallocated)
3059         {
3060             xptr = state->x;
3061             vptr = state->v;
3062             init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3063             state->natoms = tpx.natoms;
3064             state->nalloc = tpx.natoms;
3065             state->x      = xptr;
3066             state->v      = vptr;
3067         }
3068         else
3069         {
3070             init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3071         }
3072     }
3073
3074 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3075
3076     do_test(fio, tpx.bBox, state->box);
3077     do_section(fio, eitemBOX, bRead);
3078     if (tpx.bBox)
3079     {
3080         gmx_fio_ndo_rvec(fio, state->box, DIM);
3081         if (file_version >= 51)
3082         {
3083             gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3084         }
3085         else
3086         {
3087             /* We initialize box_rel after reading the inputrec */
3088             clear_mat(state->box_rel);
3089         }
3090         if (file_version >= 28)
3091         {
3092             gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3093             if (file_version < 56)
3094             {
3095                 matrix mdum;
3096                 gmx_fio_ndo_rvec(fio, mdum, DIM);
3097             }
3098         }
3099     }
3100
3101     if (state->ngtc > 0 && file_version >= 28)
3102     {
3103         real *dumv;
3104         /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3105         /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3106         /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3107         snew(dumv, state->ngtc);
3108         if (file_version < 69)
3109         {
3110             bDum = gmx_fio_ndo_real(fio, dumv, state->ngtc);
3111         }
3112         /* These used to be the Berendsen tcoupl_lambda's */
3113         bDum = gmx_fio_ndo_real(fio, dumv, state->ngtc);
3114         sfree(dumv);
3115     }
3116
3117     /* Prior to tpx version 26, the inputrec was here.
3118      * I moved it to enable partial forward-compatibility
3119      * for analysis/viewer programs.
3120      */
3121     if (file_version < 26)
3122     {
3123         do_test(fio, tpx.bIr, ir);
3124         do_section(fio, eitemIR, bRead);
3125         if (tpx.bIr)
3126         {
3127             if (ir)
3128             {
3129                 do_inputrec(fio, ir, bRead, file_version,
3130                             mtop ? &mtop->ffparams.fudgeQQ : NULL);
3131                 if (bRead && debug)
3132                 {
3133                     pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3134                 }
3135             }
3136             else
3137             {
3138                 do_inputrec(fio, &dum_ir, bRead, file_version,
3139                             mtop ? &mtop->ffparams.fudgeQQ : NULL);
3140                 if (bRead && debug)
3141                 {
3142                     pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3143                 }
3144                 done_inputrec(&dum_ir);
3145             }
3146
3147         }
3148     }
3149
3150     do_test(fio, tpx.bTop, mtop);
3151     do_section(fio, eitemTOP, bRead);
3152     if (tpx.bTop)
3153     {
3154         if (mtop)
3155         {
3156             do_mtop(fio, mtop, bRead, file_version);
3157         }
3158         else
3159         {
3160             do_mtop(fio, &dum_top, bRead, file_version);
3161             done_mtop(&dum_top, TRUE);
3162         }
3163     }
3164     do_test(fio, tpx.bX, state->x);
3165     do_section(fio, eitemX, bRead);
3166     if (tpx.bX)
3167     {
3168         if (bRead)
3169         {
3170             state->flags |= (1<<estX);
3171         }
3172         gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3173     }
3174
3175     do_test(fio, tpx.bV, state->v);
3176     do_section(fio, eitemV, bRead);
3177     if (tpx.bV)
3178     {
3179         if (bRead)
3180         {
3181             state->flags |= (1<<estV);
3182         }
3183         gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3184     }
3185
3186     do_test(fio, tpx.bF, f);
3187     do_section(fio, eitemF, bRead);
3188     if (tpx.bF)
3189     {
3190         gmx_fio_ndo_rvec(fio, f, state->natoms);
3191     }
3192
3193     /* Starting with tpx version 26, we have the inputrec
3194      * at the end of the file, so we can ignore it
3195      * if the file is never than the software (but still the
3196      * same generation - see comments at the top of this file.
3197      *
3198      *
3199      */
3200     ePBC          = -1;
3201     bPeriodicMols = FALSE;
3202     if (file_version >= 26)
3203     {
3204         do_test(fio, tpx.bIr, ir);
3205         do_section(fio, eitemIR, bRead);
3206         if (tpx.bIr)
3207         {
3208             if (file_version >= 53)
3209             {
3210                 /* Removed the pbc info from do_inputrec, since we always want it */
3211                 if (!bRead)
3212                 {
3213                     ePBC          = ir->ePBC;
3214                     bPeriodicMols = ir->bPeriodicMols;
3215                 }
3216                 gmx_fio_do_int(fio, ePBC);
3217                 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3218             }
3219             if (file_generation <= tpx_generation && ir)
3220             {
3221                 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3222                 if (bRead && debug)
3223                 {
3224                     pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3225                 }
3226                 if (file_version < 51)
3227                 {
3228                     set_box_rel(ir, state);
3229                 }
3230                 if (file_version < 53)
3231                 {
3232                     ePBC          = ir->ePBC;
3233                     bPeriodicMols = ir->bPeriodicMols;
3234                 }
3235             }
3236             if (bRead && ir && file_version >= 53)
3237             {
3238                 /* We need to do this after do_inputrec, since that initializes ir */
3239                 ir->ePBC          = ePBC;
3240                 ir->bPeriodicMols = bPeriodicMols;
3241             }
3242         }
3243     }
3244
3245     if (bRead)
3246     {
3247         if (tpx.bIr && ir)
3248         {
3249             if (state->ngtc == 0)
3250             {
3251                 /* Reading old version without tcoupl state data: set it */
3252                 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3253             }
3254             if (tpx.bTop && mtop)
3255             {
3256                 if (file_version < 57)
3257                 {
3258                     if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3259                     {
3260                         ir->eDisre = edrSimple;
3261                     }
3262                     else
3263                     {
3264                         ir->eDisre = edrNone;
3265                     }
3266                 }
3267                 set_disres_npair(mtop);
3268             }
3269         }
3270
3271         if (tpx.bTop && mtop)
3272         {
3273             gmx_mtop_finalize(mtop);
3274         }
3275
3276         if (file_version >= 57)
3277         {
3278             char *env;
3279             int   ienv;
3280             env = getenv("GMX_NOCHARGEGROUPS");
3281             if (env != NULL)
3282             {
3283                 sscanf(env, "%d", &ienv);
3284                 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3285                         ienv);
3286                 if (ienv > 0)
3287                 {
3288                     fprintf(stderr,
3289                             "Will make single atomic charge groups in non-solvent%s\n",
3290                             ienv > 1 ? " and solvent" : "");
3291                     gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3292                 }
3293                 fprintf(stderr, "\n");
3294             }
3295         }
3296     }
3297
3298     return ePBC;
3299 }
3300
3301 /************************************************************
3302  *
3303  *  The following routines are the exported ones
3304  *
3305  ************************************************************/
3306
3307 t_fileio *open_tpx(const char *fn, const char *mode)
3308 {
3309     return gmx_fio_open(fn, mode);
3310 }
3311
3312 void close_tpx(t_fileio *fio)
3313 {
3314     gmx_fio_close(fio);
3315 }
3316
3317 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3318                     int *file_version, int *file_generation)
3319 {
3320     t_fileio *fio;
3321
3322     fio = open_tpx(fn, "r");
3323     do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3324     close_tpx(fio);
3325 }
3326
3327 void write_tpx_state(const char *fn,
3328                      t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3329 {
3330     t_fileio *fio;
3331
3332     fio = open_tpx(fn, "w");
3333     do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3334     close_tpx(fio);
3335 }
3336
3337 void read_tpx_state(const char *fn,
3338                     t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3339 {
3340     t_fileio *fio;
3341
3342     fio = open_tpx(fn, "r");
3343     do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3344     close_tpx(fio);
3345 }
3346
3347 int read_tpx(const char *fn,
3348              t_inputrec *ir, matrix box, int *natoms,
3349              rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3350 {
3351     t_fileio *fio;
3352     t_state   state;
3353     int       ePBC;
3354
3355     state.x = x;
3356     state.v = v;
3357     fio     = open_tpx(fn, "r");
3358     ePBC    = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3359     close_tpx(fio);
3360     *natoms = state.natoms;
3361     if (box)
3362     {
3363         copy_mat(state.box, box);
3364     }
3365     state.x = NULL;
3366     state.v = NULL;
3367     done_state(&state);
3368
3369     return ePBC;
3370 }
3371
3372 int read_tpx_top(const char *fn,
3373                  t_inputrec *ir, matrix box, int *natoms,
3374                  rvec *x, rvec *v, rvec *f, t_topology *top)
3375 {
3376     gmx_mtop_t  mtop;
3377     t_topology *ltop;
3378     int         ePBC;
3379
3380     ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3381
3382     *top = gmx_mtop_t_to_t_topology(&mtop);
3383
3384     return ePBC;
3385 }
3386
3387 gmx_bool fn2bTPX(const char *file)
3388 {
3389     switch (fn2ftp(file))
3390     {
3391         case efTPR:
3392         case efTPB:
3393         case efTPA:
3394             return TRUE;
3395         default:
3396             return FALSE;
3397     }
3398 }
3399
3400 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3401                        rvec **x, rvec **v, matrix box, gmx_bool bMass)
3402 {
3403     t_tpxheader      header;
3404     int              natoms, i, version, generation;
3405     gmx_bool         bTop, bXNULL = FALSE;
3406     gmx_mtop_t      *mtop;
3407     t_topology      *topconv;
3408     gmx_atomprop_t   aps;
3409
3410     bTop  = fn2bTPX(infile);
3411     *ePBC = -1;
3412     if (bTop)
3413     {
3414         read_tpxheader(infile, &header, TRUE, &version, &generation);
3415         if (x)
3416         {
3417             snew(*x, header.natoms);
3418         }
3419         if (v)
3420         {
3421             snew(*v, header.natoms);
3422         }
3423         snew(mtop, 1);
3424         *ePBC = read_tpx(infile, NULL, box, &natoms,
3425                          (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3426         *top = gmx_mtop_t_to_t_topology(mtop);
3427         sfree(mtop);
3428         strcpy(title, *top->name);
3429         tpx_make_chain_identifiers(&top->atoms, &top->mols);
3430     }
3431     else
3432     {
3433         get_stx_coordnum(infile, &natoms);
3434         init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3435         if (x == NULL)
3436         {
3437             snew(x, 1);
3438             bXNULL = TRUE;
3439         }
3440         snew(*x, natoms);
3441         if (v)
3442         {
3443             snew(*v, natoms);
3444         }
3445         read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3446         if (bXNULL)
3447         {
3448             sfree(*x);
3449             sfree(x);
3450         }
3451         if (bMass)
3452         {
3453             aps = gmx_atomprop_init();
3454             for (i = 0; (i < natoms); i++)
3455             {
3456                 if (!gmx_atomprop_query(aps, epropMass,
3457                                         *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3458                                         *top->atoms.atomname[i],
3459                                         &(top->atoms.atom[i].m)))
3460                 {
3461                     if (debug)
3462                     {
3463                         fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3464                                 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3465                                 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3466                                 *top->atoms.atomname[i]);
3467                     }
3468                 }
3469             }
3470             gmx_atomprop_destroy(aps);
3471         }
3472         top->idef.ntypes = -1;
3473     }
3474
3475     return bTop;
3476 }