1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
34 * GROningen Mixture of Alchemy and Childrens' Stories
40 /* This file is completely threadsafe - keep it that way! */
42 #include <thread_mpi.h>
50 #include "gmx_fatal.h"
64 #include "mtop_util.h"
66 #define TPX_TAG_RELEASE "release"
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.
73 static const char *tpx_tag = TPX_TAG_RELEASE;
75 /* This number should be increased whenever the file format changes! */
76 static const int tpx_version = 79;
78 /* This number should only be increased when you edit the TOPOLOGY section
79 * of the tpx format. This way we can maintain forward compatibility too
80 * for all analysis tools and/or external programs that only need to
81 * know the atom/residue names, charges, and bond connectivity.
83 * It first appeared in tpx version 26, when I also moved the inputrecord
84 * to the end of the tpx file, so we can just skip it if we only
87 static const int tpx_generation = 24;
89 /* This number should be the most recent backwards incompatible version
90 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
92 static const int tpx_incompatible_version = 9;
96 /* Struct used to maintain tpx compatibility when function types are added */
98 int fvnr; /* file version number in which the function type first appeared */
99 int ftype; /* function type */
103 *The entries should be ordered in:
104 * 1. ascending file version number
105 * 2. ascending function type number
107 /*static const t_ftupd ftupd[] = {
108 { 20, F_CUBICBONDS },
112 { 22, F_DISRESVIOL },
118 { 26, F_DIHRESVIOL },
119 { 30, F_CROSS_BOND_BONDS },
120 { 30, F_CROSS_BOND_ANGLES },
121 { 30, F_UREY_BRADLEY },
122 { 30, F_POLARIZATION },
126 *The entries should be ordered in:
127 * 1. ascending function type number
128 * 2. ascending file version number
130 /* question; what is the purpose of the commented code above? */
131 static const t_ftupd ftupd[] = {
132 { 20, F_CUBICBONDS },
137 { 43, F_TABBONDSNC },
138 { 70, F_RESTRBONDS },
139 { 76, F_LINEAR_ANGLES },
140 { 30, F_CROSS_BOND_BONDS },
141 { 30, F_CROSS_BOND_ANGLES },
142 { 30, F_UREY_BRADLEY },
143 { 34, F_QUARTIC_ANGLES },
153 { 72, F_NPSOLVATION },
155 { 41, F_LJC_PAIRS_NB },
158 { 32, F_COUL_RECIP },
160 { 30, F_POLARIZATION },
162 { 22, F_DISRESVIOL },
166 { 26, F_DIHRESVIOL },
171 { 46, F_ECONSERVED },
175 { 76, F_ANHARM_POL },
178 { 79, F_DVDL_BONDED, },
179 { 79, F_DVDL_RESTRAINT },
180 { 79, F_DVDL_TEMPERATURE },
183 #define NFTUPD asize(ftupd)
185 /* Needed for backward compatibility */
188 static void _do_section(t_fileio *fio,int key,gmx_bool bRead,const char *src,
194 if (gmx_fio_getftp(fio) == efTPA) {
196 gmx_fio_write_string(fio,itemstr[key]);
197 bDbg = gmx_fio_getdebug(fio);
198 gmx_fio_setdebug(fio,FALSE);
199 gmx_fio_write_string(fio,comment_str[key]);
200 gmx_fio_setdebug(fio,bDbg);
203 if (gmx_fio_getdebug(fio))
204 fprintf(stderr,"Looking for section %s (%s, %d)",
205 itemstr[key],src,line);
208 gmx_fio_do_string(fio,buf);
209 } while ((gmx_strcasecmp(buf,itemstr[key]) != 0));
211 if (gmx_strcasecmp(buf,itemstr[key]) != 0)
212 gmx_fatal(FARGS,"\nCould not find section heading %s",itemstr[key]);
213 else if (gmx_fio_getdebug(fio))
214 fprintf(stderr," and found it\n");
219 #define do_section(fio,key,bRead) _do_section(fio,key,bRead,__FILE__,__LINE__)
221 /**************************************************************
223 * Now the higer level routines that do io of the structures and arrays
225 **************************************************************/
226 static void do_pullgrp(t_fileio *fio, t_pullgrp *pgrp, gmx_bool bRead,
232 gmx_fio_do_int(fio,pgrp->nat);
234 snew(pgrp->ind,pgrp->nat);
235 bDum=gmx_fio_ndo_int(fio,pgrp->ind,pgrp->nat);
236 gmx_fio_do_int(fio,pgrp->nweight);
238 snew(pgrp->weight,pgrp->nweight);
239 bDum=gmx_fio_ndo_real(fio,pgrp->weight,pgrp->nweight);
240 gmx_fio_do_int(fio,pgrp->pbcatom);
241 gmx_fio_do_rvec(fio,pgrp->vec);
242 gmx_fio_do_rvec(fio,pgrp->init);
243 gmx_fio_do_real(fio,pgrp->rate);
244 gmx_fio_do_real(fio,pgrp->k);
245 if (file_version >= 56) {
246 gmx_fio_do_real(fio,pgrp->kB);
252 static void do_expandedvals(t_fileio *fio,t_expanded *expand,int n_lambda, gmx_bool bRead, int file_version)
254 /* i is used in the ndo_double macro*/
260 if (file_version >= 79)
266 snew(expand->init_lambda_weights,n_lambda);
268 bDum=gmx_fio_ndo_real(fio,expand->init_lambda_weights,n_lambda);
269 gmx_fio_do_gmx_bool(fio,expand->bInit_weights);
272 gmx_fio_do_int(fio,expand->nstexpanded);
273 gmx_fio_do_int(fio,expand->elmcmove);
274 gmx_fio_do_int(fio,expand->elamstats);
275 gmx_fio_do_int(fio,expand->lmc_repeats);
276 gmx_fio_do_int(fio,expand->gibbsdeltalam);
277 gmx_fio_do_int(fio,expand->lmc_forced_nstart);
278 gmx_fio_do_int(fio,expand->lmc_seed);
279 gmx_fio_do_real(fio,expand->mc_temp);
280 gmx_fio_do_int(fio,expand->bSymmetrizedTMatrix);
281 gmx_fio_do_int(fio,expand->nstTij);
282 gmx_fio_do_int(fio,expand->minvarmin);
283 gmx_fio_do_int(fio,expand->c_range);
284 gmx_fio_do_real(fio,expand->wl_scale);
285 gmx_fio_do_real(fio,expand->wl_ratio);
286 gmx_fio_do_real(fio,expand->init_wl_delta);
287 gmx_fio_do_gmx_bool(fio,expand->bWLoneovert);
288 gmx_fio_do_int(fio,expand->elmceq);
289 gmx_fio_do_int(fio,expand->equil_steps);
290 gmx_fio_do_int(fio,expand->equil_samples);
291 gmx_fio_do_int(fio,expand->equil_n_at_lam);
292 gmx_fio_do_real(fio,expand->equil_wl_delta);
293 gmx_fio_do_real(fio,expand->equil_ratio);
297 static void do_simtempvals(t_fileio *fio,t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
302 if (file_version >= 79)
304 gmx_fio_do_int(fio,simtemp->eSimTempScale);
305 gmx_fio_do_real(fio,simtemp->simtemp_high);
306 gmx_fio_do_real(fio,simtemp->simtemp_low);
311 snew(simtemp->temperatures,n_lambda);
313 bDum=gmx_fio_ndo_real(fio,simtemp->temperatures,n_lambda);
318 static void do_fepvals(t_fileio *fio,t_lambda *fepvals,gmx_bool bRead, int file_version)
320 /* i is defined in the ndo_double macro; use g to iterate. */
326 /* free energy values */
327 if (file_version >= 79)
329 gmx_fio_do_int(fio,fepvals->init_fep_state);
330 gmx_fio_do_double(fio,fepvals->init_lambda);
331 gmx_fio_do_double(fio,fepvals->delta_lambda);
333 else if (file_version >= 59) {
334 gmx_fio_do_double(fio,fepvals->init_lambda);
335 gmx_fio_do_double(fio,fepvals->delta_lambda);
337 gmx_fio_do_real(fio,rdum);
338 fepvals->init_lambda = rdum;
339 gmx_fio_do_real(fio,rdum);
340 fepvals->delta_lambda = rdum;
342 if (file_version >= 79)
344 gmx_fio_do_int(fio,fepvals->n_lambda);
347 snew(fepvals->all_lambda,efptNR);
349 for (g=0;g<efptNR;g++)
351 if (fepvals->n_lambda > 0) {
354 snew(fepvals->all_lambda[g],fepvals->n_lambda);
356 bDum=gmx_fio_ndo_double(fio,fepvals->all_lambda[g],fepvals->n_lambda);
357 bDum=gmx_fio_ndo_int(fio,fepvals->separate_dvdl,efptNR);
359 else if (fepvals->init_lambda >= 0)
361 fepvals->separate_dvdl[efptFEP] = TRUE;
365 else if (file_version >= 64)
367 gmx_fio_do_int(fio,fepvals->n_lambda);
368 snew(fepvals->all_lambda,efptNR);
371 snew(fepvals->all_lambda[efptFEP],fepvals->n_lambda);
373 bDum=gmx_fio_ndo_double(fio,fepvals->all_lambda[efptFEP],fepvals->n_lambda);
374 if (fepvals->init_lambda >= 0)
376 fepvals->separate_dvdl[efptFEP] = TRUE;
378 /* still allocate the all_lambda array's contents. */
379 for (g=0;g<efptNR;g++)
381 if (fepvals->n_lambda > 0) {
384 snew(fepvals->all_lambda[g],fepvals->n_lambda);
391 fepvals->n_lambda = 0;
392 fepvals->all_lambda = NULL;
393 if (fepvals->init_lambda >= 0)
395 fepvals->separate_dvdl[efptFEP] = TRUE;
398 if (file_version >= 13)
400 gmx_fio_do_real(fio,fepvals->sc_alpha);
404 fepvals->sc_alpha = 0;
406 if (file_version >= 38)
408 gmx_fio_do_int(fio,fepvals->sc_power);
412 fepvals->sc_power = 2;
414 if (file_version >= 79)
416 gmx_fio_do_real(fio,fepvals->sc_r_power);
420 fepvals->sc_r_power = 6.0;
422 if (file_version >= 15)
424 gmx_fio_do_real(fio,fepvals->sc_sigma);
428 fepvals->sc_sigma = 0.3;
432 if (file_version >= 71)
434 fepvals->sc_sigma_min = fepvals->sc_sigma;
438 fepvals->sc_sigma_min = 0;
441 if (file_version >= 79)
443 gmx_fio_do_int(fio,fepvals->bScCoul);
447 fepvals->bScCoul = TRUE;
449 if (file_version >= 64) {
450 gmx_fio_do_int(fio,fepvals->nstdhdl);
452 fepvals->nstdhdl = 1;
455 if (file_version >= 73)
457 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
458 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
462 fepvals->separate_dhdl_file = esepdhdlfileYES;
463 fepvals->dhdl_derivatives = edhdlderivativesYES;
465 if (file_version >= 71)
467 gmx_fio_do_int(fio,fepvals->dh_hist_size);
468 gmx_fio_do_double(fio,fepvals->dh_hist_spacing);
472 fepvals->dh_hist_size = 0;
473 fepvals->dh_hist_spacing = 0.1;
475 if (file_version >= 79)
477 gmx_fio_do_int(fio,fepvals->bPrintEnergy);
481 fepvals->bPrintEnergy = FALSE;
485 static void do_pull(t_fileio *fio, t_pull *pull,gmx_bool bRead, int file_version)
489 gmx_fio_do_int(fio,pull->ngrp);
490 gmx_fio_do_int(fio,pull->eGeom);
491 gmx_fio_do_ivec(fio,pull->dim);
492 gmx_fio_do_real(fio,pull->cyl_r1);
493 gmx_fio_do_real(fio,pull->cyl_r0);
494 gmx_fio_do_real(fio,pull->constr_tol);
495 gmx_fio_do_int(fio,pull->nstxout);
496 gmx_fio_do_int(fio,pull->nstfout);
498 snew(pull->grp,pull->ngrp+1);
499 for(g=0; g<pull->ngrp+1; g++)
500 do_pullgrp(fio,&pull->grp[g],bRead,file_version);
504 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg,gmx_bool bRead, int file_version)
509 gmx_fio_do_int(fio,rotg->eType);
510 gmx_fio_do_int(fio,rotg->bMassW);
511 gmx_fio_do_int(fio,rotg->nat);
513 snew(rotg->ind,rotg->nat);
514 gmx_fio_ndo_int(fio,rotg->ind,rotg->nat);
516 snew(rotg->x_ref,rotg->nat);
517 gmx_fio_ndo_rvec(fio,rotg->x_ref,rotg->nat);
518 gmx_fio_do_rvec(fio,rotg->vec);
519 gmx_fio_do_rvec(fio,rotg->pivot);
520 gmx_fio_do_real(fio,rotg->rate);
521 gmx_fio_do_real(fio,rotg->k);
522 gmx_fio_do_real(fio,rotg->slab_dist);
523 gmx_fio_do_real(fio,rotg->min_gaussian);
524 gmx_fio_do_real(fio,rotg->eps);
525 gmx_fio_do_int(fio,rotg->eFittype);
526 gmx_fio_do_int(fio,rotg->PotAngle_nstep);
527 gmx_fio_do_real(fio,rotg->PotAngle_step);
530 static void do_rot(t_fileio *fio, t_rot *rot,gmx_bool bRead, int file_version)
534 gmx_fio_do_int(fio,rot->ngrp);
535 gmx_fio_do_int(fio,rot->nstrout);
536 gmx_fio_do_int(fio,rot->nstsout);
538 snew(rot->grp,rot->ngrp);
539 for(g=0; g<rot->ngrp; g++)
540 do_rotgrp(fio, &rot->grp[g],bRead,file_version);
544 static void do_inputrec(t_fileio *fio, t_inputrec *ir,gmx_bool bRead,
545 int file_version, real *fudgeQQ)
547 int i,j,k,*tmp,idum=0;
552 real zerotemptime,finish_t,init_temp,finish_temp;
554 if (file_version != tpx_version)
556 /* Give a warning about features that are not accessible */
557 fprintf(stderr,"Note: file tpx version %d, software tpx version %d\n",
558 file_version,tpx_version);
566 if (file_version == 0)
571 /* Basic inputrec stuff */
572 gmx_fio_do_int(fio,ir->eI);
573 if (file_version >= 62) {
574 gmx_fio_do_gmx_large_int(fio, ir->nsteps);
576 gmx_fio_do_int(fio,idum);
579 if(file_version > 25) {
580 if (file_version >= 62) {
581 gmx_fio_do_gmx_large_int(fio, ir->init_step);
583 gmx_fio_do_int(fio,idum);
584 ir->init_step = idum;
590 if(file_version >= 58)
591 gmx_fio_do_int(fio,ir->simulation_part);
593 ir->simulation_part=1;
595 if (file_version >= 67) {
596 gmx_fio_do_int(fio,ir->nstcalcenergy);
598 ir->nstcalcenergy = 1;
600 if (file_version < 53) {
601 /* The pbc info has been moved out of do_inputrec,
602 * since we always want it, also without reading the inputrec.
604 gmx_fio_do_int(fio,ir->ePBC);
605 if ((file_version <= 15) && (ir->ePBC == 2))
607 if (file_version >= 45) {
608 gmx_fio_do_int(fio,ir->bPeriodicMols);
612 ir->bPeriodicMols = TRUE;
614 ir->bPeriodicMols = FALSE;
618 gmx_fio_do_int(fio,ir->ns_type);
619 gmx_fio_do_int(fio,ir->nstlist);
620 gmx_fio_do_int(fio,ir->ndelta);
621 if (file_version < 41) {
622 gmx_fio_do_int(fio,idum);
623 gmx_fio_do_int(fio,idum);
625 if (file_version >= 45)
626 gmx_fio_do_real(fio,ir->rtpi);
629 gmx_fio_do_int(fio,ir->nstcomm);
630 if (file_version > 34)
631 gmx_fio_do_int(fio,ir->comm_mode);
632 else if (ir->nstcomm < 0)
633 ir->comm_mode = ecmANGULAR;
635 ir->comm_mode = ecmLINEAR;
636 ir->nstcomm = abs(ir->nstcomm);
638 if(file_version > 25)
639 gmx_fio_do_int(fio,ir->nstcheckpoint);
643 gmx_fio_do_int(fio,ir->nstcgsteep);
646 gmx_fio_do_int(fio,ir->nbfgscorr);
650 gmx_fio_do_int(fio,ir->nstlog);
651 gmx_fio_do_int(fio,ir->nstxout);
652 gmx_fio_do_int(fio,ir->nstvout);
653 gmx_fio_do_int(fio,ir->nstfout);
654 gmx_fio_do_int(fio,ir->nstenergy);
655 gmx_fio_do_int(fio,ir->nstxtcout);
656 if (file_version >= 59) {
657 gmx_fio_do_double(fio,ir->init_t);
658 gmx_fio_do_double(fio,ir->delta_t);
660 gmx_fio_do_real(fio,rdum);
662 gmx_fio_do_real(fio,rdum);
665 gmx_fio_do_real(fio,ir->xtcprec);
666 if (file_version < 19) {
667 gmx_fio_do_int(fio,idum);
668 gmx_fio_do_int(fio,idum);
670 if(file_version < 18)
671 gmx_fio_do_int(fio,idum);
672 gmx_fio_do_real(fio,ir->rlist);
673 if (file_version >= 67) {
674 gmx_fio_do_real(fio,ir->rlistlong);
676 gmx_fio_do_int(fio,ir->coulombtype);
677 if (file_version < 32 && ir->coulombtype == eelRF)
678 ir->coulombtype = eelRF_NEC;
679 gmx_fio_do_real(fio,ir->rcoulomb_switch);
680 gmx_fio_do_real(fio,ir->rcoulomb);
681 gmx_fio_do_int(fio,ir->vdwtype);
682 gmx_fio_do_real(fio,ir->rvdw_switch);
683 gmx_fio_do_real(fio,ir->rvdw);
684 if (file_version < 67) {
685 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
687 gmx_fio_do_int(fio,ir->eDispCorr);
688 gmx_fio_do_real(fio,ir->epsilon_r);
689 if (file_version >= 37) {
690 gmx_fio_do_real(fio,ir->epsilon_rf);
692 if (EEL_RF(ir->coulombtype)) {
693 ir->epsilon_rf = ir->epsilon_r;
696 ir->epsilon_rf = 1.0;
699 if (file_version >= 29)
700 gmx_fio_do_real(fio,ir->tabext);
704 if(file_version > 25) {
705 gmx_fio_do_int(fio,ir->gb_algorithm);
706 gmx_fio_do_int(fio,ir->nstgbradii);
707 gmx_fio_do_real(fio,ir->rgbradii);
708 gmx_fio_do_real(fio,ir->gb_saltconc);
709 gmx_fio_do_int(fio,ir->implicit_solvent);
711 ir->gb_algorithm=egbSTILL;
715 ir->implicit_solvent=eisNO;
719 gmx_fio_do_real(fio,ir->gb_epsilon_solvent);
720 gmx_fio_do_real(fio,ir->gb_obc_alpha);
721 gmx_fio_do_real(fio,ir->gb_obc_beta);
722 gmx_fio_do_real(fio,ir->gb_obc_gamma);
725 gmx_fio_do_real(fio,ir->gb_dielectric_offset);
726 gmx_fio_do_int(fio,ir->sa_algorithm);
730 ir->gb_dielectric_offset = 0.009;
731 ir->sa_algorithm = esaAPPROX;
733 gmx_fio_do_real(fio,ir->sa_surface_tension);
735 /* Override sa_surface_tension if it is not changed in the mpd-file */
736 if(ir->sa_surface_tension<0)
738 if(ir->gb_algorithm==egbSTILL)
740 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
742 else if(ir->gb_algorithm==egbHCT || ir->gb_algorithm==egbOBC)
744 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
751 /* Better use sensible values than insane (0.0) ones... */
752 ir->gb_epsilon_solvent = 80;
753 ir->gb_obc_alpha = 1.0;
754 ir->gb_obc_beta = 0.8;
755 ir->gb_obc_gamma = 4.85;
756 ir->sa_surface_tension = 2.092;
760 gmx_fio_do_int(fio,ir->nkx);
761 gmx_fio_do_int(fio,ir->nky);
762 gmx_fio_do_int(fio,ir->nkz);
763 gmx_fio_do_int(fio,ir->pme_order);
764 gmx_fio_do_real(fio,ir->ewald_rtol);
766 if (file_version >=24)
767 gmx_fio_do_int(fio,ir->ewald_geometry);
769 ir->ewald_geometry=eewg3D;
771 if (file_version <=17) {
772 ir->epsilon_surface=0;
773 if (file_version==17)
774 gmx_fio_do_int(fio,idum);
777 gmx_fio_do_real(fio,ir->epsilon_surface);
779 gmx_fio_do_gmx_bool(fio,ir->bOptFFT);
781 gmx_fio_do_gmx_bool(fio,ir->bContinuation);
782 gmx_fio_do_int(fio,ir->etc);
783 /* before version 18, ir->etc was a gmx_bool (ir->btc),
784 * but the values 0 and 1 still mean no and
785 * berendsen temperature coupling, respectively.
787 if (file_version >= 79) {
788 gmx_fio_do_gmx_bool(fio,ir->bPrintNHChains);
790 if (file_version >= 71)
792 gmx_fio_do_int(fio,ir->nsttcouple);
796 ir->nsttcouple = ir->nstcalcenergy;
798 if (file_version <= 15)
800 gmx_fio_do_int(fio,idum);
802 if (file_version <=17)
804 gmx_fio_do_int(fio,ir->epct);
805 if (file_version <= 15)
809 ir->epct = epctSURFACETENSION;
811 gmx_fio_do_int(fio,idum);
814 /* we have removed the NO alternative at the beginning */
818 ir->epct=epctISOTROPIC;
822 ir->epc=epcBERENDSEN;
827 gmx_fio_do_int(fio,ir->epc);
828 gmx_fio_do_int(fio,ir->epct);
830 if (file_version >= 71)
832 gmx_fio_do_int(fio,ir->nstpcouple);
836 ir->nstpcouple = ir->nstcalcenergy;
838 gmx_fio_do_real(fio,ir->tau_p);
839 if (file_version <= 15) {
840 gmx_fio_do_rvec(fio,vdum);
841 clear_mat(ir->ref_p);
843 ir->ref_p[i][i] = vdum[i];
845 gmx_fio_do_rvec(fio,ir->ref_p[XX]);
846 gmx_fio_do_rvec(fio,ir->ref_p[YY]);
847 gmx_fio_do_rvec(fio,ir->ref_p[ZZ]);
849 if (file_version <= 15) {
850 gmx_fio_do_rvec(fio,vdum);
851 clear_mat(ir->compress);
853 ir->compress[i][i] = vdum[i];
856 gmx_fio_do_rvec(fio,ir->compress[XX]);
857 gmx_fio_do_rvec(fio,ir->compress[YY]);
858 gmx_fio_do_rvec(fio,ir->compress[ZZ]);
860 if (file_version >= 47) {
861 gmx_fio_do_int(fio,ir->refcoord_scaling);
862 gmx_fio_do_rvec(fio,ir->posres_com);
863 gmx_fio_do_rvec(fio,ir->posres_comB);
865 ir->refcoord_scaling = erscNO;
866 clear_rvec(ir->posres_com);
867 clear_rvec(ir->posres_comB);
869 if((file_version > 25) && (file_version < 79))
870 gmx_fio_do_int(fio,ir->andersen_seed);
873 if(file_version < 26) {
874 gmx_fio_do_gmx_bool(fio,bSimAnn);
875 gmx_fio_do_real(fio,zerotemptime);
878 if (file_version < 37)
879 gmx_fio_do_real(fio,rdum);
881 gmx_fio_do_real(fio,ir->shake_tol);
882 if (file_version < 54)
883 gmx_fio_do_real(fio,*fudgeQQ);
885 gmx_fio_do_int(fio,ir->efep);
886 if (file_version <= 14 && ir->efep != efepNO)
890 do_fepvals(fio,ir->fepvals,bRead,file_version);
892 if (file_version >= 79)
894 gmx_fio_do_gmx_bool(fio,ir->bSimTemp);
902 ir->bSimTemp = FALSE;
906 do_simtempvals(fio,ir->simtempvals,ir->fepvals->n_lambda,bRead,file_version);
909 if (file_version >= 79)
911 gmx_fio_do_gmx_bool(fio,ir->bExpanded);
914 ir->bExpanded = TRUE;
918 ir->bExpanded = FALSE;
923 do_expandedvals(fio,ir->expandedvals,ir->fepvals->n_lambda,bRead,file_version);
925 if (file_version >= 57) {
926 gmx_fio_do_int(fio,ir->eDisre);
928 gmx_fio_do_int(fio,ir->eDisreWeighting);
929 if (file_version < 22) {
930 if (ir->eDisreWeighting == 0)
931 ir->eDisreWeighting = edrwEqual;
933 ir->eDisreWeighting = edrwConservative;
935 gmx_fio_do_gmx_bool(fio,ir->bDisreMixed);
936 gmx_fio_do_real(fio,ir->dr_fc);
937 gmx_fio_do_real(fio,ir->dr_tau);
938 gmx_fio_do_int(fio,ir->nstdisreout);
939 if (file_version >= 22) {
940 gmx_fio_do_real(fio,ir->orires_fc);
941 gmx_fio_do_real(fio,ir->orires_tau);
942 gmx_fio_do_int(fio,ir->nstorireout);
948 if(file_version >= 26 && file_version < 79) {
949 gmx_fio_do_real(fio,ir->dihre_fc);
950 if (file_version < 56)
952 gmx_fio_do_real(fio,rdum);
953 gmx_fio_do_int(fio,idum);
959 gmx_fio_do_real(fio,ir->em_stepsize);
960 gmx_fio_do_real(fio,ir->em_tol);
961 if (file_version >= 22)
962 gmx_fio_do_gmx_bool(fio,ir->bShakeSOR);
964 ir->bShakeSOR = TRUE;
965 if (file_version >= 11)
966 gmx_fio_do_int(fio,ir->niter);
969 fprintf(stderr,"Note: niter not in run input file, setting it to %d\n",
972 if (file_version >= 21)
973 gmx_fio_do_real(fio,ir->fc_stepsize);
976 gmx_fio_do_int(fio,ir->eConstrAlg);
977 gmx_fio_do_int(fio,ir->nProjOrder);
978 gmx_fio_do_real(fio,ir->LincsWarnAngle);
979 if (file_version <= 14)
980 gmx_fio_do_int(fio,idum);
981 if (file_version >=26)
982 gmx_fio_do_int(fio,ir->nLincsIter);
985 fprintf(stderr,"Note: nLincsIter not in run input file, setting it to %d\n",
988 if (file_version < 33)
989 gmx_fio_do_real(fio,bd_temp);
990 gmx_fio_do_real(fio,ir->bd_fric);
991 gmx_fio_do_int(fio,ir->ld_seed);
992 if (file_version >= 33) {
994 gmx_fio_do_rvec(fio,ir->deform[i]);
997 clear_rvec(ir->deform[i]);
999 if (file_version >= 14)
1000 gmx_fio_do_real(fio,ir->cos_accel);
1003 gmx_fio_do_int(fio,ir->userint1);
1004 gmx_fio_do_int(fio,ir->userint2);
1005 gmx_fio_do_int(fio,ir->userint3);
1006 gmx_fio_do_int(fio,ir->userint4);
1007 gmx_fio_do_real(fio,ir->userreal1);
1008 gmx_fio_do_real(fio,ir->userreal2);
1009 gmx_fio_do_real(fio,ir->userreal3);
1010 gmx_fio_do_real(fio,ir->userreal4);
1013 if (file_version >= 77) {
1014 gmx_fio_do_gmx_bool(fio,ir->bAdress);
1016 if (bRead) snew(ir->adress, 1);
1017 gmx_fio_do_int(fio,ir->adress->type);
1018 gmx_fio_do_real(fio,ir->adress->const_wf);
1019 gmx_fio_do_real(fio,ir->adress->ex_width);
1020 gmx_fio_do_real(fio,ir->adress->hy_width);
1021 gmx_fio_do_int(fio,ir->adress->icor);
1022 gmx_fio_do_int(fio,ir->adress->site);
1023 gmx_fio_do_rvec(fio,ir->adress->refs);
1024 gmx_fio_do_int(fio,ir->adress->n_tf_grps);
1025 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1026 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1027 gmx_fio_do_int(fio,ir->adress->do_hybridpairs);
1029 if (bRead)snew(ir->adress->tf_table_index,ir->adress->n_tf_grps);
1030 if (ir->adress->n_tf_grps > 0) {
1031 bDum=gmx_fio_ndo_int(fio,ir->adress->tf_table_index,ir->adress->n_tf_grps);
1033 if (bRead)snew(ir->adress->group_explicit,ir->adress->n_energy_grps);
1034 if (ir->adress->n_energy_grps > 0) {
1035 bDum=gmx_fio_ndo_int(fio, ir->adress->group_explicit,ir->adress->n_energy_grps);
1039 ir->bAdress = FALSE;
1043 if (file_version >= 48) {
1044 gmx_fio_do_int(fio,ir->ePull);
1045 if (ir->ePull != epullNO) {
1048 do_pull(fio, ir->pull,bRead,file_version);
1051 ir->ePull = epullNO;
1054 /* Enforced rotation */
1055 if (file_version >= 74) {
1056 gmx_fio_do_int(fio,ir->bRot);
1057 if (ir->bRot == TRUE) {
1060 do_rot(fio, ir->rot,bRead,file_version);
1067 gmx_fio_do_int(fio,ir->opts.ngtc);
1068 if (file_version >= 69) {
1069 gmx_fio_do_int(fio,ir->opts.nhchainlength);
1071 ir->opts.nhchainlength = 1;
1073 gmx_fio_do_int(fio,ir->opts.ngacc);
1074 gmx_fio_do_int(fio,ir->opts.ngfrz);
1075 gmx_fio_do_int(fio,ir->opts.ngener);
1078 snew(ir->opts.nrdf, ir->opts.ngtc);
1079 snew(ir->opts.ref_t, ir->opts.ngtc);
1080 snew(ir->opts.annealing, ir->opts.ngtc);
1081 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1082 snew(ir->opts.anneal_time, ir->opts.ngtc);
1083 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1084 snew(ir->opts.tau_t, ir->opts.ngtc);
1085 snew(ir->opts.nFreeze,ir->opts.ngfrz);
1086 snew(ir->opts.acc, ir->opts.ngacc);
1087 snew(ir->opts.egp_flags,ir->opts.ngener*ir->opts.ngener);
1089 if (ir->opts.ngtc > 0) {
1090 if (bRead && file_version<13) {
1091 snew(tmp,ir->opts.ngtc);
1092 bDum=gmx_fio_ndo_int(fio,tmp, ir->opts.ngtc);
1093 for(i=0; i<ir->opts.ngtc; i++)
1094 ir->opts.nrdf[i] = tmp[i];
1097 bDum=gmx_fio_ndo_real(fio,ir->opts.nrdf, ir->opts.ngtc);
1099 bDum=gmx_fio_ndo_real(fio,ir->opts.ref_t,ir->opts.ngtc);
1100 bDum=gmx_fio_ndo_real(fio,ir->opts.tau_t,ir->opts.ngtc);
1101 if (file_version<33 && ir->eI==eiBD) {
1102 for(i=0; i<ir->opts.ngtc; i++)
1103 ir->opts.tau_t[i] = bd_temp;
1106 if (ir->opts.ngfrz > 0)
1107 bDum=gmx_fio_ndo_ivec(fio,ir->opts.nFreeze,ir->opts.ngfrz);
1108 if (ir->opts.ngacc > 0)
1109 gmx_fio_ndo_rvec(fio,ir->opts.acc,ir->opts.ngacc);
1110 if (file_version >= 12)
1111 bDum=gmx_fio_ndo_int(fio,ir->opts.egp_flags,
1112 ir->opts.ngener*ir->opts.ngener);
1114 if(bRead && file_version < 26) {
1115 for(i=0;i<ir->opts.ngtc;i++) {
1117 ir->opts.annealing[i] = eannSINGLE;
1118 ir->opts.anneal_npoints[i] = 2;
1119 snew(ir->opts.anneal_time[i],2);
1120 snew(ir->opts.anneal_temp[i],2);
1121 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1122 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1123 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1124 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1125 ir->opts.anneal_time[i][0] = ir->init_t;
1126 ir->opts.anneal_time[i][1] = finish_t;
1127 ir->opts.anneal_temp[i][0] = init_temp;
1128 ir->opts.anneal_temp[i][1] = finish_temp;
1130 ir->opts.annealing[i] = eannNO;
1131 ir->opts.anneal_npoints[i] = 0;
1135 /* file version 26 or later */
1136 /* First read the lists with annealing and npoints for each group */
1137 bDum=gmx_fio_ndo_int(fio,ir->opts.annealing,ir->opts.ngtc);
1138 bDum=gmx_fio_ndo_int(fio,ir->opts.anneal_npoints,ir->opts.ngtc);
1139 for(j=0;j<(ir->opts.ngtc);j++) {
1140 k=ir->opts.anneal_npoints[j];
1142 snew(ir->opts.anneal_time[j],k);
1143 snew(ir->opts.anneal_temp[j],k);
1145 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_time[j],k);
1146 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_temp[j],k);
1150 if (file_version >= 45) {
1151 gmx_fio_do_int(fio,ir->nwall);
1152 gmx_fio_do_int(fio,ir->wall_type);
1153 if (file_version >= 50)
1154 gmx_fio_do_real(fio,ir->wall_r_linpot);
1156 ir->wall_r_linpot = -1;
1157 gmx_fio_do_int(fio,ir->wall_atomtype[0]);
1158 gmx_fio_do_int(fio,ir->wall_atomtype[1]);
1159 gmx_fio_do_real(fio,ir->wall_density[0]);
1160 gmx_fio_do_real(fio,ir->wall_density[1]);
1161 gmx_fio_do_real(fio,ir->wall_ewald_zfac);
1165 ir->wall_atomtype[0] = -1;
1166 ir->wall_atomtype[1] = -1;
1167 ir->wall_density[0] = 0;
1168 ir->wall_density[1] = 0;
1169 ir->wall_ewald_zfac = 3;
1171 /* Cosine stuff for electric fields */
1172 for(j=0; (j<DIM); j++) {
1173 gmx_fio_do_int(fio,ir->ex[j].n);
1174 gmx_fio_do_int(fio,ir->et[j].n);
1176 snew(ir->ex[j].a, ir->ex[j].n);
1177 snew(ir->ex[j].phi,ir->ex[j].n);
1178 snew(ir->et[j].a, ir->et[j].n);
1179 snew(ir->et[j].phi,ir->et[j].n);
1181 bDum=gmx_fio_ndo_real(fio,ir->ex[j].a, ir->ex[j].n);
1182 bDum=gmx_fio_ndo_real(fio,ir->ex[j].phi,ir->ex[j].n);
1183 bDum=gmx_fio_ndo_real(fio,ir->et[j].a, ir->et[j].n);
1184 bDum=gmx_fio_ndo_real(fio,ir->et[j].phi,ir->et[j].n);
1188 if(file_version>=39){
1189 gmx_fio_do_gmx_bool(fio,ir->bQMMM);
1190 gmx_fio_do_int(fio,ir->QMMMscheme);
1191 gmx_fio_do_real(fio,ir->scalefactor);
1192 gmx_fio_do_int(fio,ir->opts.ngQM);
1194 snew(ir->opts.QMmethod, ir->opts.ngQM);
1195 snew(ir->opts.QMbasis, ir->opts.ngQM);
1196 snew(ir->opts.QMcharge, ir->opts.ngQM);
1197 snew(ir->opts.QMmult, ir->opts.ngQM);
1198 snew(ir->opts.bSH, ir->opts.ngQM);
1199 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1200 snew(ir->opts.CASelectrons,ir->opts.ngQM);
1201 snew(ir->opts.SAon, ir->opts.ngQM);
1202 snew(ir->opts.SAoff, ir->opts.ngQM);
1203 snew(ir->opts.SAsteps, ir->opts.ngQM);
1204 snew(ir->opts.bOPT, ir->opts.ngQM);
1205 snew(ir->opts.bTS, ir->opts.ngQM);
1207 if (ir->opts.ngQM > 0) {
1208 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmethod,ir->opts.ngQM);
1209 bDum=gmx_fio_ndo_int(fio,ir->opts.QMbasis,ir->opts.ngQM);
1210 bDum=gmx_fio_ndo_int(fio,ir->opts.QMcharge,ir->opts.ngQM);
1211 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmult,ir->opts.ngQM);
1212 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bSH,ir->opts.ngQM);
1213 bDum=gmx_fio_ndo_int(fio,ir->opts.CASorbitals,ir->opts.ngQM);
1214 bDum=gmx_fio_ndo_int(fio,ir->opts.CASelectrons,ir->opts.ngQM);
1215 bDum=gmx_fio_ndo_real(fio,ir->opts.SAon,ir->opts.ngQM);
1216 bDum=gmx_fio_ndo_real(fio,ir->opts.SAoff,ir->opts.ngQM);
1217 bDum=gmx_fio_ndo_int(fio,ir->opts.SAsteps,ir->opts.ngQM);
1218 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bOPT,ir->opts.ngQM);
1219 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bTS,ir->opts.ngQM);
1221 /* end of QMMM stuff */
1226 static void do_harm(t_fileio *fio, t_iparams *iparams,gmx_bool bRead)
1228 gmx_fio_do_real(fio,iparams->harmonic.rA);
1229 gmx_fio_do_real(fio,iparams->harmonic.krA);
1230 gmx_fio_do_real(fio,iparams->harmonic.rB);
1231 gmx_fio_do_real(fio,iparams->harmonic.krB);
1234 void do_iparams(t_fileio *fio, t_functype ftype,t_iparams *iparams,
1235 gmx_bool bRead, int file_version)
1242 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1250 do_harm(fio, iparams,bRead);
1251 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead) {
1252 /* Correct incorrect storage of parameters */
1253 iparams->pdihs.phiB = iparams->pdihs.phiA;
1254 iparams->pdihs.cpB = iparams->pdihs.cpA;
1257 case F_LINEAR_ANGLES:
1258 gmx_fio_do_real(fio,iparams->linangle.klinA);
1259 gmx_fio_do_real(fio,iparams->linangle.aA);
1260 gmx_fio_do_real(fio,iparams->linangle.klinB);
1261 gmx_fio_do_real(fio,iparams->linangle.aB);
1264 gmx_fio_do_real(fio,iparams->fene.bm);
1265 gmx_fio_do_real(fio,iparams->fene.kb);
1268 gmx_fio_do_real(fio,iparams->restraint.lowA);
1269 gmx_fio_do_real(fio,iparams->restraint.up1A);
1270 gmx_fio_do_real(fio,iparams->restraint.up2A);
1271 gmx_fio_do_real(fio,iparams->restraint.kA);
1272 gmx_fio_do_real(fio,iparams->restraint.lowB);
1273 gmx_fio_do_real(fio,iparams->restraint.up1B);
1274 gmx_fio_do_real(fio,iparams->restraint.up2B);
1275 gmx_fio_do_real(fio,iparams->restraint.kB);
1281 gmx_fio_do_real(fio,iparams->tab.kA);
1282 gmx_fio_do_int(fio,iparams->tab.table);
1283 gmx_fio_do_real(fio,iparams->tab.kB);
1285 case F_CROSS_BOND_BONDS:
1286 gmx_fio_do_real(fio,iparams->cross_bb.r1e);
1287 gmx_fio_do_real(fio,iparams->cross_bb.r2e);
1288 gmx_fio_do_real(fio,iparams->cross_bb.krr);
1290 case F_CROSS_BOND_ANGLES:
1291 gmx_fio_do_real(fio,iparams->cross_ba.r1e);
1292 gmx_fio_do_real(fio,iparams->cross_ba.r2e);
1293 gmx_fio_do_real(fio,iparams->cross_ba.r3e);
1294 gmx_fio_do_real(fio,iparams->cross_ba.krt);
1296 case F_UREY_BRADLEY:
1297 gmx_fio_do_real(fio,iparams->u_b.thetaA);
1298 gmx_fio_do_real(fio,iparams->u_b.kthetaA);
1299 gmx_fio_do_real(fio,iparams->u_b.r13A);
1300 gmx_fio_do_real(fio,iparams->u_b.kUBA);
1301 if (file_version >= 79) {
1302 gmx_fio_do_real(fio,iparams->u_b.thetaB);
1303 gmx_fio_do_real(fio,iparams->u_b.kthetaB);
1304 gmx_fio_do_real(fio,iparams->u_b.r13B);
1305 gmx_fio_do_real(fio,iparams->u_b.kUBB);
1307 iparams->u_b.thetaB=iparams->u_b.thetaA;
1308 iparams->u_b.kthetaB=iparams->u_b.kthetaA;
1309 iparams->u_b.r13B=iparams->u_b.r13A;
1310 iparams->u_b.kUBB=iparams->u_b.kUBA;
1313 case F_QUARTIC_ANGLES:
1314 gmx_fio_do_real(fio,iparams->qangle.theta);
1315 bDum=gmx_fio_ndo_real(fio,iparams->qangle.c,5);
1318 gmx_fio_do_real(fio,iparams->bham.a);
1319 gmx_fio_do_real(fio,iparams->bham.b);
1320 gmx_fio_do_real(fio,iparams->bham.c);
1323 gmx_fio_do_real(fio,iparams->morse.b0A);
1324 gmx_fio_do_real(fio,iparams->morse.cbA);
1325 gmx_fio_do_real(fio,iparams->morse.betaA);
1326 if (file_version >= 79) {
1327 gmx_fio_do_real(fio,iparams->morse.b0B);
1328 gmx_fio_do_real(fio,iparams->morse.cbB);
1329 gmx_fio_do_real(fio,iparams->morse.betaB);
1331 iparams->morse.b0B = iparams->morse.b0A;
1332 iparams->morse.cbB = iparams->morse.cbA;
1333 iparams->morse.betaB = iparams->morse.betaA;
1337 gmx_fio_do_real(fio,iparams->cubic.b0);
1338 gmx_fio_do_real(fio,iparams->cubic.kb);
1339 gmx_fio_do_real(fio,iparams->cubic.kcub);
1343 case F_POLARIZATION:
1344 gmx_fio_do_real(fio,iparams->polarize.alpha);
1347 gmx_fio_do_real(fio,iparams->anharm_polarize.alpha);
1348 gmx_fio_do_real(fio,iparams->anharm_polarize.drcut);
1349 gmx_fio_do_real(fio,iparams->anharm_polarize.khyp);
1352 if (file_version < 31)
1353 gmx_fatal(FARGS,"Old tpr files with water_polarization not supported. Make a new.");
1354 gmx_fio_do_real(fio,iparams->wpol.al_x);
1355 gmx_fio_do_real(fio,iparams->wpol.al_y);
1356 gmx_fio_do_real(fio,iparams->wpol.al_z);
1357 gmx_fio_do_real(fio,iparams->wpol.rOH);
1358 gmx_fio_do_real(fio,iparams->wpol.rHH);
1359 gmx_fio_do_real(fio,iparams->wpol.rOD);
1362 gmx_fio_do_real(fio,iparams->thole.a);
1363 gmx_fio_do_real(fio,iparams->thole.alpha1);
1364 gmx_fio_do_real(fio,iparams->thole.alpha2);
1365 gmx_fio_do_real(fio,iparams->thole.rfac);
1368 gmx_fio_do_real(fio,iparams->lj.c6);
1369 gmx_fio_do_real(fio,iparams->lj.c12);
1372 gmx_fio_do_real(fio,iparams->lj14.c6A);
1373 gmx_fio_do_real(fio,iparams->lj14.c12A);
1374 gmx_fio_do_real(fio,iparams->lj14.c6B);
1375 gmx_fio_do_real(fio,iparams->lj14.c12B);
1378 gmx_fio_do_real(fio,iparams->ljc14.fqq);
1379 gmx_fio_do_real(fio,iparams->ljc14.qi);
1380 gmx_fio_do_real(fio,iparams->ljc14.qj);
1381 gmx_fio_do_real(fio,iparams->ljc14.c6);
1382 gmx_fio_do_real(fio,iparams->ljc14.c12);
1384 case F_LJC_PAIRS_NB:
1385 gmx_fio_do_real(fio,iparams->ljcnb.qi);
1386 gmx_fio_do_real(fio,iparams->ljcnb.qj);
1387 gmx_fio_do_real(fio,iparams->ljcnb.c6);
1388 gmx_fio_do_real(fio,iparams->ljcnb.c12);
1394 gmx_fio_do_real(fio,iparams->pdihs.phiA);
1395 gmx_fio_do_real(fio,iparams->pdihs.cpA);
1396 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42) {
1397 /* Read the incorrectly stored multiplicity */
1398 gmx_fio_do_real(fio,iparams->harmonic.rB);
1399 gmx_fio_do_real(fio,iparams->harmonic.krB);
1400 iparams->pdihs.phiB = iparams->pdihs.phiA;
1401 iparams->pdihs.cpB = iparams->pdihs.cpA;
1403 gmx_fio_do_real(fio,iparams->pdihs.phiB);
1404 gmx_fio_do_real(fio,iparams->pdihs.cpB);
1405 gmx_fio_do_int(fio,iparams->pdihs.mult);
1409 gmx_fio_do_int(fio,iparams->disres.label);
1410 gmx_fio_do_int(fio,iparams->disres.type);
1411 gmx_fio_do_real(fio,iparams->disres.low);
1412 gmx_fio_do_real(fio,iparams->disres.up1);
1413 gmx_fio_do_real(fio,iparams->disres.up2);
1414 gmx_fio_do_real(fio,iparams->disres.kfac);
1417 gmx_fio_do_int(fio,iparams->orires.ex);
1418 gmx_fio_do_int(fio,iparams->orires.label);
1419 gmx_fio_do_int(fio,iparams->orires.power);
1420 gmx_fio_do_real(fio,iparams->orires.c);
1421 gmx_fio_do_real(fio,iparams->orires.obs);
1422 gmx_fio_do_real(fio,iparams->orires.kfac);
1425 gmx_fio_do_real(fio,iparams->dihres.phiA);
1426 gmx_fio_do_real(fio,iparams->dihres.dphiA);
1427 gmx_fio_do_real(fio,iparams->dihres.kfacA);
1428 if (file_version >= 72) {
1429 gmx_fio_do_real(fio,iparams->dihres.phiB);
1430 gmx_fio_do_real(fio,iparams->dihres.dphiB);
1431 gmx_fio_do_real(fio,iparams->dihres.kfacB);
1433 iparams->dihres.phiB=iparams->dihres.phiA;
1434 iparams->dihres.dphiB=iparams->dihres.dphiA;
1435 iparams->dihres.kfacB=iparams->dihres.kfacA;
1439 gmx_fio_do_rvec(fio,iparams->posres.pos0A);
1440 gmx_fio_do_rvec(fio,iparams->posres.fcA);
1441 if (bRead && file_version < 27) {
1442 copy_rvec(iparams->posres.pos0A,iparams->posres.pos0B);
1443 copy_rvec(iparams->posres.fcA,iparams->posres.fcB);
1445 gmx_fio_do_rvec(fio,iparams->posres.pos0B);
1446 gmx_fio_do_rvec(fio,iparams->posres.fcB);
1450 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA,NR_RBDIHS);
1451 if(file_version>=25)
1452 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB,NR_RBDIHS);
1455 /* Fourier dihedrals are internally represented
1456 * as Ryckaert-Bellemans since those are faster to compute.
1458 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA, NR_RBDIHS);
1459 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB, NR_RBDIHS);
1463 gmx_fio_do_real(fio,iparams->constr.dA);
1464 gmx_fio_do_real(fio,iparams->constr.dB);
1467 gmx_fio_do_real(fio,iparams->settle.doh);
1468 gmx_fio_do_real(fio,iparams->settle.dhh);
1471 gmx_fio_do_real(fio,iparams->vsite.a);
1476 gmx_fio_do_real(fio,iparams->vsite.a);
1477 gmx_fio_do_real(fio,iparams->vsite.b);
1482 gmx_fio_do_real(fio,iparams->vsite.a);
1483 gmx_fio_do_real(fio,iparams->vsite.b);
1484 gmx_fio_do_real(fio,iparams->vsite.c);
1487 gmx_fio_do_int(fio,iparams->vsiten.n);
1488 gmx_fio_do_real(fio,iparams->vsiten.a);
1493 /* We got rid of some parameters in version 68 */
1494 if(bRead && file_version<68)
1496 gmx_fio_do_real(fio,rdum);
1497 gmx_fio_do_real(fio,rdum);
1498 gmx_fio_do_real(fio,rdum);
1499 gmx_fio_do_real(fio,rdum);
1501 gmx_fio_do_real(fio,iparams->gb.sar);
1502 gmx_fio_do_real(fio,iparams->gb.st);
1503 gmx_fio_do_real(fio,iparams->gb.pi);
1504 gmx_fio_do_real(fio,iparams->gb.gbr);
1505 gmx_fio_do_real(fio,iparams->gb.bmlt);
1508 gmx_fio_do_int(fio,iparams->cmap.cmapA);
1509 gmx_fio_do_int(fio,iparams->cmap.cmapB);
1512 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
1513 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
1516 gmx_fio_unset_comment(fio);
1519 static void do_ilist(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,int file_version,
1526 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1528 if (file_version < 44) {
1529 for(i=0; i<MAXNODES; i++)
1530 gmx_fio_do_int(fio,idum);
1532 gmx_fio_do_int(fio,ilist->nr);
1534 snew(ilist->iatoms,ilist->nr);
1535 bDum=gmx_fio_ndo_int(fio,ilist->iatoms,ilist->nr);
1537 gmx_fio_unset_comment(fio);
1540 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1541 gmx_bool bRead, int file_version)
1547 gmx_fio_do_int(fio,ffparams->atnr);
1548 if (file_version < 57) {
1549 gmx_fio_do_int(fio,idum);
1551 gmx_fio_do_int(fio,ffparams->ntypes);
1553 fprintf(debug,"ffparams->atnr = %d, ntypes = %d\n",
1554 ffparams->atnr,ffparams->ntypes);
1556 snew(ffparams->functype,ffparams->ntypes);
1557 snew(ffparams->iparams,ffparams->ntypes);
1559 /* Read/write all the function types */
1560 bDum=gmx_fio_ndo_int(fio,ffparams->functype,ffparams->ntypes);
1562 pr_ivec(debug,0,"functype",ffparams->functype,ffparams->ntypes,TRUE);
1564 if (file_version >= 66) {
1565 gmx_fio_do_double(fio,ffparams->reppow);
1567 ffparams->reppow = 12.0;
1570 if (file_version >= 57) {
1571 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1574 /* Check whether all these function types are supported by the code.
1575 * In practice the code is backwards compatible, which means that the
1576 * numbering may have to be altered from old numbering to new numbering
1578 for (i=0; (i<ffparams->ntypes); i++) {
1580 /* Loop over file versions */
1581 for (k=0; (k<NFTUPD); k++)
1582 /* Compare the read file_version to the update table */
1583 if ((file_version < ftupd[k].fvnr) &&
1584 (ffparams->functype[i] >= ftupd[k].ftype)) {
1585 ffparams->functype[i] += 1;
1587 fprintf(debug,"Incrementing function type %d to %d (due to %s)\n",
1588 i,ffparams->functype[i],
1589 interaction_function[ftupd[k].ftype].longname);
1594 do_iparams(fio, ffparams->functype[i],&ffparams->iparams[i],bRead,
1597 pr_iparams(debug,ffparams->functype[i],&ffparams->iparams[i]);
1601 static void add_settle_atoms(t_ilist *ilist)
1605 /* Settle used to only store the first atom: add the other two */
1606 srenew(ilist->iatoms,2*ilist->nr);
1607 for(i=ilist->nr/2-1; i>=0; i--)
1609 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
1610 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
1611 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
1612 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
1614 ilist->nr = 2*ilist->nr;
1617 static void do_ilists(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,
1620 int i,j,renum[F_NRE];
1621 gmx_bool bDum=TRUE,bClear;
1624 for(j=0; (j<F_NRE); j++) {
1627 for (k=0; k<NFTUPD; k++)
1628 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
1632 ilist[j].iatoms = NULL;
1634 do_ilist(fio, &ilist[j],bRead,file_version,j);
1635 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
1637 add_settle_atoms(&ilist[j]);
1641 if (bRead && gmx_debug_at)
1642 pr_ilist(debug,0,interaction_function[j].longname,
1643 functype,&ilist[j],TRUE);
1648 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams,gmx_moltype_t *molt,
1649 gmx_bool bRead, int file_version)
1651 do_ffparams(fio, ffparams,bRead,file_version);
1653 if (file_version >= 54) {
1654 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1657 do_ilists(fio, molt->ilist,bRead,file_version);
1660 static void do_block(t_fileio *fio, t_block *block,gmx_bool bRead,int file_version)
1662 int i,idum,dum_nra,*dum_a;
1665 if (file_version < 44)
1666 for(i=0; i<MAXNODES; i++)
1667 gmx_fio_do_int(fio,idum);
1668 gmx_fio_do_int(fio,block->nr);
1669 if (file_version < 51)
1670 gmx_fio_do_int(fio,dum_nra);
1672 block->nalloc_index = block->nr+1;
1673 snew(block->index,block->nalloc_index);
1675 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1677 if (file_version < 51 && dum_nra > 0) {
1678 snew(dum_a,dum_nra);
1679 bDum=gmx_fio_ndo_int(fio,dum_a,dum_nra);
1684 static void do_blocka(t_fileio *fio, t_blocka *block,gmx_bool bRead,
1690 if (file_version < 44)
1691 for(i=0; i<MAXNODES; i++)
1692 gmx_fio_do_int(fio,idum);
1693 gmx_fio_do_int(fio,block->nr);
1694 gmx_fio_do_int(fio,block->nra);
1696 block->nalloc_index = block->nr+1;
1697 snew(block->index,block->nalloc_index);
1698 block->nalloc_a = block->nra;
1699 snew(block->a,block->nalloc_a);
1701 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1702 bDum=gmx_fio_ndo_int(fio,block->a,block->nra);
1705 static void do_atom(t_fileio *fio, t_atom *atom,int ngrp,gmx_bool bRead,
1706 int file_version, gmx_groups_t *groups,int atnr)
1710 gmx_fio_do_real(fio,atom->m);
1711 gmx_fio_do_real(fio,atom->q);
1712 gmx_fio_do_real(fio,atom->mB);
1713 gmx_fio_do_real(fio,atom->qB);
1714 gmx_fio_do_ushort(fio, atom->type);
1715 gmx_fio_do_ushort(fio, atom->typeB);
1716 gmx_fio_do_int(fio,atom->ptype);
1717 gmx_fio_do_int(fio,atom->resind);
1718 if (file_version >= 52)
1719 gmx_fio_do_int(fio,atom->atomnumber);
1721 atom->atomnumber = NOTSET;
1722 if (file_version < 23)
1724 else if (file_version < 39)
1729 if (file_version < 57) {
1730 unsigned char uchar[egcNR];
1731 gmx_fio_ndo_uchar(fio,uchar,myngrp);
1732 for(i=myngrp; (i<ngrp); i++) {
1735 /* Copy the old data format to the groups struct */
1736 for(i=0; i<ngrp; i++) {
1737 groups->grpnr[i][atnr] = uchar[i];
1742 static void do_grps(t_fileio *fio, int ngrp,t_grps grps[],gmx_bool bRead,
1748 if (file_version < 23)
1750 else if (file_version < 39)
1755 for(j=0; (j<ngrp); j++) {
1757 gmx_fio_do_int(fio,grps[j].nr);
1759 snew(grps[j].nm_ind,grps[j].nr);
1760 bDum=gmx_fio_ndo_int(fio,grps[j].nm_ind,grps[j].nr);
1764 snew(grps[j].nm_ind,grps[j].nr);
1769 static void do_symstr(t_fileio *fio, char ***nm,gmx_bool bRead,t_symtab *symtab)
1774 gmx_fio_do_int(fio,ls);
1775 *nm = get_symtab_handle(symtab,ls);
1778 ls = lookup_symtab(symtab,*nm);
1779 gmx_fio_do_int(fio,ls);
1783 static void do_strstr(t_fileio *fio, int nstr,char ***nm,gmx_bool bRead,
1788 for (j=0; (j<nstr); j++)
1789 do_symstr(fio, &(nm[j]),bRead,symtab);
1792 static void do_resinfo(t_fileio *fio, int n,t_resinfo *ri,gmx_bool bRead,
1793 t_symtab *symtab, int file_version)
1797 for (j=0; (j<n); j++) {
1798 do_symstr(fio, &(ri[j].name),bRead,symtab);
1799 if (file_version >= 63) {
1800 gmx_fio_do_int(fio,ri[j].nr);
1801 gmx_fio_do_uchar(fio, ri[j].ic);
1809 static void do_atoms(t_fileio *fio, t_atoms *atoms,gmx_bool bRead,t_symtab *symtab,
1811 gmx_groups_t *groups)
1815 gmx_fio_do_int(fio,atoms->nr);
1816 gmx_fio_do_int(fio,atoms->nres);
1817 if (file_version < 57) {
1818 gmx_fio_do_int(fio,groups->ngrpname);
1819 for(i=0; i<egcNR; i++) {
1820 groups->ngrpnr[i] = atoms->nr;
1821 snew(groups->grpnr[i],groups->ngrpnr[i]);
1825 snew(atoms->atom,atoms->nr);
1826 snew(atoms->atomname,atoms->nr);
1827 snew(atoms->atomtype,atoms->nr);
1828 snew(atoms->atomtypeB,atoms->nr);
1829 snew(atoms->resinfo,atoms->nres);
1830 if (file_version < 57) {
1831 snew(groups->grpname,groups->ngrpname);
1833 atoms->pdbinfo = NULL;
1835 for(i=0; (i<atoms->nr); i++) {
1836 do_atom(fio, &atoms->atom[i],egcNR,bRead, file_version,groups,i);
1838 do_strstr(fio, atoms->nr,atoms->atomname,bRead,symtab);
1839 if (bRead && (file_version <= 20)) {
1840 for(i=0; i<atoms->nr; i++) {
1841 atoms->atomtype[i] = put_symtab(symtab,"?");
1842 atoms->atomtypeB[i] = put_symtab(symtab,"?");
1845 do_strstr(fio, atoms->nr,atoms->atomtype,bRead,symtab);
1846 do_strstr(fio, atoms->nr,atoms->atomtypeB,bRead,symtab);
1848 do_resinfo(fio, atoms->nres,atoms->resinfo,bRead,symtab,file_version);
1850 if (file_version < 57) {
1851 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1853 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1857 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
1858 gmx_bool bRead,t_symtab *symtab,
1864 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1865 gmx_fio_do_int(fio,groups->ngrpname);
1867 snew(groups->grpname,groups->ngrpname);
1869 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1870 for(g=0; g<egcNR; g++) {
1871 gmx_fio_do_int(fio,groups->ngrpnr[g]);
1872 if (groups->ngrpnr[g] == 0) {
1874 groups->grpnr[g] = NULL;
1878 snew(groups->grpnr[g],groups->ngrpnr[g]);
1880 bDum=gmx_fio_ndo_uchar(fio, groups->grpnr[g],groups->ngrpnr[g]);
1885 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes,gmx_bool bRead,
1886 t_symtab *symtab,int file_version)
1889 gmx_bool bDum = TRUE;
1891 if (file_version > 25) {
1892 gmx_fio_do_int(fio,atomtypes->nr);
1895 snew(atomtypes->radius,j);
1896 snew(atomtypes->vol,j);
1897 snew(atomtypes->surftens,j);
1898 snew(atomtypes->atomnumber,j);
1899 snew(atomtypes->gb_radius,j);
1900 snew(atomtypes->S_hct,j);
1902 bDum=gmx_fio_ndo_real(fio,atomtypes->radius,j);
1903 bDum=gmx_fio_ndo_real(fio,atomtypes->vol,j);
1904 bDum=gmx_fio_ndo_real(fio,atomtypes->surftens,j);
1905 if(file_version >= 40)
1907 bDum=gmx_fio_ndo_int(fio,atomtypes->atomnumber,j);
1909 if(file_version >= 60)
1911 bDum=gmx_fio_ndo_real(fio,atomtypes->gb_radius,j);
1912 bDum=gmx_fio_ndo_real(fio,atomtypes->S_hct,j);
1915 /* File versions prior to 26 cannot do GBSA,
1916 * so they dont use this structure
1919 atomtypes->radius = NULL;
1920 atomtypes->vol = NULL;
1921 atomtypes->surftens = NULL;
1922 atomtypes->atomnumber = NULL;
1923 atomtypes->gb_radius = NULL;
1924 atomtypes->S_hct = NULL;
1928 static void do_symtab(t_fileio *fio, t_symtab *symtab,gmx_bool bRead)
1934 gmx_fio_do_int(fio,symtab->nr);
1937 snew(symtab->symbuf,1);
1938 symbuf = symtab->symbuf;
1939 symbuf->bufsize = nr;
1940 snew(symbuf->buf,nr);
1941 for (i=0; (i<nr); i++) {
1942 gmx_fio_do_string(fio,buf);
1943 symbuf->buf[i]=strdup(buf);
1947 symbuf = symtab->symbuf;
1948 while (symbuf!=NULL) {
1949 for (i=0; (i<symbuf->bufsize) && (i<nr); i++)
1950 gmx_fio_do_string(fio,symbuf->buf[i]);
1952 symbuf=symbuf->next;
1955 gmx_fatal(FARGS,"nr of symtab strings left: %d",nr);
1959 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
1961 int i,j,ngrid,gs,nelem;
1963 gmx_fio_do_int(fio,cmap_grid->ngrid);
1964 gmx_fio_do_int(fio,cmap_grid->grid_spacing);
1966 ngrid = cmap_grid->ngrid;
1967 gs = cmap_grid->grid_spacing;
1972 snew(cmap_grid->cmapdata,ngrid);
1974 for(i=0;i<cmap_grid->ngrid;i++)
1976 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
1980 for(i=0;i<cmap_grid->ngrid;i++)
1982 for(j=0;j<nelem;j++)
1984 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4]);
1985 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+1]);
1986 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+2]);
1987 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+3]);
1993 void tpx_make_chain_identifiers(t_atoms *atoms,t_block *mols)
1999 /* We always assign a new chain number, but save the chain id characters
2000 * for larger molecules.
2002 #define CHAIN_MIN_ATOMS 15
2006 for(m=0; m<mols->nr; m++)
2009 a1=mols->index[m+1];
2010 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2019 for(a=a0; a<a1; a++)
2021 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2022 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2027 /* Blank out the chain id if there was only one chain */
2030 for(r=0; r<atoms->nres; r++)
2032 atoms->resinfo[r].chainid = ' ';
2037 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt,gmx_bool bRead,
2038 t_symtab *symtab, int file_version,
2039 gmx_groups_t *groups)
2043 if (file_version >= 57) {
2044 do_symstr(fio, &(molt->name),bRead,symtab);
2047 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2049 if (bRead && gmx_debug_at) {
2050 pr_atoms(debug,0,"atoms",&molt->atoms,TRUE);
2053 if (file_version >= 57) {
2054 do_ilists(fio, molt->ilist,bRead,file_version);
2056 do_block(fio, &molt->cgs,bRead,file_version);
2057 if (bRead && gmx_debug_at) {
2058 pr_block(debug,0,"cgs",&molt->cgs,TRUE);
2062 /* This used to be in the atoms struct */
2063 do_blocka(fio, &molt->excls, bRead, file_version);
2066 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb,gmx_bool bRead,
2071 gmx_fio_do_int(fio,molb->type);
2072 gmx_fio_do_int(fio,molb->nmol);
2073 gmx_fio_do_int(fio,molb->natoms_mol);
2074 /* Position restraint coordinates */
2075 gmx_fio_do_int(fio,molb->nposres_xA);
2076 if (molb->nposres_xA > 0) {
2078 snew(molb->posres_xA,molb->nposres_xA);
2080 gmx_fio_ndo_rvec(fio,molb->posres_xA,molb->nposres_xA);
2082 gmx_fio_do_int(fio,molb->nposres_xB);
2083 if (molb->nposres_xB > 0) {
2085 snew(molb->posres_xB,molb->nposres_xB);
2087 gmx_fio_ndo_rvec(fio,molb->posres_xB,molb->nposres_xB);
2092 static t_block mtop_mols(gmx_mtop_t *mtop)
2098 for(mb=0; mb<mtop->nmolblock; mb++) {
2099 mols.nr += mtop->molblock[mb].nmol;
2101 mols.nalloc_index = mols.nr + 1;
2102 snew(mols.index,mols.nalloc_index);
2107 for(mb=0; mb<mtop->nmolblock; mb++) {
2108 for(mol=0; mol<mtop->molblock[mb].nmol; mol++) {
2109 a += mtop->molblock[mb].natoms_mol;
2118 static void add_posres_molblock(gmx_mtop_t *mtop)
2123 gmx_molblock_t *molb;
2126 il = &mtop->moltype[0].ilist[F_POSRES];
2132 for(i=0; i<il->nr; i+=2) {
2133 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2134 am = max(am,il->iatoms[i+1]);
2135 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2136 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2137 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ]) {
2141 /* Make the posres coordinate block end at a molecule end */
2143 while(am >= mtop->mols.index[mol+1]) {
2146 molb = &mtop->molblock[0];
2147 molb->nposres_xA = mtop->mols.index[mol+1];
2148 snew(molb->posres_xA,molb->nposres_xA);
2150 molb->nposres_xB = molb->nposres_xA;
2151 snew(molb->posres_xB,molb->nposres_xB);
2153 molb->nposres_xB = 0;
2155 for(i=0; i<il->nr; i+=2) {
2156 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2157 a = il->iatoms[i+1];
2158 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2159 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2160 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2162 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2163 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2164 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2169 static void set_disres_npair(gmx_mtop_t *mtop)
2176 ip = mtop->ffparams.iparams;
2178 for(mt=0; mt<mtop->nmoltype; mt++) {
2179 il = &mtop->moltype[mt].ilist[F_DISRES];
2183 for(i=0; i<il->nr; i+=3) {
2185 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label) {
2186 ip[a[i]].disres.npair = npair;
2194 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop,gmx_bool bRead,
2202 do_symtab(fio, &(mtop->symtab),bRead);
2204 pr_symtab(debug,0,"symtab",&mtop->symtab);
2206 do_symstr(fio, &(mtop->name),bRead,&(mtop->symtab));
2208 if (file_version >= 57) {
2209 do_ffparams(fio, &mtop->ffparams,bRead,file_version);
2211 gmx_fio_do_int(fio,mtop->nmoltype);
2216 snew(mtop->moltype,mtop->nmoltype);
2217 if (file_version < 57) {
2218 mtop->moltype[0].name = mtop->name;
2221 for(mt=0; mt<mtop->nmoltype; mt++) {
2222 do_moltype(fio, &mtop->moltype[mt],bRead,&mtop->symtab,file_version,
2226 if (file_version >= 57) {
2227 gmx_fio_do_int(fio,mtop->nmolblock);
2229 mtop->nmolblock = 1;
2232 snew(mtop->molblock,mtop->nmolblock);
2234 if (file_version >= 57) {
2235 for(mb=0; mb<mtop->nmolblock; mb++) {
2236 do_molblock(fio, &mtop->molblock[mb],bRead,file_version);
2238 gmx_fio_do_int(fio,mtop->natoms);
2240 mtop->molblock[0].type = 0;
2241 mtop->molblock[0].nmol = 1;
2242 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
2243 mtop->molblock[0].nposres_xA = 0;
2244 mtop->molblock[0].nposres_xB = 0;
2247 do_atomtypes (fio, &(mtop->atomtypes),bRead,&(mtop->symtab), file_version);
2249 pr_atomtypes(debug,0,"atomtypes",&mtop->atomtypes,TRUE);
2251 if (file_version < 57) {
2252 /* Debug statements are inside do_idef */
2253 do_idef (fio, &mtop->ffparams,&mtop->moltype[0],bRead,file_version);
2254 mtop->natoms = mtop->moltype[0].atoms.nr;
2257 if(file_version >= 65)
2259 do_cmap(fio, &mtop->ffparams.cmap_grid,bRead);
2263 mtop->ffparams.cmap_grid.ngrid = 0;
2264 mtop->ffparams.cmap_grid.grid_spacing = 0;
2265 mtop->ffparams.cmap_grid.cmapdata = NULL;
2268 if (file_version >= 57) {
2269 do_groups(fio, &mtop->groups,bRead,&(mtop->symtab),file_version);
2272 if (file_version < 57) {
2273 do_block(fio, &mtop->moltype[0].cgs,bRead,file_version);
2274 if (bRead && gmx_debug_at) {
2275 pr_block(debug,0,"cgs",&mtop->moltype[0].cgs,TRUE);
2277 do_block(fio, &mtop->mols,bRead,file_version);
2278 /* Add the posres coordinates to the molblock */
2279 add_posres_molblock(mtop);
2282 if (file_version >= 57) {
2283 mtop->mols = mtop_mols(mtop);
2286 pr_block(debug,0,"mols",&mtop->mols,TRUE);
2290 if (file_version < 51) {
2291 /* Here used to be the shake blocks */
2292 do_blocka(fio, &dumb,bRead,file_version);
2300 close_symtab(&(mtop->symtab));
2304 /* If TopOnlyOK is TRUE then we can read even future versions
2305 * of tpx files, provided the file_generation hasn't changed.
2306 * If it is FALSE, we need the inputrecord too, and bail out
2307 * if the file is newer than the program.
2309 * The version and generation if the topology (see top of this file)
2310 * are returned in the two last arguments.
2312 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2314 static void do_tpxheader(t_fileio *fio,gmx_bool bRead,t_tpxheader *tpx,
2315 gmx_bool TopOnlyOK, int *file_version,
2316 int *file_generation)
2319 char file_tag[STRLEN];
2326 gmx_fio_checktype(fio);
2327 gmx_fio_setdebug(fio,bDebugMode());
2329 /* NEW! XDR tpb file */
2330 precision = sizeof(real);
2332 gmx_fio_do_string(fio,buf);
2333 if (strncmp(buf,"VERSION",7))
2334 gmx_fatal(FARGS,"Can not read file %s,\n"
2335 " this file is from a Gromacs version which is older than 2.0\n"
2336 " Make a new one with grompp or use a gro or pdb file, if possible",
2337 gmx_fio_getname(fio));
2338 gmx_fio_do_int(fio,precision);
2339 bDouble = (precision == sizeof(double));
2340 if ((precision != sizeof(float)) && !bDouble)
2341 gmx_fatal(FARGS,"Unknown precision in file %s: real is %d bytes "
2342 "instead of %d or %d",
2343 gmx_fio_getname(fio),precision,sizeof(float),sizeof(double));
2344 gmx_fio_setprecision(fio,bDouble);
2345 fprintf(stderr,"Reading file %s, %s (%s precision)\n",
2346 gmx_fio_getname(fio),buf,bDouble ? "double" : "single");
2349 gmx_fio_write_string(fio,GromacsVersion());
2350 bDouble = (precision == sizeof(double));
2351 gmx_fio_setprecision(fio,bDouble);
2352 gmx_fio_do_int(fio,precision);
2354 sprintf(file_tag,"%s",tpx_tag);
2355 fgen = tpx_generation;
2358 /* Check versions! */
2359 gmx_fio_do_int(fio,fver);
2363 gmx_fio_do_string(fio,file_tag);
2369 /* Versions before 77 don't have the tag, set it to release */
2370 sprintf(file_tag,"%s",TPX_TAG_RELEASE);
2373 if (strcmp(file_tag,tpx_tag) != 0)
2375 fprintf(stderr,"Note: file tpx tag '%s', software tpx tag '%s'\n",
2378 /* We only support reading tpx files with the same tag as the code
2379 * or tpx files with the release tag and with lower version number.
2381 if (!(strcmp(file_tag,TPX_TAG_RELEASE) == 0 && fver < tpx_version))
2383 gmx_fatal(FARGS,"tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2384 gmx_fio_getname(fio),fver,file_tag,
2385 tpx_version,tpx_tag);
2392 gmx_fio_do_int(fio,fgen);
2399 if (file_version != NULL)
2401 *file_version = fver;
2403 if (file_generation != NULL)
2405 *file_generation = fgen;
2409 if ((fver <= tpx_incompatible_version) ||
2410 ((fver > tpx_version) && !TopOnlyOK) ||
2411 (fgen > tpx_generation))
2412 gmx_fatal(FARGS,"reading tpx file (%s) version %d with version %d program",
2413 gmx_fio_getname(fio),fver,tpx_version);
2415 do_section(fio,eitemHEADER,bRead);
2416 gmx_fio_do_int(fio,tpx->natoms);
2418 gmx_fio_do_int(fio,tpx->ngtc);
2422 gmx_fio_do_int(fio,idum);
2423 gmx_fio_do_real(fio,rdum);
2425 /*a better decision will eventually (5.0 or later) need to be made
2426 on how to treat the alchemical state of the system, which can now
2427 vary through a simulation, and cannot be completely described
2428 though a single lambda variable, or even a single state
2429 index. Eventually, should probably be a vector. MRS*/
2432 gmx_fio_do_int(fio,tpx->fep_state);
2434 gmx_fio_do_real(fio,tpx->lambda);
2435 gmx_fio_do_int(fio,tpx->bIr);
2436 gmx_fio_do_int(fio,tpx->bTop);
2437 gmx_fio_do_int(fio,tpx->bX);
2438 gmx_fio_do_int(fio,tpx->bV);
2439 gmx_fio_do_int(fio,tpx->bF);
2440 gmx_fio_do_int(fio,tpx->bBox);
2442 if((fgen > tpx_generation)) {
2443 /* This can only happen if TopOnlyOK=TRUE */
2448 static int do_tpx(t_fileio *fio, gmx_bool bRead,
2449 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop,
2450 gmx_bool bXVallocated)
2455 gmx_bool TopOnlyOK,bDum=TRUE;
2456 int file_version,file_generation;
2460 gmx_bool bPeriodicMols;
2463 tpx.natoms = state->natoms;
2464 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
2465 tpx.fep_state = state->fep_state;
2466 tpx.lambda = state->lambda[efptFEP];
2467 tpx.bIr = (ir != NULL);
2468 tpx.bTop = (mtop != NULL);
2469 tpx.bX = (state->x != NULL);
2470 tpx.bV = (state->v != NULL);
2471 tpx.bF = (f != NULL);
2475 TopOnlyOK = (ir==NULL);
2477 do_tpxheader(fio,bRead,&tpx,TopOnlyOK,&file_version,&file_generation);
2481 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
2482 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
2486 init_state(state,0,tpx.ngtc,0,0,0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
2487 state->natoms = tpx.natoms;
2488 state->nalloc = tpx.natoms;
2492 init_state(state,tpx.natoms,tpx.ngtc,0,0,0); /* nose-hoover chains */
2496 #define do_test(fio,b,p) if (bRead && (p!=NULL) && !b) gmx_fatal(FARGS,"No %s in %s",#p,gmx_fio_getname(fio))
2498 do_test(fio,tpx.bBox,state->box);
2499 do_section(fio,eitemBOX,bRead);
2501 gmx_fio_ndo_rvec(fio,state->box,DIM);
2502 if (file_version >= 51) {
2503 gmx_fio_ndo_rvec(fio,state->box_rel,DIM);
2505 /* We initialize box_rel after reading the inputrec */
2506 clear_mat(state->box_rel);
2508 if (file_version >= 28) {
2509 gmx_fio_ndo_rvec(fio,state->boxv,DIM);
2510 if (file_version < 56) {
2512 gmx_fio_ndo_rvec(fio,mdum,DIM);
2517 if (state->ngtc > 0 && file_version >= 28) {
2519 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
2520 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
2521 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
2522 snew(dumv,state->ngtc);
2523 if (file_version < 69) {
2524 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2526 /* These used to be the Berendsen tcoupl_lambda's */
2527 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2531 /* Prior to tpx version 26, the inputrec was here.
2532 * I moved it to enable partial forward-compatibility
2533 * for analysis/viewer programs.
2535 if(file_version<26) {
2536 do_test(fio,tpx.bIr,ir);
2537 do_section(fio,eitemIR,bRead);
2540 do_inputrec(fio, ir,bRead,file_version,
2541 mtop ? &mtop->ffparams.fudgeQQ : NULL);
2543 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2546 do_inputrec(fio, &dum_ir,bRead,file_version,
2547 mtop ? &mtop->ffparams.fudgeQQ :NULL);
2549 pr_inputrec(debug,0,"inputrec",&dum_ir,FALSE);
2550 done_inputrec(&dum_ir);
2556 do_test(fio,tpx.bTop,mtop);
2557 do_section(fio,eitemTOP,bRead);
2560 do_mtop(fio,mtop,bRead, file_version);
2562 do_mtop(fio,&dum_top,bRead,file_version);
2563 done_mtop(&dum_top,TRUE);
2566 do_test(fio,tpx.bX,state->x);
2567 do_section(fio,eitemX,bRead);
2570 state->flags |= (1<<estX);
2572 gmx_fio_ndo_rvec(fio,state->x,state->natoms);
2575 do_test(fio,tpx.bV,state->v);
2576 do_section(fio,eitemV,bRead);
2579 state->flags |= (1<<estV);
2581 gmx_fio_ndo_rvec(fio,state->v,state->natoms);
2584 do_test(fio,tpx.bF,f);
2585 do_section(fio,eitemF,bRead);
2586 if (tpx.bF) gmx_fio_ndo_rvec(fio,f,state->natoms);
2588 /* Starting with tpx version 26, we have the inputrec
2589 * at the end of the file, so we can ignore it
2590 * if the file is never than the software (but still the
2591 * same generation - see comments at the top of this file.
2596 bPeriodicMols = FALSE;
2597 if (file_version >= 26) {
2598 do_test(fio,tpx.bIr,ir);
2599 do_section(fio,eitemIR,bRead);
2601 if (file_version >= 53) {
2602 /* Removed the pbc info from do_inputrec, since we always want it */
2605 bPeriodicMols = ir->bPeriodicMols;
2607 gmx_fio_do_int(fio,ePBC);
2608 gmx_fio_do_gmx_bool(fio,bPeriodicMols);
2610 if (file_generation <= tpx_generation && ir) {
2611 do_inputrec(fio, ir,bRead,file_version,mtop ? &mtop->ffparams.fudgeQQ : NULL);
2613 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2614 if (file_version < 51)
2615 set_box_rel(ir,state);
2616 if (file_version < 53) {
2618 bPeriodicMols = ir->bPeriodicMols;
2621 if (bRead && ir && file_version >= 53) {
2622 /* We need to do this after do_inputrec, since that initializes ir */
2624 ir->bPeriodicMols = bPeriodicMols;
2633 if (state->ngtc == 0)
2635 /* Reading old version without tcoupl state data: set it */
2636 init_gtc_state(state,ir->opts.ngtc,0,ir->opts.nhchainlength);
2638 if (tpx.bTop && mtop)
2640 if (file_version < 57)
2642 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
2644 ir->eDisre = edrSimple;
2648 ir->eDisre = edrNone;
2651 set_disres_npair(mtop);
2655 if (tpx.bTop && mtop)
2657 gmx_mtop_finalize(mtop);
2660 if (file_version >= 57)
2664 env = getenv("GMX_NOCHARGEGROUPS");
2667 sscanf(env,"%d",&ienv);
2668 fprintf(stderr,"\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
2673 "Will make single atomic charge groups in non-solvent%s\n",
2674 ienv > 1 ? " and solvent" : "");
2675 gmx_mtop_make_atomic_charge_groups(mtop,ienv==1);
2677 fprintf(stderr,"\n");
2685 /************************************************************
2687 * The following routines are the exported ones
2689 ************************************************************/
2691 t_fileio *open_tpx(const char *fn,const char *mode)
2693 return gmx_fio_open(fn,mode);
2696 void close_tpx(t_fileio *fio)
2701 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
2702 int *file_version, int *file_generation)
2706 fio = open_tpx(fn,"r");
2707 do_tpxheader(fio,TRUE,tpx,TopOnlyOK,file_version,file_generation);
2711 void write_tpx_state(const char *fn,
2712 t_inputrec *ir,t_state *state,gmx_mtop_t *mtop)
2716 fio = open_tpx(fn,"w");
2717 do_tpx(fio,FALSE,ir,state,NULL,mtop,FALSE);
2721 void read_tpx_state(const char *fn,
2722 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop)
2726 fio = open_tpx(fn,"r");
2727 do_tpx(fio,TRUE,ir,state,f,mtop,FALSE);
2731 int read_tpx(const char *fn,
2732 t_inputrec *ir, matrix box,int *natoms,
2733 rvec *x,rvec *v,rvec *f,gmx_mtop_t *mtop)
2741 fio = open_tpx(fn,"r");
2742 ePBC = do_tpx(fio,TRUE,ir,&state,f,mtop,TRUE);
2744 *natoms = state.natoms;
2746 copy_mat(state.box,box);
2754 int read_tpx_top(const char *fn,
2755 t_inputrec *ir, matrix box,int *natoms,
2756 rvec *x,rvec *v,rvec *f,t_topology *top)
2762 ePBC = read_tpx(fn,ir,box,natoms,x,v,f,&mtop);
2764 *top = gmx_mtop_t_to_t_topology(&mtop);
2769 gmx_bool fn2bTPX(const char *file)
2771 switch (fn2ftp(file)) {
2781 gmx_bool read_tps_conf(const char *infile,char *title,t_topology *top,int *ePBC,
2782 rvec **x,rvec **v,matrix box,gmx_bool bMass)
2785 int natoms,i,version,generation;
2786 gmx_bool bTop,bXNULL=FALSE;
2788 t_topology *topconv;
2791 bTop = fn2bTPX(infile);
2794 read_tpxheader(infile,&header,TRUE,&version,&generation);
2796 snew(*x,header.natoms);
2798 snew(*v,header.natoms);
2800 *ePBC = read_tpx(infile,NULL,box,&natoms,
2801 (x==NULL) ? NULL : *x,(v==NULL) ? NULL : *v,NULL,mtop);
2802 *top = gmx_mtop_t_to_t_topology(mtop);
2804 strcpy(title,*top->name);
2805 tpx_make_chain_identifiers(&top->atoms,&top->mols);
2808 get_stx_coordnum(infile,&natoms);
2809 init_t_atoms(&top->atoms,natoms,(fn2ftp(infile) == efPDB));
2818 read_stx_conf(infile,title,&top->atoms,*x,(v==NULL) ? NULL : *v,ePBC,box);
2825 aps = gmx_atomprop_init();
2826 for(i=0; (i<natoms); i++)
2827 if (!gmx_atomprop_query(aps,epropMass,
2828 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2829 *top->atoms.atomname[i],
2830 &(top->atoms.atom[i].m))) {
2832 fprintf(debug,"Can not find mass for atom %s %d %s, setting to 1\n",
2833 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2834 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
2835 *top->atoms.atomname[i]);
2837 gmx_atomprop_destroy(aps);
2839 top->idef.ntypes=-1;