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