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