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 = 92;
78 /* This number should only be increased when you edit the TOPOLOGY section
79 * or the HEADER of the tpx format.
80 * This way we can maintain forward compatibility too for all analysis tools
81 * and/or external programs that only need to know the atom/residue names,
82 * charges, and bond connectivity.
84 * It first appeared in tpx version 26, when I also moved the inputrecord
85 * to the end of the tpx file, so we can just skip it if we only
88 static const int tpx_generation = 25;
90 /* This number should be the most recent backwards incompatible version
91 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
93 static const int tpx_incompatible_version = 9;
97 /* Struct used to maintain tpx compatibility when function types are added */
99 int fvnr; /* file version number in which the function type first appeared */
100 int ftype; /* function type */
104 *The entries should be ordered in:
105 * 1. ascending file version number
106 * 2. ascending function type number
108 /*static const t_ftupd ftupd[] = {
109 { 20, F_CUBICBONDS },
113 { 22, F_DISRESVIOL },
119 { 26, F_DIHRESVIOL },
120 { 30, F_CROSS_BOND_BONDS },
121 { 30, F_CROSS_BOND_ANGLES },
122 { 30, F_UREY_BRADLEY },
123 { 30, F_POLARIZATION },
127 *The entries should be ordered in:
128 * 1. ascending function type number
129 * 2. ascending file version number
131 /* question; what is the purpose of the commented code above? */
132 static const t_ftupd ftupd[] = {
133 { 20, F_CUBICBONDS },
138 { 43, F_TABBONDSNC },
139 { 70, F_RESTRBONDS },
140 { 76, F_LINEAR_ANGLES },
141 { 30, F_CROSS_BOND_BONDS },
142 { 30, F_CROSS_BOND_ANGLES },
143 { 30, F_UREY_BRADLEY },
144 { 34, F_QUARTIC_ANGLES },
154 { 72, F_NPSOLVATION },
156 { 41, F_LJC_PAIRS_NB },
159 { 32, F_COUL_RECIP },
161 { 30, F_POLARIZATION },
164 { 22, F_DISRESVIOL },
168 { 26, F_DIHRESVIOL },
173 { 46, F_ECONSERVED },
174 { 69, F_VTEMP_NOLONGERUSED},
177 { 76, F_ANHARM_POL },
180 { 79, F_DVDL_BONDED, },
181 { 79, F_DVDL_RESTRAINT },
182 { 79, F_DVDL_TEMPERATURE },
185 #define NFTUPD asize(ftupd)
187 /* Needed for backward compatibility */
190 static void _do_section(t_fileio *fio,int key,gmx_bool bRead,const char *src,
196 if (gmx_fio_getftp(fio) == efTPA) {
198 gmx_fio_write_string(fio,itemstr[key]);
199 bDbg = gmx_fio_getdebug(fio);
200 gmx_fio_setdebug(fio,FALSE);
201 gmx_fio_write_string(fio,comment_str[key]);
202 gmx_fio_setdebug(fio,bDbg);
205 if (gmx_fio_getdebug(fio))
206 fprintf(stderr,"Looking for section %s (%s, %d)",
207 itemstr[key],src,line);
210 gmx_fio_do_string(fio,buf);
211 } while ((gmx_strcasecmp(buf,itemstr[key]) != 0));
213 if (gmx_strcasecmp(buf,itemstr[key]) != 0)
214 gmx_fatal(FARGS,"\nCould not find section heading %s",itemstr[key]);
215 else if (gmx_fio_getdebug(fio))
216 fprintf(stderr," and found it\n");
221 #define do_section(fio,key,bRead) _do_section(fio,key,bRead,__FILE__,__LINE__)
223 /**************************************************************
225 * Now the higer level routines that do io of the structures and arrays
227 **************************************************************/
228 static void do_pullgrp(t_fileio *fio, t_pullgrp *pgrp, gmx_bool bRead,
234 gmx_fio_do_int(fio,pgrp->nat);
236 snew(pgrp->ind,pgrp->nat);
237 bDum=gmx_fio_ndo_int(fio,pgrp->ind,pgrp->nat);
238 gmx_fio_do_int(fio,pgrp->nweight);
240 snew(pgrp->weight,pgrp->nweight);
241 bDum=gmx_fio_ndo_real(fio,pgrp->weight,pgrp->nweight);
242 gmx_fio_do_int(fio,pgrp->pbcatom);
243 gmx_fio_do_rvec(fio,pgrp->vec);
244 gmx_fio_do_rvec(fio,pgrp->init);
245 gmx_fio_do_real(fio,pgrp->rate);
246 gmx_fio_do_real(fio,pgrp->k);
247 if (file_version >= 56) {
248 gmx_fio_do_real(fio,pgrp->kB);
254 static void do_expandedvals(t_fileio *fio,t_expanded *expand,t_lambda *fepvals, gmx_bool bRead, int file_version)
256 /* i is used in the ndo_double macro*/
261 int n_lambda=fepvals->n_lambda;
263 /* reset the lambda calculation window */
264 fepvals->lambda_start_n = 0;
265 fepvals->lambda_stop_n = n_lambda;
266 if (file_version >= 79)
272 snew(expand->init_lambda_weights,n_lambda);
274 bDum=gmx_fio_ndo_real(fio,expand->init_lambda_weights,n_lambda);
275 gmx_fio_do_gmx_bool(fio,expand->bInit_weights);
278 gmx_fio_do_int(fio,expand->nstexpanded);
279 gmx_fio_do_int(fio,expand->elmcmove);
280 gmx_fio_do_int(fio,expand->elamstats);
281 gmx_fio_do_int(fio,expand->lmc_repeats);
282 gmx_fio_do_int(fio,expand->gibbsdeltalam);
283 gmx_fio_do_int(fio,expand->lmc_forced_nstart);
284 gmx_fio_do_int(fio,expand->lmc_seed);
285 gmx_fio_do_real(fio,expand->mc_temp);
286 gmx_fio_do_int(fio,expand->bSymmetrizedTMatrix);
287 gmx_fio_do_int(fio,expand->nstTij);
288 gmx_fio_do_int(fio,expand->minvarmin);
289 gmx_fio_do_int(fio,expand->c_range);
290 gmx_fio_do_real(fio,expand->wl_scale);
291 gmx_fio_do_real(fio,expand->wl_ratio);
292 gmx_fio_do_real(fio,expand->init_wl_delta);
293 gmx_fio_do_gmx_bool(fio,expand->bWLoneovert);
294 gmx_fio_do_int(fio,expand->elmceq);
295 gmx_fio_do_int(fio,expand->equil_steps);
296 gmx_fio_do_int(fio,expand->equil_samples);
297 gmx_fio_do_int(fio,expand->equil_n_at_lam);
298 gmx_fio_do_real(fio,expand->equil_wl_delta);
299 gmx_fio_do_real(fio,expand->equil_ratio);
303 static void do_simtempvals(t_fileio *fio,t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
308 if (file_version >= 79)
310 gmx_fio_do_int(fio,simtemp->eSimTempScale);
311 gmx_fio_do_real(fio,simtemp->simtemp_high);
312 gmx_fio_do_real(fio,simtemp->simtemp_low);
317 snew(simtemp->temperatures,n_lambda);
319 bDum=gmx_fio_ndo_real(fio,simtemp->temperatures,n_lambda);
324 static void do_fepvals(t_fileio *fio,t_lambda *fepvals,gmx_bool bRead, int file_version)
326 /* i is defined in the ndo_double macro; use g to iterate. */
332 /* free energy values */
334 if (file_version >= 79)
336 gmx_fio_do_int(fio,fepvals->init_fep_state);
337 gmx_fio_do_double(fio,fepvals->init_lambda);
338 gmx_fio_do_double(fio,fepvals->delta_lambda);
340 else if (file_version >= 59) {
341 gmx_fio_do_double(fio,fepvals->init_lambda);
342 gmx_fio_do_double(fio,fepvals->delta_lambda);
344 gmx_fio_do_real(fio,rdum);
345 fepvals->init_lambda = rdum;
346 gmx_fio_do_real(fio,rdum);
347 fepvals->delta_lambda = rdum;
349 if (file_version >= 79)
351 gmx_fio_do_int(fio,fepvals->n_lambda);
354 snew(fepvals->all_lambda,efptNR);
356 for (g=0;g<efptNR;g++)
358 if (fepvals->n_lambda > 0) {
361 snew(fepvals->all_lambda[g],fepvals->n_lambda);
363 bDum=gmx_fio_ndo_double(fio,fepvals->all_lambda[g],fepvals->n_lambda);
364 bDum=gmx_fio_ndo_int(fio,fepvals->separate_dvdl,efptNR);
366 else if (fepvals->init_lambda >= 0)
368 fepvals->separate_dvdl[efptFEP] = TRUE;
372 else if (file_version >= 64)
374 gmx_fio_do_int(fio,fepvals->n_lambda);
379 snew(fepvals->all_lambda,efptNR);
380 /* still allocate the all_lambda array's contents. */
381 for(g=0;g<efptNR;g++)
383 if (fepvals->n_lambda > 0) {
384 snew(fepvals->all_lambda[g],fepvals->n_lambda);
388 bDum=gmx_fio_ndo_double(fio,fepvals->all_lambda[efptFEP],
390 if (fepvals->init_lambda >= 0)
394 fepvals->separate_dvdl[efptFEP] = TRUE;
398 /* copy the contents of the efptFEP lambda component to all
399 the other components */
400 for(g=0;g<efptNR;g++)
402 for(h=0;h<fepvals->n_lambda;h++)
406 fepvals->all_lambda[g][h] =
407 fepvals->all_lambda[efptFEP][h];
416 fepvals->n_lambda = 0;
417 fepvals->all_lambda = NULL;
418 if (fepvals->init_lambda >= 0)
420 fepvals->separate_dvdl[efptFEP] = TRUE;
423 if (file_version >= 13)
425 gmx_fio_do_real(fio,fepvals->sc_alpha);
429 fepvals->sc_alpha = 0;
431 if (file_version >= 38)
433 gmx_fio_do_int(fio,fepvals->sc_power);
437 fepvals->sc_power = 2;
439 if (file_version >= 79)
441 gmx_fio_do_real(fio,fepvals->sc_r_power);
445 fepvals->sc_r_power = 6.0;
447 if (file_version >= 15)
449 gmx_fio_do_real(fio,fepvals->sc_sigma);
453 fepvals->sc_sigma = 0.3;
457 if (file_version >= 71)
459 fepvals->sc_sigma_min = fepvals->sc_sigma;
463 fepvals->sc_sigma_min = 0;
466 if (file_version >= 79)
468 gmx_fio_do_int(fio,fepvals->bScCoul);
472 fepvals->bScCoul = TRUE;
474 if (file_version >= 64) {
475 gmx_fio_do_int(fio,fepvals->nstdhdl);
477 fepvals->nstdhdl = 1;
480 if (file_version >= 73)
482 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
483 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
487 fepvals->separate_dhdl_file = esepdhdlfileYES;
488 fepvals->dhdl_derivatives = edhdlderivativesYES;
490 if (file_version >= 71)
492 gmx_fio_do_int(fio,fepvals->dh_hist_size);
493 gmx_fio_do_double(fio,fepvals->dh_hist_spacing);
497 fepvals->dh_hist_size = 0;
498 fepvals->dh_hist_spacing = 0.1;
500 if (file_version >= 79)
502 gmx_fio_do_int(fio,fepvals->bPrintEnergy);
506 fepvals->bPrintEnergy = FALSE;
509 /* handle lambda_neighbors */
510 if ((file_version >= 83 && file_version < 90) || file_version >= 92 )
512 gmx_fio_do_int(fio,fepvals->lambda_neighbors);
513 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state>=0) &&
514 (fepvals->init_lambda < 0) )
516 fepvals->lambda_start_n = (fepvals->init_fep_state -
517 fepvals->lambda_neighbors);
518 fepvals->lambda_stop_n = (fepvals->init_fep_state +
519 fepvals->lambda_neighbors + 1);
520 if (fepvals->lambda_start_n < 0)
522 fepvals->lambda_start_n = 0;;
524 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
526 fepvals->lambda_stop_n = fepvals->n_lambda;
531 fepvals->lambda_start_n = 0;
532 fepvals->lambda_stop_n = fepvals->n_lambda;
537 fepvals->lambda_start_n = 0;
538 fepvals->lambda_stop_n = fepvals->n_lambda;
542 static void do_pull(t_fileio *fio, t_pull *pull,gmx_bool bRead, int file_version)
546 gmx_fio_do_int(fio,pull->ngrp);
547 gmx_fio_do_int(fio,pull->eGeom);
548 gmx_fio_do_ivec(fio,pull->dim);
549 gmx_fio_do_real(fio,pull->cyl_r1);
550 gmx_fio_do_real(fio,pull->cyl_r0);
551 gmx_fio_do_real(fio,pull->constr_tol);
552 gmx_fio_do_int(fio,pull->nstxout);
553 gmx_fio_do_int(fio,pull->nstfout);
555 snew(pull->grp,pull->ngrp+1);
556 for(g=0; g<pull->ngrp+1; g++)
557 do_pullgrp(fio,&pull->grp[g],bRead,file_version);
561 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg,gmx_bool bRead, int file_version)
566 gmx_fio_do_int(fio,rotg->eType);
567 gmx_fio_do_int(fio,rotg->bMassW);
568 gmx_fio_do_int(fio,rotg->nat);
570 snew(rotg->ind,rotg->nat);
571 gmx_fio_ndo_int(fio,rotg->ind,rotg->nat);
573 snew(rotg->x_ref,rotg->nat);
574 gmx_fio_ndo_rvec(fio,rotg->x_ref,rotg->nat);
575 gmx_fio_do_rvec(fio,rotg->vec);
576 gmx_fio_do_rvec(fio,rotg->pivot);
577 gmx_fio_do_real(fio,rotg->rate);
578 gmx_fio_do_real(fio,rotg->k);
579 gmx_fio_do_real(fio,rotg->slab_dist);
580 gmx_fio_do_real(fio,rotg->min_gaussian);
581 gmx_fio_do_real(fio,rotg->eps);
582 gmx_fio_do_int(fio,rotg->eFittype);
583 gmx_fio_do_int(fio,rotg->PotAngle_nstep);
584 gmx_fio_do_real(fio,rotg->PotAngle_step);
587 static void do_rot(t_fileio *fio, t_rot *rot,gmx_bool bRead, int file_version)
591 gmx_fio_do_int(fio,rot->ngrp);
592 gmx_fio_do_int(fio,rot->nstrout);
593 gmx_fio_do_int(fio,rot->nstsout);
595 snew(rot->grp,rot->ngrp);
596 for(g=0; g<rot->ngrp; g++)
597 do_rotgrp(fio, &rot->grp[g],bRead,file_version);
601 static void do_inputrec(t_fileio *fio, t_inputrec *ir,gmx_bool bRead,
602 int file_version, real *fudgeQQ)
604 int i,j,k,*tmp,idum=0;
609 real zerotemptime,finish_t,init_temp,finish_temp;
611 if (file_version != tpx_version)
613 /* Give a warning about features that are not accessible */
614 fprintf(stderr,"Note: file tpx version %d, software tpx version %d\n",
615 file_version,tpx_version);
623 if (file_version == 0)
628 /* Basic inputrec stuff */
629 gmx_fio_do_int(fio,ir->eI);
630 if (file_version >= 62) {
631 gmx_fio_do_gmx_large_int(fio, ir->nsteps);
633 gmx_fio_do_int(fio,idum);
636 if(file_version > 25) {
637 if (file_version >= 62) {
638 gmx_fio_do_gmx_large_int(fio, ir->init_step);
640 gmx_fio_do_int(fio,idum);
641 ir->init_step = idum;
647 if(file_version >= 58)
648 gmx_fio_do_int(fio,ir->simulation_part);
650 ir->simulation_part=1;
652 if (file_version >= 67) {
653 gmx_fio_do_int(fio,ir->nstcalcenergy);
655 ir->nstcalcenergy = 1;
657 if (file_version < 53) {
658 /* The pbc info has been moved out of do_inputrec,
659 * since we always want it, also without reading the inputrec.
661 gmx_fio_do_int(fio,ir->ePBC);
662 if ((file_version <= 15) && (ir->ePBC == 2))
664 if (file_version >= 45) {
665 gmx_fio_do_int(fio,ir->bPeriodicMols);
669 ir->bPeriodicMols = TRUE;
671 ir->bPeriodicMols = FALSE;
675 if (file_version >= 81)
677 gmx_fio_do_int(fio,ir->cutoff_scheme);
681 ir->cutoff_scheme = ecutsGROUP;
683 gmx_fio_do_int(fio,ir->ns_type);
684 gmx_fio_do_int(fio,ir->nstlist);
685 gmx_fio_do_int(fio,ir->ndelta);
686 if (file_version < 41) {
687 gmx_fio_do_int(fio,idum);
688 gmx_fio_do_int(fio,idum);
690 if (file_version >= 45)
691 gmx_fio_do_real(fio,ir->rtpi);
694 gmx_fio_do_int(fio,ir->nstcomm);
695 if (file_version > 34)
696 gmx_fio_do_int(fio,ir->comm_mode);
697 else if (ir->nstcomm < 0)
698 ir->comm_mode = ecmANGULAR;
700 ir->comm_mode = ecmLINEAR;
701 ir->nstcomm = abs(ir->nstcomm);
703 if(file_version > 25)
704 gmx_fio_do_int(fio,ir->nstcheckpoint);
708 gmx_fio_do_int(fio,ir->nstcgsteep);
711 gmx_fio_do_int(fio,ir->nbfgscorr);
715 gmx_fio_do_int(fio,ir->nstlog);
716 gmx_fio_do_int(fio,ir->nstxout);
717 gmx_fio_do_int(fio,ir->nstvout);
718 gmx_fio_do_int(fio,ir->nstfout);
719 gmx_fio_do_int(fio,ir->nstenergy);
720 gmx_fio_do_int(fio,ir->nstxtcout);
721 if (file_version >= 59) {
722 gmx_fio_do_double(fio,ir->init_t);
723 gmx_fio_do_double(fio,ir->delta_t);
725 gmx_fio_do_real(fio,rdum);
727 gmx_fio_do_real(fio,rdum);
730 gmx_fio_do_real(fio,ir->xtcprec);
731 if (file_version < 19) {
732 gmx_fio_do_int(fio,idum);
733 gmx_fio_do_int(fio,idum);
735 if(file_version < 18)
736 gmx_fio_do_int(fio,idum);
737 if (file_version >= 81) {
738 gmx_fio_do_real(fio,ir->verletbuf_drift);
740 ir->verletbuf_drift = 0;
742 gmx_fio_do_real(fio,ir->rlist);
743 if (file_version >= 67) {
744 gmx_fio_do_real(fio,ir->rlistlong);
746 if(file_version >= 82 && file_version != 90)
748 gmx_fio_do_int(fio,ir->nstcalclr);
752 /* Calculate at NS steps */
753 ir->nstcalclr = ir->nstlist;
755 gmx_fio_do_int(fio,ir->coulombtype);
756 if (file_version < 32 && ir->coulombtype == eelRF)
757 ir->coulombtype = eelRF_NEC;
758 if (file_version >= 81)
760 gmx_fio_do_int(fio,ir->coulomb_modifier);
764 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
766 gmx_fio_do_real(fio,ir->rcoulomb_switch);
767 gmx_fio_do_real(fio,ir->rcoulomb);
768 gmx_fio_do_int(fio,ir->vdwtype);
769 if (file_version >= 81)
771 gmx_fio_do_int(fio,ir->vdw_modifier);
775 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
777 gmx_fio_do_real(fio,ir->rvdw_switch);
778 gmx_fio_do_real(fio,ir->rvdw);
779 if (file_version < 67) {
780 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
782 gmx_fio_do_int(fio,ir->eDispCorr);
783 gmx_fio_do_real(fio,ir->epsilon_r);
784 if (file_version >= 37) {
785 gmx_fio_do_real(fio,ir->epsilon_rf);
787 if (EEL_RF(ir->coulombtype)) {
788 ir->epsilon_rf = ir->epsilon_r;
791 ir->epsilon_rf = 1.0;
794 if (file_version >= 29)
795 gmx_fio_do_real(fio,ir->tabext);
799 if(file_version > 25) {
800 gmx_fio_do_int(fio,ir->gb_algorithm);
801 gmx_fio_do_int(fio,ir->nstgbradii);
802 gmx_fio_do_real(fio,ir->rgbradii);
803 gmx_fio_do_real(fio,ir->gb_saltconc);
804 gmx_fio_do_int(fio,ir->implicit_solvent);
806 ir->gb_algorithm=egbSTILL;
810 ir->implicit_solvent=eisNO;
814 gmx_fio_do_real(fio,ir->gb_epsilon_solvent);
815 gmx_fio_do_real(fio,ir->gb_obc_alpha);
816 gmx_fio_do_real(fio,ir->gb_obc_beta);
817 gmx_fio_do_real(fio,ir->gb_obc_gamma);
820 gmx_fio_do_real(fio,ir->gb_dielectric_offset);
821 gmx_fio_do_int(fio,ir->sa_algorithm);
825 ir->gb_dielectric_offset = 0.009;
826 ir->sa_algorithm = esaAPPROX;
828 gmx_fio_do_real(fio,ir->sa_surface_tension);
830 /* Override sa_surface_tension if it is not changed in the mpd-file */
831 if(ir->sa_surface_tension<0)
833 if(ir->gb_algorithm==egbSTILL)
835 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
837 else if(ir->gb_algorithm==egbHCT || ir->gb_algorithm==egbOBC)
839 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
846 /* Better use sensible values than insane (0.0) ones... */
847 ir->gb_epsilon_solvent = 80;
848 ir->gb_obc_alpha = 1.0;
849 ir->gb_obc_beta = 0.8;
850 ir->gb_obc_gamma = 4.85;
851 ir->sa_surface_tension = 2.092;
855 if (file_version >= 81)
857 gmx_fio_do_real(fio,ir->fourier_spacing);
861 ir->fourier_spacing = 0.0;
863 gmx_fio_do_int(fio,ir->nkx);
864 gmx_fio_do_int(fio,ir->nky);
865 gmx_fio_do_int(fio,ir->nkz);
866 gmx_fio_do_int(fio,ir->pme_order);
867 gmx_fio_do_real(fio,ir->ewald_rtol);
869 if (file_version >=24)
870 gmx_fio_do_int(fio,ir->ewald_geometry);
872 ir->ewald_geometry=eewg3D;
874 if (file_version <=17) {
875 ir->epsilon_surface=0;
876 if (file_version==17)
877 gmx_fio_do_int(fio,idum);
880 gmx_fio_do_real(fio,ir->epsilon_surface);
882 gmx_fio_do_gmx_bool(fio,ir->bOptFFT);
884 gmx_fio_do_gmx_bool(fio,ir->bContinuation);
885 gmx_fio_do_int(fio,ir->etc);
886 /* before version 18, ir->etc was a gmx_bool (ir->btc),
887 * but the values 0 and 1 still mean no and
888 * berendsen temperature coupling, respectively.
890 if (file_version >= 79) {
891 gmx_fio_do_gmx_bool(fio,ir->bPrintNHChains);
893 if (file_version >= 71)
895 gmx_fio_do_int(fio,ir->nsttcouple);
899 ir->nsttcouple = ir->nstcalcenergy;
901 if (file_version <= 15)
903 gmx_fio_do_int(fio,idum);
905 if (file_version <=17)
907 gmx_fio_do_int(fio,ir->epct);
908 if (file_version <= 15)
912 ir->epct = epctSURFACETENSION;
914 gmx_fio_do_int(fio,idum);
917 /* we have removed the NO alternative at the beginning */
921 ir->epct=epctISOTROPIC;
925 ir->epc=epcBERENDSEN;
930 gmx_fio_do_int(fio,ir->epc);
931 gmx_fio_do_int(fio,ir->epct);
933 if (file_version >= 71)
935 gmx_fio_do_int(fio,ir->nstpcouple);
939 ir->nstpcouple = ir->nstcalcenergy;
941 gmx_fio_do_real(fio,ir->tau_p);
942 if (file_version <= 15) {
943 gmx_fio_do_rvec(fio,vdum);
944 clear_mat(ir->ref_p);
946 ir->ref_p[i][i] = vdum[i];
948 gmx_fio_do_rvec(fio,ir->ref_p[XX]);
949 gmx_fio_do_rvec(fio,ir->ref_p[YY]);
950 gmx_fio_do_rvec(fio,ir->ref_p[ZZ]);
952 if (file_version <= 15) {
953 gmx_fio_do_rvec(fio,vdum);
954 clear_mat(ir->compress);
956 ir->compress[i][i] = vdum[i];
959 gmx_fio_do_rvec(fio,ir->compress[XX]);
960 gmx_fio_do_rvec(fio,ir->compress[YY]);
961 gmx_fio_do_rvec(fio,ir->compress[ZZ]);
963 if (file_version >= 47) {
964 gmx_fio_do_int(fio,ir->refcoord_scaling);
965 gmx_fio_do_rvec(fio,ir->posres_com);
966 gmx_fio_do_rvec(fio,ir->posres_comB);
968 ir->refcoord_scaling = erscNO;
969 clear_rvec(ir->posres_com);
970 clear_rvec(ir->posres_comB);
972 if((file_version > 25) && (file_version < 79))
973 gmx_fio_do_int(fio,ir->andersen_seed);
976 if(file_version < 26) {
977 gmx_fio_do_gmx_bool(fio,bSimAnn);
978 gmx_fio_do_real(fio,zerotemptime);
981 if (file_version < 37)
982 gmx_fio_do_real(fio,rdum);
984 gmx_fio_do_real(fio,ir->shake_tol);
985 if (file_version < 54)
986 gmx_fio_do_real(fio,*fudgeQQ);
988 gmx_fio_do_int(fio,ir->efep);
989 if (file_version <= 14 && ir->efep != efepNO)
993 do_fepvals(fio,ir->fepvals,bRead,file_version);
995 if (file_version >= 79)
997 gmx_fio_do_gmx_bool(fio,ir->bSimTemp);
1000 ir->bSimTemp = TRUE;
1005 ir->bSimTemp = FALSE;
1009 do_simtempvals(fio,ir->simtempvals,ir->fepvals->n_lambda,bRead,file_version);
1012 if (file_version >= 79)
1014 gmx_fio_do_gmx_bool(fio,ir->bExpanded);
1017 ir->bExpanded = TRUE;
1021 ir->bExpanded = FALSE;
1026 do_expandedvals(fio,ir->expandedvals,ir->fepvals,bRead,file_version);
1028 if (file_version >= 57) {
1029 gmx_fio_do_int(fio,ir->eDisre);
1031 gmx_fio_do_int(fio,ir->eDisreWeighting);
1032 if (file_version < 22) {
1033 if (ir->eDisreWeighting == 0)
1034 ir->eDisreWeighting = edrwEqual;
1036 ir->eDisreWeighting = edrwConservative;
1038 gmx_fio_do_gmx_bool(fio,ir->bDisreMixed);
1039 gmx_fio_do_real(fio,ir->dr_fc);
1040 gmx_fio_do_real(fio,ir->dr_tau);
1041 gmx_fio_do_int(fio,ir->nstdisreout);
1042 if (file_version >= 22) {
1043 gmx_fio_do_real(fio,ir->orires_fc);
1044 gmx_fio_do_real(fio,ir->orires_tau);
1045 gmx_fio_do_int(fio,ir->nstorireout);
1049 ir->nstorireout = 0;
1051 if(file_version >= 26 && file_version < 79) {
1052 gmx_fio_do_real(fio,ir->dihre_fc);
1053 if (file_version < 56)
1055 gmx_fio_do_real(fio,rdum);
1056 gmx_fio_do_int(fio,idum);
1062 gmx_fio_do_real(fio,ir->em_stepsize);
1063 gmx_fio_do_real(fio,ir->em_tol);
1064 if (file_version >= 22)
1065 gmx_fio_do_gmx_bool(fio,ir->bShakeSOR);
1067 ir->bShakeSOR = TRUE;
1068 if (file_version >= 11)
1069 gmx_fio_do_int(fio,ir->niter);
1072 fprintf(stderr,"Note: niter not in run input file, setting it to %d\n",
1075 if (file_version >= 21)
1076 gmx_fio_do_real(fio,ir->fc_stepsize);
1078 ir->fc_stepsize = 0;
1079 gmx_fio_do_int(fio,ir->eConstrAlg);
1080 gmx_fio_do_int(fio,ir->nProjOrder);
1081 gmx_fio_do_real(fio,ir->LincsWarnAngle);
1082 if (file_version <= 14)
1083 gmx_fio_do_int(fio,idum);
1084 if (file_version >=26)
1085 gmx_fio_do_int(fio,ir->nLincsIter);
1088 fprintf(stderr,"Note: nLincsIter not in run input file, setting it to %d\n",
1091 if (file_version < 33)
1092 gmx_fio_do_real(fio,bd_temp);
1093 gmx_fio_do_real(fio,ir->bd_fric);
1094 gmx_fio_do_int(fio,ir->ld_seed);
1095 if (file_version >= 33) {
1096 for(i=0; i<DIM; i++)
1097 gmx_fio_do_rvec(fio,ir->deform[i]);
1099 for(i=0; i<DIM; i++)
1100 clear_rvec(ir->deform[i]);
1102 if (file_version >= 14)
1103 gmx_fio_do_real(fio,ir->cos_accel);
1106 gmx_fio_do_int(fio,ir->userint1);
1107 gmx_fio_do_int(fio,ir->userint2);
1108 gmx_fio_do_int(fio,ir->userint3);
1109 gmx_fio_do_int(fio,ir->userint4);
1110 gmx_fio_do_real(fio,ir->userreal1);
1111 gmx_fio_do_real(fio,ir->userreal2);
1112 gmx_fio_do_real(fio,ir->userreal3);
1113 gmx_fio_do_real(fio,ir->userreal4);
1116 if (file_version >= 77) {
1117 gmx_fio_do_gmx_bool(fio,ir->bAdress);
1119 if (bRead) snew(ir->adress, 1);
1120 gmx_fio_do_int(fio,ir->adress->type);
1121 gmx_fio_do_real(fio,ir->adress->const_wf);
1122 gmx_fio_do_real(fio,ir->adress->ex_width);
1123 gmx_fio_do_real(fio,ir->adress->hy_width);
1124 gmx_fio_do_int(fio,ir->adress->icor);
1125 gmx_fio_do_int(fio,ir->adress->site);
1126 gmx_fio_do_rvec(fio,ir->adress->refs);
1127 gmx_fio_do_int(fio,ir->adress->n_tf_grps);
1128 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1129 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1130 gmx_fio_do_int(fio,ir->adress->do_hybridpairs);
1132 if (bRead)snew(ir->adress->tf_table_index,ir->adress->n_tf_grps);
1133 if (ir->adress->n_tf_grps > 0) {
1134 bDum=gmx_fio_ndo_int(fio,ir->adress->tf_table_index,ir->adress->n_tf_grps);
1136 if (bRead)snew(ir->adress->group_explicit,ir->adress->n_energy_grps);
1137 if (ir->adress->n_energy_grps > 0) {
1138 bDum=gmx_fio_ndo_int(fio, ir->adress->group_explicit,ir->adress->n_energy_grps);
1142 ir->bAdress = FALSE;
1146 if (file_version >= 48) {
1147 gmx_fio_do_int(fio,ir->ePull);
1148 if (ir->ePull != epullNO) {
1151 do_pull(fio, ir->pull,bRead,file_version);
1154 ir->ePull = epullNO;
1157 /* Enforced rotation */
1158 if (file_version >= 74) {
1159 gmx_fio_do_int(fio,ir->bRot);
1160 if (ir->bRot == TRUE) {
1163 do_rot(fio, ir->rot,bRead,file_version);
1170 gmx_fio_do_int(fio,ir->opts.ngtc);
1171 if (file_version >= 69) {
1172 gmx_fio_do_int(fio,ir->opts.nhchainlength);
1174 ir->opts.nhchainlength = 1;
1176 gmx_fio_do_int(fio,ir->opts.ngacc);
1177 gmx_fio_do_int(fio,ir->opts.ngfrz);
1178 gmx_fio_do_int(fio,ir->opts.ngener);
1181 snew(ir->opts.nrdf, ir->opts.ngtc);
1182 snew(ir->opts.ref_t, ir->opts.ngtc);
1183 snew(ir->opts.annealing, ir->opts.ngtc);
1184 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1185 snew(ir->opts.anneal_time, ir->opts.ngtc);
1186 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1187 snew(ir->opts.tau_t, ir->opts.ngtc);
1188 snew(ir->opts.nFreeze,ir->opts.ngfrz);
1189 snew(ir->opts.acc, ir->opts.ngacc);
1190 snew(ir->opts.egp_flags,ir->opts.ngener*ir->opts.ngener);
1192 if (ir->opts.ngtc > 0) {
1193 if (bRead && file_version<13) {
1194 snew(tmp,ir->opts.ngtc);
1195 bDum=gmx_fio_ndo_int(fio,tmp, ir->opts.ngtc);
1196 for(i=0; i<ir->opts.ngtc; i++)
1197 ir->opts.nrdf[i] = tmp[i];
1200 bDum=gmx_fio_ndo_real(fio,ir->opts.nrdf, ir->opts.ngtc);
1202 bDum=gmx_fio_ndo_real(fio,ir->opts.ref_t,ir->opts.ngtc);
1203 bDum=gmx_fio_ndo_real(fio,ir->opts.tau_t,ir->opts.ngtc);
1204 if (file_version<33 && ir->eI==eiBD) {
1205 for(i=0; i<ir->opts.ngtc; i++)
1206 ir->opts.tau_t[i] = bd_temp;
1209 if (ir->opts.ngfrz > 0)
1210 bDum=gmx_fio_ndo_ivec(fio,ir->opts.nFreeze,ir->opts.ngfrz);
1211 if (ir->opts.ngacc > 0)
1212 gmx_fio_ndo_rvec(fio,ir->opts.acc,ir->opts.ngacc);
1213 if (file_version >= 12)
1214 bDum=gmx_fio_ndo_int(fio,ir->opts.egp_flags,
1215 ir->opts.ngener*ir->opts.ngener);
1217 if(bRead && file_version < 26) {
1218 for(i=0;i<ir->opts.ngtc;i++) {
1220 ir->opts.annealing[i] = eannSINGLE;
1221 ir->opts.anneal_npoints[i] = 2;
1222 snew(ir->opts.anneal_time[i],2);
1223 snew(ir->opts.anneal_temp[i],2);
1224 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1225 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1226 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1227 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1228 ir->opts.anneal_time[i][0] = ir->init_t;
1229 ir->opts.anneal_time[i][1] = finish_t;
1230 ir->opts.anneal_temp[i][0] = init_temp;
1231 ir->opts.anneal_temp[i][1] = finish_temp;
1233 ir->opts.annealing[i] = eannNO;
1234 ir->opts.anneal_npoints[i] = 0;
1238 /* file version 26 or later */
1239 /* First read the lists with annealing and npoints for each group */
1240 bDum=gmx_fio_ndo_int(fio,ir->opts.annealing,ir->opts.ngtc);
1241 bDum=gmx_fio_ndo_int(fio,ir->opts.anneal_npoints,ir->opts.ngtc);
1242 for(j=0;j<(ir->opts.ngtc);j++) {
1243 k=ir->opts.anneal_npoints[j];
1245 snew(ir->opts.anneal_time[j],k);
1246 snew(ir->opts.anneal_temp[j],k);
1248 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_time[j],k);
1249 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_temp[j],k);
1253 if (file_version >= 45) {
1254 gmx_fio_do_int(fio,ir->nwall);
1255 gmx_fio_do_int(fio,ir->wall_type);
1256 if (file_version >= 50)
1257 gmx_fio_do_real(fio,ir->wall_r_linpot);
1259 ir->wall_r_linpot = -1;
1260 gmx_fio_do_int(fio,ir->wall_atomtype[0]);
1261 gmx_fio_do_int(fio,ir->wall_atomtype[1]);
1262 gmx_fio_do_real(fio,ir->wall_density[0]);
1263 gmx_fio_do_real(fio,ir->wall_density[1]);
1264 gmx_fio_do_real(fio,ir->wall_ewald_zfac);
1268 ir->wall_atomtype[0] = -1;
1269 ir->wall_atomtype[1] = -1;
1270 ir->wall_density[0] = 0;
1271 ir->wall_density[1] = 0;
1272 ir->wall_ewald_zfac = 3;
1274 /* Cosine stuff for electric fields */
1275 for(j=0; (j<DIM); j++) {
1276 gmx_fio_do_int(fio,ir->ex[j].n);
1277 gmx_fio_do_int(fio,ir->et[j].n);
1279 snew(ir->ex[j].a, ir->ex[j].n);
1280 snew(ir->ex[j].phi,ir->ex[j].n);
1281 snew(ir->et[j].a, ir->et[j].n);
1282 snew(ir->et[j].phi,ir->et[j].n);
1284 bDum=gmx_fio_ndo_real(fio,ir->ex[j].a, ir->ex[j].n);
1285 bDum=gmx_fio_ndo_real(fio,ir->ex[j].phi,ir->ex[j].n);
1286 bDum=gmx_fio_ndo_real(fio,ir->et[j].a, ir->et[j].n);
1287 bDum=gmx_fio_ndo_real(fio,ir->et[j].phi,ir->et[j].n);
1291 if(file_version>=39){
1292 gmx_fio_do_gmx_bool(fio,ir->bQMMM);
1293 gmx_fio_do_int(fio,ir->QMMMscheme);
1294 gmx_fio_do_real(fio,ir->scalefactor);
1295 gmx_fio_do_int(fio,ir->opts.ngQM);
1297 snew(ir->opts.QMmethod, ir->opts.ngQM);
1298 snew(ir->opts.QMbasis, ir->opts.ngQM);
1299 snew(ir->opts.QMcharge, ir->opts.ngQM);
1300 snew(ir->opts.QMmult, ir->opts.ngQM);
1301 snew(ir->opts.bSH, ir->opts.ngQM);
1302 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1303 snew(ir->opts.CASelectrons,ir->opts.ngQM);
1304 snew(ir->opts.SAon, ir->opts.ngQM);
1305 snew(ir->opts.SAoff, ir->opts.ngQM);
1306 snew(ir->opts.SAsteps, ir->opts.ngQM);
1307 snew(ir->opts.bOPT, ir->opts.ngQM);
1308 snew(ir->opts.bTS, ir->opts.ngQM);
1310 if (ir->opts.ngQM > 0) {
1311 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmethod,ir->opts.ngQM);
1312 bDum=gmx_fio_ndo_int(fio,ir->opts.QMbasis,ir->opts.ngQM);
1313 bDum=gmx_fio_ndo_int(fio,ir->opts.QMcharge,ir->opts.ngQM);
1314 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmult,ir->opts.ngQM);
1315 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bSH,ir->opts.ngQM);
1316 bDum=gmx_fio_ndo_int(fio,ir->opts.CASorbitals,ir->opts.ngQM);
1317 bDum=gmx_fio_ndo_int(fio,ir->opts.CASelectrons,ir->opts.ngQM);
1318 bDum=gmx_fio_ndo_real(fio,ir->opts.SAon,ir->opts.ngQM);
1319 bDum=gmx_fio_ndo_real(fio,ir->opts.SAoff,ir->opts.ngQM);
1320 bDum=gmx_fio_ndo_int(fio,ir->opts.SAsteps,ir->opts.ngQM);
1321 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bOPT,ir->opts.ngQM);
1322 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bTS,ir->opts.ngQM);
1324 /* end of QMMM stuff */
1329 static void do_harm(t_fileio *fio, t_iparams *iparams,gmx_bool bRead)
1331 gmx_fio_do_real(fio,iparams->harmonic.rA);
1332 gmx_fio_do_real(fio,iparams->harmonic.krA);
1333 gmx_fio_do_real(fio,iparams->harmonic.rB);
1334 gmx_fio_do_real(fio,iparams->harmonic.krB);
1337 void do_iparams(t_fileio *fio, t_functype ftype,t_iparams *iparams,
1338 gmx_bool bRead, int file_version)
1345 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1353 do_harm(fio, iparams,bRead);
1354 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead) {
1355 /* Correct incorrect storage of parameters */
1356 iparams->pdihs.phiB = iparams->pdihs.phiA;
1357 iparams->pdihs.cpB = iparams->pdihs.cpA;
1360 case F_LINEAR_ANGLES:
1361 gmx_fio_do_real(fio,iparams->linangle.klinA);
1362 gmx_fio_do_real(fio,iparams->linangle.aA);
1363 gmx_fio_do_real(fio,iparams->linangle.klinB);
1364 gmx_fio_do_real(fio,iparams->linangle.aB);
1367 gmx_fio_do_real(fio,iparams->fene.bm);
1368 gmx_fio_do_real(fio,iparams->fene.kb);
1371 gmx_fio_do_real(fio,iparams->restraint.lowA);
1372 gmx_fio_do_real(fio,iparams->restraint.up1A);
1373 gmx_fio_do_real(fio,iparams->restraint.up2A);
1374 gmx_fio_do_real(fio,iparams->restraint.kA);
1375 gmx_fio_do_real(fio,iparams->restraint.lowB);
1376 gmx_fio_do_real(fio,iparams->restraint.up1B);
1377 gmx_fio_do_real(fio,iparams->restraint.up2B);
1378 gmx_fio_do_real(fio,iparams->restraint.kB);
1384 gmx_fio_do_real(fio,iparams->tab.kA);
1385 gmx_fio_do_int(fio,iparams->tab.table);
1386 gmx_fio_do_real(fio,iparams->tab.kB);
1388 case F_CROSS_BOND_BONDS:
1389 gmx_fio_do_real(fio,iparams->cross_bb.r1e);
1390 gmx_fio_do_real(fio,iparams->cross_bb.r2e);
1391 gmx_fio_do_real(fio,iparams->cross_bb.krr);
1393 case F_CROSS_BOND_ANGLES:
1394 gmx_fio_do_real(fio,iparams->cross_ba.r1e);
1395 gmx_fio_do_real(fio,iparams->cross_ba.r2e);
1396 gmx_fio_do_real(fio,iparams->cross_ba.r3e);
1397 gmx_fio_do_real(fio,iparams->cross_ba.krt);
1399 case F_UREY_BRADLEY:
1400 gmx_fio_do_real(fio,iparams->u_b.thetaA);
1401 gmx_fio_do_real(fio,iparams->u_b.kthetaA);
1402 gmx_fio_do_real(fio,iparams->u_b.r13A);
1403 gmx_fio_do_real(fio,iparams->u_b.kUBA);
1404 if (file_version >= 79) {
1405 gmx_fio_do_real(fio,iparams->u_b.thetaB);
1406 gmx_fio_do_real(fio,iparams->u_b.kthetaB);
1407 gmx_fio_do_real(fio,iparams->u_b.r13B);
1408 gmx_fio_do_real(fio,iparams->u_b.kUBB);
1410 iparams->u_b.thetaB=iparams->u_b.thetaA;
1411 iparams->u_b.kthetaB=iparams->u_b.kthetaA;
1412 iparams->u_b.r13B=iparams->u_b.r13A;
1413 iparams->u_b.kUBB=iparams->u_b.kUBA;
1416 case F_QUARTIC_ANGLES:
1417 gmx_fio_do_real(fio,iparams->qangle.theta);
1418 bDum=gmx_fio_ndo_real(fio,iparams->qangle.c,5);
1421 gmx_fio_do_real(fio,iparams->bham.a);
1422 gmx_fio_do_real(fio,iparams->bham.b);
1423 gmx_fio_do_real(fio,iparams->bham.c);
1426 gmx_fio_do_real(fio,iparams->morse.b0A);
1427 gmx_fio_do_real(fio,iparams->morse.cbA);
1428 gmx_fio_do_real(fio,iparams->morse.betaA);
1429 if (file_version >= 79) {
1430 gmx_fio_do_real(fio,iparams->morse.b0B);
1431 gmx_fio_do_real(fio,iparams->morse.cbB);
1432 gmx_fio_do_real(fio,iparams->morse.betaB);
1434 iparams->morse.b0B = iparams->morse.b0A;
1435 iparams->morse.cbB = iparams->morse.cbA;
1436 iparams->morse.betaB = iparams->morse.betaA;
1440 gmx_fio_do_real(fio,iparams->cubic.b0);
1441 gmx_fio_do_real(fio,iparams->cubic.kb);
1442 gmx_fio_do_real(fio,iparams->cubic.kcub);
1446 case F_POLARIZATION:
1447 gmx_fio_do_real(fio,iparams->polarize.alpha);
1450 gmx_fio_do_real(fio,iparams->anharm_polarize.alpha);
1451 gmx_fio_do_real(fio,iparams->anharm_polarize.drcut);
1452 gmx_fio_do_real(fio,iparams->anharm_polarize.khyp);
1455 if (file_version < 31)
1456 gmx_fatal(FARGS,"Old tpr files with water_polarization not supported. Make a new.");
1457 gmx_fio_do_real(fio,iparams->wpol.al_x);
1458 gmx_fio_do_real(fio,iparams->wpol.al_y);
1459 gmx_fio_do_real(fio,iparams->wpol.al_z);
1460 gmx_fio_do_real(fio,iparams->wpol.rOH);
1461 gmx_fio_do_real(fio,iparams->wpol.rHH);
1462 gmx_fio_do_real(fio,iparams->wpol.rOD);
1465 gmx_fio_do_real(fio,iparams->thole.a);
1466 gmx_fio_do_real(fio,iparams->thole.alpha1);
1467 gmx_fio_do_real(fio,iparams->thole.alpha2);
1468 gmx_fio_do_real(fio,iparams->thole.rfac);
1471 gmx_fio_do_real(fio,iparams->lj.c6);
1472 gmx_fio_do_real(fio,iparams->lj.c12);
1475 gmx_fio_do_real(fio,iparams->lj14.c6A);
1476 gmx_fio_do_real(fio,iparams->lj14.c12A);
1477 gmx_fio_do_real(fio,iparams->lj14.c6B);
1478 gmx_fio_do_real(fio,iparams->lj14.c12B);
1481 gmx_fio_do_real(fio,iparams->ljc14.fqq);
1482 gmx_fio_do_real(fio,iparams->ljc14.qi);
1483 gmx_fio_do_real(fio,iparams->ljc14.qj);
1484 gmx_fio_do_real(fio,iparams->ljc14.c6);
1485 gmx_fio_do_real(fio,iparams->ljc14.c12);
1487 case F_LJC_PAIRS_NB:
1488 gmx_fio_do_real(fio,iparams->ljcnb.qi);
1489 gmx_fio_do_real(fio,iparams->ljcnb.qj);
1490 gmx_fio_do_real(fio,iparams->ljcnb.c6);
1491 gmx_fio_do_real(fio,iparams->ljcnb.c12);
1497 gmx_fio_do_real(fio,iparams->pdihs.phiA);
1498 gmx_fio_do_real(fio,iparams->pdihs.cpA);
1499 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42) {
1500 /* Read the incorrectly stored multiplicity */
1501 gmx_fio_do_real(fio,iparams->harmonic.rB);
1502 gmx_fio_do_real(fio,iparams->harmonic.krB);
1503 iparams->pdihs.phiB = iparams->pdihs.phiA;
1504 iparams->pdihs.cpB = iparams->pdihs.cpA;
1506 gmx_fio_do_real(fio,iparams->pdihs.phiB);
1507 gmx_fio_do_real(fio,iparams->pdihs.cpB);
1508 gmx_fio_do_int(fio,iparams->pdihs.mult);
1512 gmx_fio_do_int(fio,iparams->disres.label);
1513 gmx_fio_do_int(fio,iparams->disres.type);
1514 gmx_fio_do_real(fio,iparams->disres.low);
1515 gmx_fio_do_real(fio,iparams->disres.up1);
1516 gmx_fio_do_real(fio,iparams->disres.up2);
1517 gmx_fio_do_real(fio,iparams->disres.kfac);
1520 gmx_fio_do_int(fio,iparams->orires.ex);
1521 gmx_fio_do_int(fio,iparams->orires.label);
1522 gmx_fio_do_int(fio,iparams->orires.power);
1523 gmx_fio_do_real(fio,iparams->orires.c);
1524 gmx_fio_do_real(fio,iparams->orires.obs);
1525 gmx_fio_do_real(fio,iparams->orires.kfac);
1528 if ( file_version < 72) {
1529 gmx_fio_do_int(fio,idum);
1530 gmx_fio_do_int(fio,idum);
1532 gmx_fio_do_real(fio,iparams->dihres.phiA);
1533 gmx_fio_do_real(fio,iparams->dihres.dphiA);
1534 gmx_fio_do_real(fio,iparams->dihres.kfacA);
1535 if (file_version >= 72) {
1536 gmx_fio_do_real(fio,iparams->dihres.phiB);
1537 gmx_fio_do_real(fio,iparams->dihres.dphiB);
1538 gmx_fio_do_real(fio,iparams->dihres.kfacB);
1540 iparams->dihres.phiB=iparams->dihres.phiA;
1541 iparams->dihres.dphiB=iparams->dihres.dphiA;
1542 iparams->dihres.kfacB=iparams->dihres.kfacA;
1546 gmx_fio_do_rvec(fio,iparams->posres.pos0A);
1547 gmx_fio_do_rvec(fio,iparams->posres.fcA);
1548 if (bRead && file_version < 27) {
1549 copy_rvec(iparams->posres.pos0A,iparams->posres.pos0B);
1550 copy_rvec(iparams->posres.fcA,iparams->posres.fcB);
1552 gmx_fio_do_rvec(fio,iparams->posres.pos0B);
1553 gmx_fio_do_rvec(fio,iparams->posres.fcB);
1557 gmx_fio_do_int(fio,iparams->fbposres.geom);
1558 gmx_fio_do_rvec(fio,iparams->fbposres.pos0);
1559 gmx_fio_do_real(fio,iparams->fbposres.r);
1560 gmx_fio_do_real(fio,iparams->fbposres.k);
1563 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA,NR_RBDIHS);
1564 if(file_version>=25)
1565 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB,NR_RBDIHS);
1568 /* Fourier dihedrals are internally represented
1569 * as Ryckaert-Bellemans since those are faster to compute.
1571 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA, NR_RBDIHS);
1572 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB, NR_RBDIHS);
1576 gmx_fio_do_real(fio,iparams->constr.dA);
1577 gmx_fio_do_real(fio,iparams->constr.dB);
1580 gmx_fio_do_real(fio,iparams->settle.doh);
1581 gmx_fio_do_real(fio,iparams->settle.dhh);
1584 gmx_fio_do_real(fio,iparams->vsite.a);
1589 gmx_fio_do_real(fio,iparams->vsite.a);
1590 gmx_fio_do_real(fio,iparams->vsite.b);
1595 gmx_fio_do_real(fio,iparams->vsite.a);
1596 gmx_fio_do_real(fio,iparams->vsite.b);
1597 gmx_fio_do_real(fio,iparams->vsite.c);
1600 gmx_fio_do_int(fio,iparams->vsiten.n);
1601 gmx_fio_do_real(fio,iparams->vsiten.a);
1606 /* We got rid of some parameters in version 68 */
1607 if(bRead && file_version<68)
1609 gmx_fio_do_real(fio,rdum);
1610 gmx_fio_do_real(fio,rdum);
1611 gmx_fio_do_real(fio,rdum);
1612 gmx_fio_do_real(fio,rdum);
1614 gmx_fio_do_real(fio,iparams->gb.sar);
1615 gmx_fio_do_real(fio,iparams->gb.st);
1616 gmx_fio_do_real(fio,iparams->gb.pi);
1617 gmx_fio_do_real(fio,iparams->gb.gbr);
1618 gmx_fio_do_real(fio,iparams->gb.bmlt);
1621 gmx_fio_do_int(fio,iparams->cmap.cmapA);
1622 gmx_fio_do_int(fio,iparams->cmap.cmapB);
1625 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
1626 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
1629 gmx_fio_unset_comment(fio);
1632 static void do_ilist(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,int file_version,
1639 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1641 if (file_version < 44) {
1642 for(i=0; i<MAXNODES; i++)
1643 gmx_fio_do_int(fio,idum);
1645 gmx_fio_do_int(fio,ilist->nr);
1647 snew(ilist->iatoms,ilist->nr);
1648 bDum=gmx_fio_ndo_int(fio,ilist->iatoms,ilist->nr);
1650 gmx_fio_unset_comment(fio);
1653 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1654 gmx_bool bRead, int file_version)
1660 gmx_fio_do_int(fio,ffparams->atnr);
1661 if (file_version < 57) {
1662 gmx_fio_do_int(fio,idum);
1664 gmx_fio_do_int(fio,ffparams->ntypes);
1666 fprintf(debug,"ffparams->atnr = %d, ntypes = %d\n",
1667 ffparams->atnr,ffparams->ntypes);
1669 snew(ffparams->functype,ffparams->ntypes);
1670 snew(ffparams->iparams,ffparams->ntypes);
1672 /* Read/write all the function types */
1673 bDum=gmx_fio_ndo_int(fio,ffparams->functype,ffparams->ntypes);
1675 pr_ivec(debug,0,"functype",ffparams->functype,ffparams->ntypes,TRUE);
1677 if (file_version >= 66) {
1678 gmx_fio_do_double(fio,ffparams->reppow);
1680 ffparams->reppow = 12.0;
1683 if (file_version >= 57) {
1684 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1687 /* Check whether all these function types are supported by the code.
1688 * In practice the code is backwards compatible, which means that the
1689 * numbering may have to be altered from old numbering to new numbering
1691 for (i=0; (i<ffparams->ntypes); i++) {
1693 /* Loop over file versions */
1694 for (k=0; (k<NFTUPD); k++)
1695 /* Compare the read file_version to the update table */
1696 if ((file_version < ftupd[k].fvnr) &&
1697 (ffparams->functype[i] >= ftupd[k].ftype)) {
1698 ffparams->functype[i] += 1;
1700 fprintf(debug,"Incrementing function type %d to %d (due to %s)\n",
1701 i,ffparams->functype[i],
1702 interaction_function[ftupd[k].ftype].longname);
1707 do_iparams(fio, ffparams->functype[i],&ffparams->iparams[i],bRead,
1710 pr_iparams(debug,ffparams->functype[i],&ffparams->iparams[i]);
1714 static void add_settle_atoms(t_ilist *ilist)
1718 /* Settle used to only store the first atom: add the other two */
1719 srenew(ilist->iatoms,2*ilist->nr);
1720 for(i=ilist->nr/2-1; i>=0; i--)
1722 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
1723 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
1724 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
1725 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
1727 ilist->nr = 2*ilist->nr;
1730 static void do_ilists(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,
1733 int i,j,renum[F_NRE];
1734 gmx_bool bDum=TRUE,bClear;
1737 for(j=0; (j<F_NRE); j++) {
1740 for (k=0; k<NFTUPD; k++)
1741 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
1745 ilist[j].iatoms = NULL;
1747 do_ilist(fio, &ilist[j],bRead,file_version,j);
1748 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
1750 add_settle_atoms(&ilist[j]);
1754 if (bRead && gmx_debug_at)
1755 pr_ilist(debug,0,interaction_function[j].longname,
1756 functype,&ilist[j],TRUE);
1761 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams,gmx_moltype_t *molt,
1762 gmx_bool bRead, int file_version)
1764 do_ffparams(fio, ffparams,bRead,file_version);
1766 if (file_version >= 54) {
1767 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1770 do_ilists(fio, molt->ilist,bRead,file_version);
1773 static void do_block(t_fileio *fio, t_block *block,gmx_bool bRead,int file_version)
1775 int i,idum,dum_nra,*dum_a;
1778 if (file_version < 44)
1779 for(i=0; i<MAXNODES; i++)
1780 gmx_fio_do_int(fio,idum);
1781 gmx_fio_do_int(fio,block->nr);
1782 if (file_version < 51)
1783 gmx_fio_do_int(fio,dum_nra);
1785 block->nalloc_index = block->nr+1;
1786 snew(block->index,block->nalloc_index);
1788 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1790 if (file_version < 51 && dum_nra > 0) {
1791 snew(dum_a,dum_nra);
1792 bDum=gmx_fio_ndo_int(fio,dum_a,dum_nra);
1797 static void do_blocka(t_fileio *fio, t_blocka *block,gmx_bool bRead,
1803 if (file_version < 44)
1804 for(i=0; i<MAXNODES; i++)
1805 gmx_fio_do_int(fio,idum);
1806 gmx_fio_do_int(fio,block->nr);
1807 gmx_fio_do_int(fio,block->nra);
1809 block->nalloc_index = block->nr+1;
1810 snew(block->index,block->nalloc_index);
1811 block->nalloc_a = block->nra;
1812 snew(block->a,block->nalloc_a);
1814 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1815 bDum=gmx_fio_ndo_int(fio,block->a,block->nra);
1818 static void do_atom(t_fileio *fio, t_atom *atom,int ngrp,gmx_bool bRead,
1819 int file_version, gmx_groups_t *groups,int atnr)
1823 gmx_fio_do_real(fio,atom->m);
1824 gmx_fio_do_real(fio,atom->q);
1825 gmx_fio_do_real(fio,atom->mB);
1826 gmx_fio_do_real(fio,atom->qB);
1827 gmx_fio_do_ushort(fio, atom->type);
1828 gmx_fio_do_ushort(fio, atom->typeB);
1829 gmx_fio_do_int(fio,atom->ptype);
1830 gmx_fio_do_int(fio,atom->resind);
1831 if (file_version >= 52)
1832 gmx_fio_do_int(fio,atom->atomnumber);
1834 atom->atomnumber = NOTSET;
1835 if (file_version < 23)
1837 else if (file_version < 39)
1842 if (file_version < 57) {
1843 unsigned char uchar[egcNR];
1844 gmx_fio_ndo_uchar(fio,uchar,myngrp);
1845 for(i=myngrp; (i<ngrp); i++) {
1848 /* Copy the old data format to the groups struct */
1849 for(i=0; i<ngrp; i++) {
1850 groups->grpnr[i][atnr] = uchar[i];
1855 static void do_grps(t_fileio *fio, int ngrp,t_grps grps[],gmx_bool bRead,
1861 if (file_version < 23)
1863 else if (file_version < 39)
1868 for(j=0; (j<ngrp); j++) {
1870 gmx_fio_do_int(fio,grps[j].nr);
1872 snew(grps[j].nm_ind,grps[j].nr);
1873 bDum=gmx_fio_ndo_int(fio,grps[j].nm_ind,grps[j].nr);
1877 snew(grps[j].nm_ind,grps[j].nr);
1882 static void do_symstr(t_fileio *fio, char ***nm,gmx_bool bRead,t_symtab *symtab)
1887 gmx_fio_do_int(fio,ls);
1888 *nm = get_symtab_handle(symtab,ls);
1891 ls = lookup_symtab(symtab,*nm);
1892 gmx_fio_do_int(fio,ls);
1896 static void do_strstr(t_fileio *fio, int nstr,char ***nm,gmx_bool bRead,
1901 for (j=0; (j<nstr); j++)
1902 do_symstr(fio, &(nm[j]),bRead,symtab);
1905 static void do_resinfo(t_fileio *fio, int n,t_resinfo *ri,gmx_bool bRead,
1906 t_symtab *symtab, int file_version)
1910 for (j=0; (j<n); j++) {
1911 do_symstr(fio, &(ri[j].name),bRead,symtab);
1912 if (file_version >= 63) {
1913 gmx_fio_do_int(fio,ri[j].nr);
1914 gmx_fio_do_uchar(fio, ri[j].ic);
1922 static void do_atoms(t_fileio *fio, t_atoms *atoms,gmx_bool bRead,t_symtab *symtab,
1924 gmx_groups_t *groups)
1928 gmx_fio_do_int(fio,atoms->nr);
1929 gmx_fio_do_int(fio,atoms->nres);
1930 if (file_version < 57) {
1931 gmx_fio_do_int(fio,groups->ngrpname);
1932 for(i=0; i<egcNR; i++) {
1933 groups->ngrpnr[i] = atoms->nr;
1934 snew(groups->grpnr[i],groups->ngrpnr[i]);
1938 snew(atoms->atom,atoms->nr);
1939 snew(atoms->atomname,atoms->nr);
1940 snew(atoms->atomtype,atoms->nr);
1941 snew(atoms->atomtypeB,atoms->nr);
1942 snew(atoms->resinfo,atoms->nres);
1943 if (file_version < 57) {
1944 snew(groups->grpname,groups->ngrpname);
1946 atoms->pdbinfo = NULL;
1948 for(i=0; (i<atoms->nr); i++) {
1949 do_atom(fio, &atoms->atom[i],egcNR,bRead, file_version,groups,i);
1951 do_strstr(fio, atoms->nr,atoms->atomname,bRead,symtab);
1952 if (bRead && (file_version <= 20)) {
1953 for(i=0; i<atoms->nr; i++) {
1954 atoms->atomtype[i] = put_symtab(symtab,"?");
1955 atoms->atomtypeB[i] = put_symtab(symtab,"?");
1958 do_strstr(fio, atoms->nr,atoms->atomtype,bRead,symtab);
1959 do_strstr(fio, atoms->nr,atoms->atomtypeB,bRead,symtab);
1961 do_resinfo(fio, atoms->nres,atoms->resinfo,bRead,symtab,file_version);
1963 if (file_version < 57) {
1964 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1966 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1970 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
1971 gmx_bool bRead,t_symtab *symtab,
1977 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1978 gmx_fio_do_int(fio,groups->ngrpname);
1980 snew(groups->grpname,groups->ngrpname);
1982 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1983 for(g=0; g<egcNR; g++) {
1984 gmx_fio_do_int(fio,groups->ngrpnr[g]);
1985 if (groups->ngrpnr[g] == 0) {
1987 groups->grpnr[g] = NULL;
1991 snew(groups->grpnr[g],groups->ngrpnr[g]);
1993 bDum=gmx_fio_ndo_uchar(fio, groups->grpnr[g],groups->ngrpnr[g]);
1998 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes,gmx_bool bRead,
1999 t_symtab *symtab,int file_version)
2002 gmx_bool bDum = TRUE;
2004 if (file_version > 25) {
2005 gmx_fio_do_int(fio,atomtypes->nr);
2008 snew(atomtypes->radius,j);
2009 snew(atomtypes->vol,j);
2010 snew(atomtypes->surftens,j);
2011 snew(atomtypes->atomnumber,j);
2012 snew(atomtypes->gb_radius,j);
2013 snew(atomtypes->S_hct,j);
2015 bDum=gmx_fio_ndo_real(fio,atomtypes->radius,j);
2016 bDum=gmx_fio_ndo_real(fio,atomtypes->vol,j);
2017 bDum=gmx_fio_ndo_real(fio,atomtypes->surftens,j);
2018 if(file_version >= 40)
2020 bDum=gmx_fio_ndo_int(fio,atomtypes->atomnumber,j);
2022 if(file_version >= 60)
2024 bDum=gmx_fio_ndo_real(fio,atomtypes->gb_radius,j);
2025 bDum=gmx_fio_ndo_real(fio,atomtypes->S_hct,j);
2028 /* File versions prior to 26 cannot do GBSA,
2029 * so they dont use this structure
2032 atomtypes->radius = NULL;
2033 atomtypes->vol = NULL;
2034 atomtypes->surftens = NULL;
2035 atomtypes->atomnumber = NULL;
2036 atomtypes->gb_radius = NULL;
2037 atomtypes->S_hct = NULL;
2041 static void do_symtab(t_fileio *fio, t_symtab *symtab,gmx_bool bRead)
2047 gmx_fio_do_int(fio,symtab->nr);
2050 snew(symtab->symbuf,1);
2051 symbuf = symtab->symbuf;
2052 symbuf->bufsize = nr;
2053 snew(symbuf->buf,nr);
2054 for (i=0; (i<nr); i++) {
2055 gmx_fio_do_string(fio,buf);
2056 symbuf->buf[i]=strdup(buf);
2060 symbuf = symtab->symbuf;
2061 while (symbuf!=NULL) {
2062 for (i=0; (i<symbuf->bufsize) && (i<nr); i++)
2063 gmx_fio_do_string(fio,symbuf->buf[i]);
2065 symbuf=symbuf->next;
2068 gmx_fatal(FARGS,"nr of symtab strings left: %d",nr);
2072 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2074 int i,j,ngrid,gs,nelem;
2076 gmx_fio_do_int(fio,cmap_grid->ngrid);
2077 gmx_fio_do_int(fio,cmap_grid->grid_spacing);
2079 ngrid = cmap_grid->ngrid;
2080 gs = cmap_grid->grid_spacing;
2085 snew(cmap_grid->cmapdata,ngrid);
2087 for(i=0;i<cmap_grid->ngrid;i++)
2089 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
2093 for(i=0;i<cmap_grid->ngrid;i++)
2095 for(j=0;j<nelem;j++)
2097 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4]);
2098 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+1]);
2099 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+2]);
2100 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+3]);
2106 void tpx_make_chain_identifiers(t_atoms *atoms,t_block *mols)
2112 /* We always assign a new chain number, but save the chain id characters
2113 * for larger molecules.
2115 #define CHAIN_MIN_ATOMS 15
2119 for(m=0; m<mols->nr; m++)
2122 a1=mols->index[m+1];
2123 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2132 for(a=a0; a<a1; a++)
2134 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2135 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2140 /* Blank out the chain id if there was only one chain */
2143 for(r=0; r<atoms->nres; r++)
2145 atoms->resinfo[r].chainid = ' ';
2150 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt,gmx_bool bRead,
2151 t_symtab *symtab, int file_version,
2152 gmx_groups_t *groups)
2156 if (file_version >= 57) {
2157 do_symstr(fio, &(molt->name),bRead,symtab);
2160 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2162 if (bRead && gmx_debug_at) {
2163 pr_atoms(debug,0,"atoms",&molt->atoms,TRUE);
2166 if (file_version >= 57) {
2167 do_ilists(fio, molt->ilist,bRead,file_version);
2169 do_block(fio, &molt->cgs,bRead,file_version);
2170 if (bRead && gmx_debug_at) {
2171 pr_block(debug,0,"cgs",&molt->cgs,TRUE);
2175 /* This used to be in the atoms struct */
2176 do_blocka(fio, &molt->excls, bRead, file_version);
2179 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb,gmx_bool bRead,
2184 gmx_fio_do_int(fio,molb->type);
2185 gmx_fio_do_int(fio,molb->nmol);
2186 gmx_fio_do_int(fio,molb->natoms_mol);
2187 /* Position restraint coordinates */
2188 gmx_fio_do_int(fio,molb->nposres_xA);
2189 if (molb->nposres_xA > 0) {
2191 snew(molb->posres_xA,molb->nposres_xA);
2193 gmx_fio_ndo_rvec(fio,molb->posres_xA,molb->nposres_xA);
2195 gmx_fio_do_int(fio,molb->nposres_xB);
2196 if (molb->nposres_xB > 0) {
2198 snew(molb->posres_xB,molb->nposres_xB);
2200 gmx_fio_ndo_rvec(fio,molb->posres_xB,molb->nposres_xB);
2205 static t_block mtop_mols(gmx_mtop_t *mtop)
2211 for(mb=0; mb<mtop->nmolblock; mb++) {
2212 mols.nr += mtop->molblock[mb].nmol;
2214 mols.nalloc_index = mols.nr + 1;
2215 snew(mols.index,mols.nalloc_index);
2220 for(mb=0; mb<mtop->nmolblock; mb++) {
2221 for(mol=0; mol<mtop->molblock[mb].nmol; mol++) {
2222 a += mtop->molblock[mb].natoms_mol;
2231 static void add_posres_molblock(gmx_mtop_t *mtop)
2236 gmx_molblock_t *molb;
2239 /* posres reference positions are stored in ip->posres (if present) and
2240 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2241 posres.pos0A are identical to fbposres.pos0. */
2242 il = &mtop->moltype[0].ilist[F_POSRES];
2243 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2244 if (il->nr == 0 && ilfb->nr == 0) {
2249 for(i=0; i<il->nr; i+=2) {
2250 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2251 am = max(am,il->iatoms[i+1]);
2252 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2253 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2254 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ]) {
2258 /* This loop is required if we have only flat-bottomed posres:
2260 - bFE == FALSE (no B-state for flat-bottomed posres) */
2263 for(i=0; i<ilfb->nr; i+=2) {
2264 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2265 am = max(am,ilfb->iatoms[i+1]);
2268 /* Make the posres coordinate block end at a molecule end */
2270 while(am >= mtop->mols.index[mol+1]) {
2273 molb = &mtop->molblock[0];
2274 molb->nposres_xA = mtop->mols.index[mol+1];
2275 snew(molb->posres_xA,molb->nposres_xA);
2277 molb->nposres_xB = molb->nposres_xA;
2278 snew(molb->posres_xB,molb->nposres_xB);
2280 molb->nposres_xB = 0;
2282 for(i=0; i<il->nr; i+=2) {
2283 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2284 a = il->iatoms[i+1];
2285 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2286 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2287 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2289 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2290 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2291 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2296 /* If only flat-bottomed posres are present, take reference pos from them.
2297 Here: bFE == FALSE */
2298 for(i=0; i<ilfb->nr; i+=2)
2300 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2301 a = ilfb->iatoms[i+1];
2302 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2303 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2304 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2309 static void set_disres_npair(gmx_mtop_t *mtop)
2316 ip = mtop->ffparams.iparams;
2318 for(mt=0; mt<mtop->nmoltype; mt++) {
2319 il = &mtop->moltype[mt].ilist[F_DISRES];
2323 for(i=0; i<il->nr; i+=3) {
2325 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label) {
2326 ip[a[i]].disres.npair = npair;
2334 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop,gmx_bool bRead,
2342 do_symtab(fio, &(mtop->symtab),bRead);
2344 pr_symtab(debug,0,"symtab",&mtop->symtab);
2346 do_symstr(fio, &(mtop->name),bRead,&(mtop->symtab));
2348 if (file_version >= 57) {
2349 do_ffparams(fio, &mtop->ffparams,bRead,file_version);
2351 gmx_fio_do_int(fio,mtop->nmoltype);
2356 snew(mtop->moltype,mtop->nmoltype);
2357 if (file_version < 57) {
2358 mtop->moltype[0].name = mtop->name;
2361 for(mt=0; mt<mtop->nmoltype; mt++) {
2362 do_moltype(fio, &mtop->moltype[mt],bRead,&mtop->symtab,file_version,
2366 if (file_version >= 57) {
2367 gmx_fio_do_int(fio,mtop->nmolblock);
2369 mtop->nmolblock = 1;
2372 snew(mtop->molblock,mtop->nmolblock);
2374 if (file_version >= 57) {
2375 for(mb=0; mb<mtop->nmolblock; mb++) {
2376 do_molblock(fio, &mtop->molblock[mb],bRead,file_version);
2378 gmx_fio_do_int(fio,mtop->natoms);
2380 mtop->molblock[0].type = 0;
2381 mtop->molblock[0].nmol = 1;
2382 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
2383 mtop->molblock[0].nposres_xA = 0;
2384 mtop->molblock[0].nposres_xB = 0;
2387 do_atomtypes (fio, &(mtop->atomtypes),bRead,&(mtop->symtab), file_version);
2389 pr_atomtypes(debug,0,"atomtypes",&mtop->atomtypes,TRUE);
2391 if (file_version < 57) {
2392 /* Debug statements are inside do_idef */
2393 do_idef (fio, &mtop->ffparams,&mtop->moltype[0],bRead,file_version);
2394 mtop->natoms = mtop->moltype[0].atoms.nr;
2397 if(file_version >= 65)
2399 do_cmap(fio, &mtop->ffparams.cmap_grid,bRead);
2403 mtop->ffparams.cmap_grid.ngrid = 0;
2404 mtop->ffparams.cmap_grid.grid_spacing = 0;
2405 mtop->ffparams.cmap_grid.cmapdata = NULL;
2408 if (file_version >= 57) {
2409 do_groups(fio, &mtop->groups,bRead,&(mtop->symtab),file_version);
2412 if (file_version < 57) {
2413 do_block(fio, &mtop->moltype[0].cgs,bRead,file_version);
2414 if (bRead && gmx_debug_at) {
2415 pr_block(debug,0,"cgs",&mtop->moltype[0].cgs,TRUE);
2417 do_block(fio, &mtop->mols,bRead,file_version);
2418 /* Add the posres coordinates to the molblock */
2419 add_posres_molblock(mtop);
2422 if (file_version >= 57) {
2423 mtop->mols = mtop_mols(mtop);
2426 pr_block(debug,0,"mols",&mtop->mols,TRUE);
2430 if (file_version < 51) {
2431 /* Here used to be the shake blocks */
2432 do_blocka(fio, &dumb,bRead,file_version);
2440 close_symtab(&(mtop->symtab));
2444 /* If TopOnlyOK is TRUE then we can read even future versions
2445 * of tpx files, provided the file_generation hasn't changed.
2446 * If it is FALSE, we need the inputrecord too, and bail out
2447 * if the file is newer than the program.
2449 * The version and generation if the topology (see top of this file)
2450 * are returned in the two last arguments.
2452 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2454 static void do_tpxheader(t_fileio *fio,gmx_bool bRead,t_tpxheader *tpx,
2455 gmx_bool TopOnlyOK, int *file_version,
2456 int *file_generation)
2459 char file_tag[STRLEN];
2466 gmx_fio_checktype(fio);
2467 gmx_fio_setdebug(fio,bDebugMode());
2469 /* NEW! XDR tpb file */
2470 precision = sizeof(real);
2472 gmx_fio_do_string(fio,buf);
2473 if (strncmp(buf,"VERSION",7))
2474 gmx_fatal(FARGS,"Can not read file %s,\n"
2475 " this file is from a Gromacs version which is older than 2.0\n"
2476 " Make a new one with grompp or use a gro or pdb file, if possible",
2477 gmx_fio_getname(fio));
2478 gmx_fio_do_int(fio,precision);
2479 bDouble = (precision == sizeof(double));
2480 if ((precision != sizeof(float)) && !bDouble)
2481 gmx_fatal(FARGS,"Unknown precision in file %s: real is %d bytes "
2482 "instead of %d or %d",
2483 gmx_fio_getname(fio),precision,sizeof(float),sizeof(double));
2484 gmx_fio_setprecision(fio,bDouble);
2485 fprintf(stderr,"Reading file %s, %s (%s precision)\n",
2486 gmx_fio_getname(fio),buf,bDouble ? "double" : "single");
2489 gmx_fio_write_string(fio,GromacsVersion());
2490 bDouble = (precision == sizeof(double));
2491 gmx_fio_setprecision(fio,bDouble);
2492 gmx_fio_do_int(fio,precision);
2494 sprintf(file_tag,"%s",tpx_tag);
2495 fgen = tpx_generation;
2498 /* Check versions! */
2499 gmx_fio_do_int(fio,fver);
2501 /* This is for backward compatibility with development versions 77-79
2502 * where the tag was, mistakenly, placed before the generation,
2503 * which would cause a segv instead of a proper error message
2504 * when reading the topology only from tpx with <77 code.
2506 if (fver >= 77 && fver <= 79)
2508 gmx_fio_do_string(fio,file_tag);
2513 gmx_fio_do_int(fio,fgen);
2522 gmx_fio_do_string(fio,file_tag);
2528 /* Versions before 77 don't have the tag, set it to release */
2529 sprintf(file_tag,"%s",TPX_TAG_RELEASE);
2532 if (strcmp(file_tag,tpx_tag) != 0)
2534 fprintf(stderr,"Note: file tpx tag '%s', software tpx tag '%s'\n",
2537 /* We only support reading tpx files with the same tag as the code
2538 * or tpx files with the release tag and with lower version number.
2540 if (!strcmp(file_tag,TPX_TAG_RELEASE) == 0 && fver < tpx_version)
2542 gmx_fatal(FARGS,"tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2543 gmx_fio_getname(fio),fver,file_tag,
2544 tpx_version,tpx_tag);
2549 if (file_version != NULL)
2551 *file_version = fver;
2553 if (file_generation != NULL)
2555 *file_generation = fgen;
2559 if ((fver <= tpx_incompatible_version) ||
2560 ((fver > tpx_version) && !TopOnlyOK) ||
2561 (fgen > tpx_generation))
2562 gmx_fatal(FARGS,"reading tpx file (%s) version %d with version %d program",
2563 gmx_fio_getname(fio),fver,tpx_version);
2565 do_section(fio,eitemHEADER,bRead);
2566 gmx_fio_do_int(fio,tpx->natoms);
2568 gmx_fio_do_int(fio,tpx->ngtc);
2572 gmx_fio_do_int(fio,idum);
2573 gmx_fio_do_real(fio,rdum);
2575 /*a better decision will eventually (5.0 or later) need to be made
2576 on how to treat the alchemical state of the system, which can now
2577 vary through a simulation, and cannot be completely described
2578 though a single lambda variable, or even a single state
2579 index. Eventually, should probably be a vector. MRS*/
2582 gmx_fio_do_int(fio,tpx->fep_state);
2584 gmx_fio_do_real(fio,tpx->lambda);
2585 gmx_fio_do_int(fio,tpx->bIr);
2586 gmx_fio_do_int(fio,tpx->bTop);
2587 gmx_fio_do_int(fio,tpx->bX);
2588 gmx_fio_do_int(fio,tpx->bV);
2589 gmx_fio_do_int(fio,tpx->bF);
2590 gmx_fio_do_int(fio,tpx->bBox);
2592 if((fgen > tpx_generation)) {
2593 /* This can only happen if TopOnlyOK=TRUE */
2598 static int do_tpx(t_fileio *fio, gmx_bool bRead,
2599 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop,
2600 gmx_bool bXVallocated)
2605 gmx_bool TopOnlyOK,bDum=TRUE;
2606 int file_version,file_generation;
2610 gmx_bool bPeriodicMols;
2613 tpx.natoms = state->natoms;
2614 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
2615 tpx.fep_state = state->fep_state;
2616 tpx.lambda = state->lambda[efptFEP];
2617 tpx.bIr = (ir != NULL);
2618 tpx.bTop = (mtop != NULL);
2619 tpx.bX = (state->x != NULL);
2620 tpx.bV = (state->v != NULL);
2621 tpx.bF = (f != NULL);
2625 TopOnlyOK = (ir==NULL);
2627 do_tpxheader(fio,bRead,&tpx,TopOnlyOK,&file_version,&file_generation);
2631 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
2632 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
2636 init_state(state,0,tpx.ngtc,0,0,0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
2637 state->natoms = tpx.natoms;
2638 state->nalloc = tpx.natoms;
2642 init_state(state,tpx.natoms,tpx.ngtc,0,0,0); /* nose-hoover chains */
2646 #define do_test(fio,b,p) if (bRead && (p!=NULL) && !b) gmx_fatal(FARGS,"No %s in %s",#p,gmx_fio_getname(fio))
2648 do_test(fio,tpx.bBox,state->box);
2649 do_section(fio,eitemBOX,bRead);
2651 gmx_fio_ndo_rvec(fio,state->box,DIM);
2652 if (file_version >= 51) {
2653 gmx_fio_ndo_rvec(fio,state->box_rel,DIM);
2655 /* We initialize box_rel after reading the inputrec */
2656 clear_mat(state->box_rel);
2658 if (file_version >= 28) {
2659 gmx_fio_ndo_rvec(fio,state->boxv,DIM);
2660 if (file_version < 56) {
2662 gmx_fio_ndo_rvec(fio,mdum,DIM);
2667 if (state->ngtc > 0 && file_version >= 28) {
2669 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
2670 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
2671 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
2672 snew(dumv,state->ngtc);
2673 if (file_version < 69) {
2674 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2676 /* These used to be the Berendsen tcoupl_lambda's */
2677 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2681 /* Prior to tpx version 26, the inputrec was here.
2682 * I moved it to enable partial forward-compatibility
2683 * for analysis/viewer programs.
2685 if(file_version<26) {
2686 do_test(fio,tpx.bIr,ir);
2687 do_section(fio,eitemIR,bRead);
2690 do_inputrec(fio, ir,bRead,file_version,
2691 mtop ? &mtop->ffparams.fudgeQQ : NULL);
2693 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2696 do_inputrec(fio, &dum_ir,bRead,file_version,
2697 mtop ? &mtop->ffparams.fudgeQQ :NULL);
2699 pr_inputrec(debug,0,"inputrec",&dum_ir,FALSE);
2700 done_inputrec(&dum_ir);
2706 do_test(fio,tpx.bTop,mtop);
2707 do_section(fio,eitemTOP,bRead);
2709 int mtop_file_version = file_version;
2710 /*allow reading of Gromacs 4.6 files*/
2711 if (mtop_file_version>80 && mtop_file_version<90)
2713 mtop_file_version = 79;
2716 do_mtop(fio,mtop,bRead, mtop_file_version);
2718 do_mtop(fio,&dum_top,bRead,mtop_file_version);
2719 done_mtop(&dum_top,TRUE);
2722 do_test(fio,tpx.bX,state->x);
2723 do_section(fio,eitemX,bRead);
2726 state->flags |= (1<<estX);
2728 gmx_fio_ndo_rvec(fio,state->x,state->natoms);
2731 do_test(fio,tpx.bV,state->v);
2732 do_section(fio,eitemV,bRead);
2735 state->flags |= (1<<estV);
2737 gmx_fio_ndo_rvec(fio,state->v,state->natoms);
2740 do_test(fio,tpx.bF,f);
2741 do_section(fio,eitemF,bRead);
2742 if (tpx.bF) gmx_fio_ndo_rvec(fio,f,state->natoms);
2744 /* Starting with tpx version 26, we have the inputrec
2745 * at the end of the file, so we can ignore it
2746 * if the file is never than the software (but still the
2747 * same generation - see comments at the top of this file.
2752 bPeriodicMols = FALSE;
2753 if (file_version >= 26) {
2754 do_test(fio,tpx.bIr,ir);
2755 do_section(fio,eitemIR,bRead);
2757 if (file_version >= 53) {
2758 /* Removed the pbc info from do_inputrec, since we always want it */
2761 bPeriodicMols = ir->bPeriodicMols;
2763 gmx_fio_do_int(fio,ePBC);
2764 gmx_fio_do_gmx_bool(fio,bPeriodicMols);
2766 if (file_generation <= tpx_generation && ir) {
2767 do_inputrec(fio, ir,bRead,file_version,mtop ? &mtop->ffparams.fudgeQQ : NULL);
2769 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2770 if (file_version < 51)
2771 set_box_rel(ir,state);
2772 if (file_version < 53) {
2774 bPeriodicMols = ir->bPeriodicMols;
2777 if (bRead && ir && file_version >= 53) {
2778 /* We need to do this after do_inputrec, since that initializes ir */
2780 ir->bPeriodicMols = bPeriodicMols;
2789 if (state->ngtc == 0)
2791 /* Reading old version without tcoupl state data: set it */
2792 init_gtc_state(state,ir->opts.ngtc,0,ir->opts.nhchainlength);
2794 if (tpx.bTop && mtop)
2796 if (file_version < 57)
2798 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
2800 ir->eDisre = edrSimple;
2804 ir->eDisre = edrNone;
2807 set_disres_npair(mtop);
2811 if (tpx.bTop && mtop)
2813 gmx_mtop_finalize(mtop);
2816 if (file_version >= 57)
2820 env = getenv("GMX_NOCHARGEGROUPS");
2823 sscanf(env,"%d",&ienv);
2824 fprintf(stderr,"\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
2829 "Will make single atomic charge groups in non-solvent%s\n",
2830 ienv > 1 ? " and solvent" : "");
2831 gmx_mtop_make_atomic_charge_groups(mtop,ienv==1);
2833 fprintf(stderr,"\n");
2841 /************************************************************
2843 * The following routines are the exported ones
2845 ************************************************************/
2847 t_fileio *open_tpx(const char *fn,const char *mode)
2849 return gmx_fio_open(fn,mode);
2852 void close_tpx(t_fileio *fio)
2857 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
2858 int *file_version, int *file_generation)
2862 fio = open_tpx(fn,"r");
2863 do_tpxheader(fio,TRUE,tpx,TopOnlyOK,file_version,file_generation);
2867 void write_tpx_state(const char *fn,
2868 t_inputrec *ir,t_state *state,gmx_mtop_t *mtop)
2872 fio = open_tpx(fn,"w");
2873 do_tpx(fio,FALSE,ir,state,NULL,mtop,FALSE);
2877 void read_tpx_state(const char *fn,
2878 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop)
2882 fio = open_tpx(fn,"r");
2883 do_tpx(fio,TRUE,ir,state,f,mtop,FALSE);
2887 int read_tpx(const char *fn,
2888 t_inputrec *ir, matrix box,int *natoms,
2889 rvec *x,rvec *v,rvec *f,gmx_mtop_t *mtop)
2897 fio = open_tpx(fn,"r");
2898 ePBC = do_tpx(fio,TRUE,ir,&state,f,mtop,TRUE);
2900 *natoms = state.natoms;
2902 copy_mat(state.box,box);
2910 int read_tpx_top(const char *fn,
2911 t_inputrec *ir, matrix box,int *natoms,
2912 rvec *x,rvec *v,rvec *f,t_topology *top)
2918 ePBC = read_tpx(fn,ir,box,natoms,x,v,f,&mtop);
2920 *top = gmx_mtop_t_to_t_topology(&mtop);
2925 gmx_bool fn2bTPX(const char *file)
2927 switch (fn2ftp(file)) {
2937 gmx_bool read_tps_conf(const char *infile,char *title,t_topology *top,int *ePBC,
2938 rvec **x,rvec **v,matrix box,gmx_bool bMass)
2941 int natoms,i,version,generation;
2942 gmx_bool bTop,bXNULL=FALSE;
2944 t_topology *topconv;
2947 bTop = fn2bTPX(infile);
2950 read_tpxheader(infile,&header,TRUE,&version,&generation);
2952 snew(*x,header.natoms);
2954 snew(*v,header.natoms);
2956 *ePBC = read_tpx(infile,NULL,box,&natoms,
2957 (x==NULL) ? NULL : *x,(v==NULL) ? NULL : *v,NULL,mtop);
2958 *top = gmx_mtop_t_to_t_topology(mtop);
2960 strcpy(title,*top->name);
2961 tpx_make_chain_identifiers(&top->atoms,&top->mols);
2964 get_stx_coordnum(infile,&natoms);
2965 init_t_atoms(&top->atoms,natoms,(fn2ftp(infile) == efPDB));
2974 read_stx_conf(infile,title,&top->atoms,*x,(v==NULL) ? NULL : *v,ePBC,box);
2981 aps = gmx_atomprop_init();
2982 for(i=0; (i<natoms); i++)
2983 if (!gmx_atomprop_query(aps,epropMass,
2984 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2985 *top->atoms.atomname[i],
2986 &(top->atoms.atom[i].m))) {
2988 fprintf(debug,"Can not find mass for atom %s %d %s, setting to 1\n",
2989 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2990 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
2991 *top->atoms.atomname[i]);
2993 gmx_atomprop_destroy(aps);
2995 top->idef.ntypes=-1;