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 = 80;
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 },
163 { 22, F_DISRESVIOL },
167 { 26, F_DIHRESVIOL },
172 { 46, F_ECONSERVED },
176 { 76, F_ANHARM_POL },
179 { 79, F_DVDL_BONDED, },
180 { 79, F_DVDL_RESTRAINT },
181 { 79, F_DVDL_TEMPERATURE },
184 #define NFTUPD asize(ftupd)
186 /* Needed for backward compatibility */
189 static void _do_section(t_fileio *fio,int key,gmx_bool bRead,const char *src,
195 if (gmx_fio_getftp(fio) == efTPA) {
197 gmx_fio_write_string(fio,itemstr[key]);
198 bDbg = gmx_fio_getdebug(fio);
199 gmx_fio_setdebug(fio,FALSE);
200 gmx_fio_write_string(fio,comment_str[key]);
201 gmx_fio_setdebug(fio,bDbg);
204 if (gmx_fio_getdebug(fio))
205 fprintf(stderr,"Looking for section %s (%s, %d)",
206 itemstr[key],src,line);
209 gmx_fio_do_string(fio,buf);
210 } while ((gmx_strcasecmp(buf,itemstr[key]) != 0));
212 if (gmx_strcasecmp(buf,itemstr[key]) != 0)
213 gmx_fatal(FARGS,"\nCould not find section heading %s",itemstr[key]);
214 else if (gmx_fio_getdebug(fio))
215 fprintf(stderr," and found it\n");
220 #define do_section(fio,key,bRead) _do_section(fio,key,bRead,__FILE__,__LINE__)
222 /**************************************************************
224 * Now the higer level routines that do io of the structures and arrays
226 **************************************************************/
227 static void do_pullgrp(t_fileio *fio, t_pullgrp *pgrp, gmx_bool bRead,
233 gmx_fio_do_int(fio,pgrp->nat);
235 snew(pgrp->ind,pgrp->nat);
236 bDum=gmx_fio_ndo_int(fio,pgrp->ind,pgrp->nat);
237 gmx_fio_do_int(fio,pgrp->nweight);
239 snew(pgrp->weight,pgrp->nweight);
240 bDum=gmx_fio_ndo_real(fio,pgrp->weight,pgrp->nweight);
241 gmx_fio_do_int(fio,pgrp->pbcatom);
242 gmx_fio_do_rvec(fio,pgrp->vec);
243 gmx_fio_do_rvec(fio,pgrp->init);
244 gmx_fio_do_real(fio,pgrp->rate);
245 gmx_fio_do_real(fio,pgrp->k);
246 if (file_version >= 56) {
247 gmx_fio_do_real(fio,pgrp->kB);
253 static void do_expandedvals(t_fileio *fio,t_expanded *expand,int n_lambda, gmx_bool bRead, int file_version)
255 /* i is used in the ndo_double macro*/
261 if (file_version >= 79)
267 snew(expand->init_lambda_weights,n_lambda);
269 bDum=gmx_fio_ndo_real(fio,expand->init_lambda_weights,n_lambda);
270 gmx_fio_do_gmx_bool(fio,expand->bInit_weights);
273 gmx_fio_do_int(fio,expand->nstexpanded);
274 gmx_fio_do_int(fio,expand->elmcmove);
275 gmx_fio_do_int(fio,expand->elamstats);
276 gmx_fio_do_int(fio,expand->lmc_repeats);
277 gmx_fio_do_int(fio,expand->gibbsdeltalam);
278 gmx_fio_do_int(fio,expand->lmc_forced_nstart);
279 gmx_fio_do_int(fio,expand->lmc_seed);
280 gmx_fio_do_real(fio,expand->mc_temp);
281 gmx_fio_do_int(fio,expand->bSymmetrizedTMatrix);
282 gmx_fio_do_int(fio,expand->nstTij);
283 gmx_fio_do_int(fio,expand->minvarmin);
284 gmx_fio_do_int(fio,expand->c_range);
285 gmx_fio_do_real(fio,expand->wl_scale);
286 gmx_fio_do_real(fio,expand->wl_ratio);
287 gmx_fio_do_real(fio,expand->init_wl_delta);
288 gmx_fio_do_gmx_bool(fio,expand->bWLoneovert);
289 gmx_fio_do_int(fio,expand->elmceq);
290 gmx_fio_do_int(fio,expand->equil_steps);
291 gmx_fio_do_int(fio,expand->equil_samples);
292 gmx_fio_do_int(fio,expand->equil_n_at_lam);
293 gmx_fio_do_real(fio,expand->equil_wl_delta);
294 gmx_fio_do_real(fio,expand->equil_ratio);
298 static void do_simtempvals(t_fileio *fio,t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
303 if (file_version >= 79)
305 gmx_fio_do_int(fio,simtemp->eSimTempScale);
306 gmx_fio_do_real(fio,simtemp->simtemp_high);
307 gmx_fio_do_real(fio,simtemp->simtemp_low);
312 snew(simtemp->temperatures,n_lambda);
314 bDum=gmx_fio_ndo_real(fio,simtemp->temperatures,n_lambda);
319 static void do_fepvals(t_fileio *fio,t_lambda *fepvals,gmx_bool bRead, int file_version)
321 /* i is defined in the ndo_double macro; use g to iterate. */
327 /* free energy values */
328 if (file_version >= 79)
330 gmx_fio_do_int(fio,fepvals->init_fep_state);
331 gmx_fio_do_double(fio,fepvals->init_lambda);
332 gmx_fio_do_double(fio,fepvals->delta_lambda);
334 else if (file_version >= 59) {
335 gmx_fio_do_double(fio,fepvals->init_lambda);
336 gmx_fio_do_double(fio,fepvals->delta_lambda);
338 gmx_fio_do_real(fio,rdum);
339 fepvals->init_lambda = rdum;
340 gmx_fio_do_real(fio,rdum);
341 fepvals->delta_lambda = rdum;
343 if (file_version >= 79)
345 gmx_fio_do_int(fio,fepvals->n_lambda);
348 snew(fepvals->all_lambda,efptNR);
350 for (g=0;g<efptNR;g++)
352 if (fepvals->n_lambda > 0) {
355 snew(fepvals->all_lambda[g],fepvals->n_lambda);
357 bDum=gmx_fio_ndo_double(fio,fepvals->all_lambda[g],fepvals->n_lambda);
358 bDum=gmx_fio_ndo_int(fio,fepvals->separate_dvdl,efptNR);
360 else if (fepvals->init_lambda >= 0)
362 fepvals->separate_dvdl[efptFEP] = TRUE;
366 else if (file_version >= 64)
368 gmx_fio_do_int(fio,fepvals->n_lambda);
369 snew(fepvals->all_lambda,efptNR);
372 snew(fepvals->all_lambda[efptFEP],fepvals->n_lambda);
374 bDum=gmx_fio_ndo_double(fio,fepvals->all_lambda[efptFEP],fepvals->n_lambda);
375 if (fepvals->init_lambda >= 0)
377 fepvals->separate_dvdl[efptFEP] = TRUE;
382 fepvals->n_lambda = 0;
383 fepvals->all_lambda = NULL;
384 if (fepvals->init_lambda >= 0)
386 fepvals->separate_dvdl[efptFEP] = TRUE;
389 if (file_version >= 13)
391 gmx_fio_do_real(fio,fepvals->sc_alpha);
395 fepvals->sc_alpha = 0;
397 if (file_version >= 38)
399 gmx_fio_do_int(fio,fepvals->sc_power);
403 fepvals->sc_power = 2;
405 if (file_version >= 79)
407 gmx_fio_do_real(fio,fepvals->sc_r_power);
411 fepvals->sc_r_power = 6.0;
413 if (file_version >= 15)
415 gmx_fio_do_real(fio,fepvals->sc_sigma);
419 fepvals->sc_sigma = 0.3;
423 if (file_version >= 71)
425 fepvals->sc_sigma_min = fepvals->sc_sigma;
429 fepvals->sc_sigma_min = 0;
432 if (file_version >= 79)
434 gmx_fio_do_int(fio,fepvals->bScCoul);
438 fepvals->bScCoul = TRUE;
440 if (file_version >= 64) {
441 gmx_fio_do_int(fio,fepvals->nstdhdl);
443 fepvals->nstdhdl = 1;
446 if (file_version >= 73)
448 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
449 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
453 fepvals->separate_dhdl_file = esepdhdlfileYES;
454 fepvals->dhdl_derivatives = edhdlderivativesYES;
456 if (file_version >= 71)
458 gmx_fio_do_int(fio,fepvals->dh_hist_size);
459 gmx_fio_do_double(fio,fepvals->dh_hist_spacing);
463 fepvals->dh_hist_size = 0;
464 fepvals->dh_hist_spacing = 0.1;
466 if (file_version >= 79)
468 gmx_fio_do_int(fio,fepvals->bPrintEnergy);
472 fepvals->bPrintEnergy = FALSE;
476 static void do_pull(t_fileio *fio, t_pull *pull,gmx_bool bRead, int file_version)
480 gmx_fio_do_int(fio,pull->ngrp);
481 gmx_fio_do_int(fio,pull->eGeom);
482 gmx_fio_do_ivec(fio,pull->dim);
483 gmx_fio_do_real(fio,pull->cyl_r1);
484 gmx_fio_do_real(fio,pull->cyl_r0);
485 gmx_fio_do_real(fio,pull->constr_tol);
486 gmx_fio_do_int(fio,pull->nstxout);
487 gmx_fio_do_int(fio,pull->nstfout);
489 snew(pull->grp,pull->ngrp+1);
490 for(g=0; g<pull->ngrp+1; g++)
491 do_pullgrp(fio,&pull->grp[g],bRead,file_version);
495 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg,gmx_bool bRead, int file_version)
500 gmx_fio_do_int(fio,rotg->eType);
501 gmx_fio_do_int(fio,rotg->bMassW);
502 gmx_fio_do_int(fio,rotg->nat);
504 snew(rotg->ind,rotg->nat);
505 gmx_fio_ndo_int(fio,rotg->ind,rotg->nat);
507 snew(rotg->x_ref,rotg->nat);
508 gmx_fio_ndo_rvec(fio,rotg->x_ref,rotg->nat);
509 gmx_fio_do_rvec(fio,rotg->vec);
510 gmx_fio_do_rvec(fio,rotg->pivot);
511 gmx_fio_do_real(fio,rotg->rate);
512 gmx_fio_do_real(fio,rotg->k);
513 gmx_fio_do_real(fio,rotg->slab_dist);
514 gmx_fio_do_real(fio,rotg->min_gaussian);
515 gmx_fio_do_real(fio,rotg->eps);
516 gmx_fio_do_int(fio,rotg->eFittype);
517 gmx_fio_do_int(fio,rotg->PotAngle_nstep);
518 gmx_fio_do_real(fio,rotg->PotAngle_step);
521 static void do_rot(t_fileio *fio, t_rot *rot,gmx_bool bRead, int file_version)
525 gmx_fio_do_int(fio,rot->ngrp);
526 gmx_fio_do_int(fio,rot->nstrout);
527 gmx_fio_do_int(fio,rot->nstsout);
529 snew(rot->grp,rot->ngrp);
530 for(g=0; g<rot->ngrp; g++)
531 do_rotgrp(fio, &rot->grp[g],bRead,file_version);
535 static void do_inputrec(t_fileio *fio, t_inputrec *ir,gmx_bool bRead,
536 int file_version, real *fudgeQQ)
538 int i,j,k,*tmp,idum=0;
543 real zerotemptime,finish_t,init_temp,finish_temp;
545 if (file_version != tpx_version)
547 /* Give a warning about features that are not accessible */
548 fprintf(stderr,"Note: file tpx version %d, software tpx version %d\n",
549 file_version,tpx_version);
557 if (file_version == 0)
562 /* Basic inputrec stuff */
563 gmx_fio_do_int(fio,ir->eI);
564 if (file_version >= 62) {
565 gmx_fio_do_gmx_large_int(fio, ir->nsteps);
567 gmx_fio_do_int(fio,idum);
570 if(file_version > 25) {
571 if (file_version >= 62) {
572 gmx_fio_do_gmx_large_int(fio, ir->init_step);
574 gmx_fio_do_int(fio,idum);
575 ir->init_step = idum;
581 if(file_version >= 58)
582 gmx_fio_do_int(fio,ir->simulation_part);
584 ir->simulation_part=1;
586 if (file_version >= 67) {
587 gmx_fio_do_int(fio,ir->nstcalcenergy);
589 ir->nstcalcenergy = 1;
591 if (file_version < 53) {
592 /* The pbc info has been moved out of do_inputrec,
593 * since we always want it, also without reading the inputrec.
595 gmx_fio_do_int(fio,ir->ePBC);
596 if ((file_version <= 15) && (ir->ePBC == 2))
598 if (file_version >= 45) {
599 gmx_fio_do_int(fio,ir->bPeriodicMols);
603 ir->bPeriodicMols = TRUE;
605 ir->bPeriodicMols = FALSE;
609 gmx_fio_do_int(fio,ir->ns_type);
610 gmx_fio_do_int(fio,ir->nstlist);
611 gmx_fio_do_int(fio,ir->ndelta);
612 if (file_version < 41) {
613 gmx_fio_do_int(fio,idum);
614 gmx_fio_do_int(fio,idum);
616 if (file_version >= 45)
617 gmx_fio_do_real(fio,ir->rtpi);
620 gmx_fio_do_int(fio,ir->nstcomm);
621 if (file_version > 34)
622 gmx_fio_do_int(fio,ir->comm_mode);
623 else if (ir->nstcomm < 0)
624 ir->comm_mode = ecmANGULAR;
626 ir->comm_mode = ecmLINEAR;
627 ir->nstcomm = abs(ir->nstcomm);
629 if(file_version > 25)
630 gmx_fio_do_int(fio,ir->nstcheckpoint);
634 gmx_fio_do_int(fio,ir->nstcgsteep);
637 gmx_fio_do_int(fio,ir->nbfgscorr);
641 gmx_fio_do_int(fio,ir->nstlog);
642 gmx_fio_do_int(fio,ir->nstxout);
643 gmx_fio_do_int(fio,ir->nstvout);
644 gmx_fio_do_int(fio,ir->nstfout);
645 gmx_fio_do_int(fio,ir->nstenergy);
646 gmx_fio_do_int(fio,ir->nstxtcout);
647 if (file_version >= 59) {
648 gmx_fio_do_double(fio,ir->init_t);
649 gmx_fio_do_double(fio,ir->delta_t);
651 gmx_fio_do_real(fio,rdum);
653 gmx_fio_do_real(fio,rdum);
656 gmx_fio_do_real(fio,ir->xtcprec);
657 if (file_version < 19) {
658 gmx_fio_do_int(fio,idum);
659 gmx_fio_do_int(fio,idum);
661 if(file_version < 18)
662 gmx_fio_do_int(fio,idum);
663 gmx_fio_do_real(fio,ir->rlist);
664 if (file_version >= 67) {
665 gmx_fio_do_real(fio,ir->rlistlong);
667 gmx_fio_do_int(fio,ir->coulombtype);
668 if (file_version < 32 && ir->coulombtype == eelRF)
669 ir->coulombtype = eelRF_NEC;
670 gmx_fio_do_real(fio,ir->rcoulomb_switch);
671 gmx_fio_do_real(fio,ir->rcoulomb);
672 gmx_fio_do_int(fio,ir->vdwtype);
673 gmx_fio_do_real(fio,ir->rvdw_switch);
674 gmx_fio_do_real(fio,ir->rvdw);
675 if (file_version < 67) {
676 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
678 gmx_fio_do_int(fio,ir->eDispCorr);
679 gmx_fio_do_real(fio,ir->epsilon_r);
680 if (file_version >= 37) {
681 gmx_fio_do_real(fio,ir->epsilon_rf);
683 if (EEL_RF(ir->coulombtype)) {
684 ir->epsilon_rf = ir->epsilon_r;
687 ir->epsilon_rf = 1.0;
690 if (file_version >= 29)
691 gmx_fio_do_real(fio,ir->tabext);
695 if(file_version > 25) {
696 gmx_fio_do_int(fio,ir->gb_algorithm);
697 gmx_fio_do_int(fio,ir->nstgbradii);
698 gmx_fio_do_real(fio,ir->rgbradii);
699 gmx_fio_do_real(fio,ir->gb_saltconc);
700 gmx_fio_do_int(fio,ir->implicit_solvent);
702 ir->gb_algorithm=egbSTILL;
706 ir->implicit_solvent=eisNO;
710 gmx_fio_do_real(fio,ir->gb_epsilon_solvent);
711 gmx_fio_do_real(fio,ir->gb_obc_alpha);
712 gmx_fio_do_real(fio,ir->gb_obc_beta);
713 gmx_fio_do_real(fio,ir->gb_obc_gamma);
716 gmx_fio_do_real(fio,ir->gb_dielectric_offset);
717 gmx_fio_do_int(fio,ir->sa_algorithm);
721 ir->gb_dielectric_offset = 0.009;
722 ir->sa_algorithm = esaAPPROX;
724 gmx_fio_do_real(fio,ir->sa_surface_tension);
726 /* Override sa_surface_tension if it is not changed in the mpd-file */
727 if(ir->sa_surface_tension<0)
729 if(ir->gb_algorithm==egbSTILL)
731 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
733 else if(ir->gb_algorithm==egbHCT || ir->gb_algorithm==egbOBC)
735 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
742 /* Better use sensible values than insane (0.0) ones... */
743 ir->gb_epsilon_solvent = 80;
744 ir->gb_obc_alpha = 1.0;
745 ir->gb_obc_beta = 0.8;
746 ir->gb_obc_gamma = 4.85;
747 ir->sa_surface_tension = 2.092;
751 gmx_fio_do_int(fio,ir->nkx);
752 gmx_fio_do_int(fio,ir->nky);
753 gmx_fio_do_int(fio,ir->nkz);
754 gmx_fio_do_int(fio,ir->pme_order);
755 gmx_fio_do_real(fio,ir->ewald_rtol);
757 if (file_version >=24)
758 gmx_fio_do_int(fio,ir->ewald_geometry);
760 ir->ewald_geometry=eewg3D;
762 if (file_version <=17) {
763 ir->epsilon_surface=0;
764 if (file_version==17)
765 gmx_fio_do_int(fio,idum);
768 gmx_fio_do_real(fio,ir->epsilon_surface);
770 gmx_fio_do_gmx_bool(fio,ir->bOptFFT);
772 gmx_fio_do_gmx_bool(fio,ir->bContinuation);
773 gmx_fio_do_int(fio,ir->etc);
774 /* before version 18, ir->etc was a gmx_bool (ir->btc),
775 * but the values 0 and 1 still mean no and
776 * berendsen temperature coupling, respectively.
778 if (file_version >= 79) {
779 gmx_fio_do_gmx_bool(fio,ir->bPrintNHChains);
781 if (file_version >= 71)
783 gmx_fio_do_int(fio,ir->nsttcouple);
787 ir->nsttcouple = ir->nstcalcenergy;
789 if (file_version <= 15)
791 gmx_fio_do_int(fio,idum);
793 if (file_version <=17)
795 gmx_fio_do_int(fio,ir->epct);
796 if (file_version <= 15)
800 ir->epct = epctSURFACETENSION;
802 gmx_fio_do_int(fio,idum);
805 /* we have removed the NO alternative at the beginning */
809 ir->epct=epctISOTROPIC;
813 ir->epc=epcBERENDSEN;
818 gmx_fio_do_int(fio,ir->epc);
819 gmx_fio_do_int(fio,ir->epct);
821 if (file_version >= 71)
823 gmx_fio_do_int(fio,ir->nstpcouple);
827 ir->nstpcouple = ir->nstcalcenergy;
829 gmx_fio_do_real(fio,ir->tau_p);
830 if (file_version <= 15) {
831 gmx_fio_do_rvec(fio,vdum);
832 clear_mat(ir->ref_p);
834 ir->ref_p[i][i] = vdum[i];
836 gmx_fio_do_rvec(fio,ir->ref_p[XX]);
837 gmx_fio_do_rvec(fio,ir->ref_p[YY]);
838 gmx_fio_do_rvec(fio,ir->ref_p[ZZ]);
840 if (file_version <= 15) {
841 gmx_fio_do_rvec(fio,vdum);
842 clear_mat(ir->compress);
844 ir->compress[i][i] = vdum[i];
847 gmx_fio_do_rvec(fio,ir->compress[XX]);
848 gmx_fio_do_rvec(fio,ir->compress[YY]);
849 gmx_fio_do_rvec(fio,ir->compress[ZZ]);
851 if (file_version >= 47) {
852 gmx_fio_do_int(fio,ir->refcoord_scaling);
853 gmx_fio_do_rvec(fio,ir->posres_com);
854 gmx_fio_do_rvec(fio,ir->posres_comB);
856 ir->refcoord_scaling = erscNO;
857 clear_rvec(ir->posres_com);
858 clear_rvec(ir->posres_comB);
860 if((file_version > 25) && (file_version < 79))
861 gmx_fio_do_int(fio,ir->andersen_seed);
864 if(file_version < 26) {
865 gmx_fio_do_gmx_bool(fio,bSimAnn);
866 gmx_fio_do_real(fio,zerotemptime);
869 if (file_version < 37)
870 gmx_fio_do_real(fio,rdum);
872 gmx_fio_do_real(fio,ir->shake_tol);
873 if (file_version < 54)
874 gmx_fio_do_real(fio,*fudgeQQ);
876 gmx_fio_do_int(fio,ir->efep);
877 if (file_version <= 14 && ir->efep != efepNO)
881 do_fepvals(fio,ir->fepvals,bRead,file_version);
883 if (file_version >= 79)
885 gmx_fio_do_gmx_bool(fio,ir->bSimTemp);
893 ir->bSimTemp = FALSE;
897 do_simtempvals(fio,ir->simtempvals,ir->fepvals->n_lambda,bRead,file_version);
900 if (file_version >= 79)
902 gmx_fio_do_gmx_bool(fio,ir->bExpanded);
905 ir->bExpanded = TRUE;
909 ir->bExpanded = FALSE;
914 do_expandedvals(fio,ir->expandedvals,ir->fepvals->n_lambda,bRead,file_version);
916 if (file_version >= 57) {
917 gmx_fio_do_int(fio,ir->eDisre);
919 gmx_fio_do_int(fio,ir->eDisreWeighting);
920 if (file_version < 22) {
921 if (ir->eDisreWeighting == 0)
922 ir->eDisreWeighting = edrwEqual;
924 ir->eDisreWeighting = edrwConservative;
926 gmx_fio_do_gmx_bool(fio,ir->bDisreMixed);
927 gmx_fio_do_real(fio,ir->dr_fc);
928 gmx_fio_do_real(fio,ir->dr_tau);
929 gmx_fio_do_int(fio,ir->nstdisreout);
930 if (file_version >= 22) {
931 gmx_fio_do_real(fio,ir->orires_fc);
932 gmx_fio_do_real(fio,ir->orires_tau);
933 gmx_fio_do_int(fio,ir->nstorireout);
939 if(file_version >= 26 && file_version < 79) {
940 gmx_fio_do_real(fio,ir->dihre_fc);
941 if (file_version < 56)
943 gmx_fio_do_real(fio,rdum);
944 gmx_fio_do_int(fio,idum);
950 gmx_fio_do_real(fio,ir->em_stepsize);
951 gmx_fio_do_real(fio,ir->em_tol);
952 if (file_version >= 22)
953 gmx_fio_do_gmx_bool(fio,ir->bShakeSOR);
955 ir->bShakeSOR = TRUE;
956 if (file_version >= 11)
957 gmx_fio_do_int(fio,ir->niter);
960 fprintf(stderr,"Note: niter not in run input file, setting it to %d\n",
963 if (file_version >= 21)
964 gmx_fio_do_real(fio,ir->fc_stepsize);
967 gmx_fio_do_int(fio,ir->eConstrAlg);
968 gmx_fio_do_int(fio,ir->nProjOrder);
969 gmx_fio_do_real(fio,ir->LincsWarnAngle);
970 if (file_version <= 14)
971 gmx_fio_do_int(fio,idum);
972 if (file_version >=26)
973 gmx_fio_do_int(fio,ir->nLincsIter);
976 fprintf(stderr,"Note: nLincsIter not in run input file, setting it to %d\n",
979 if (file_version < 33)
980 gmx_fio_do_real(fio,bd_temp);
981 gmx_fio_do_real(fio,ir->bd_fric);
982 gmx_fio_do_int(fio,ir->ld_seed);
983 if (file_version >= 33) {
985 gmx_fio_do_rvec(fio,ir->deform[i]);
988 clear_rvec(ir->deform[i]);
990 if (file_version >= 14)
991 gmx_fio_do_real(fio,ir->cos_accel);
994 gmx_fio_do_int(fio,ir->userint1);
995 gmx_fio_do_int(fio,ir->userint2);
996 gmx_fio_do_int(fio,ir->userint3);
997 gmx_fio_do_int(fio,ir->userint4);
998 gmx_fio_do_real(fio,ir->userreal1);
999 gmx_fio_do_real(fio,ir->userreal2);
1000 gmx_fio_do_real(fio,ir->userreal3);
1001 gmx_fio_do_real(fio,ir->userreal4);
1004 if (file_version >= 77) {
1005 gmx_fio_do_gmx_bool(fio,ir->bAdress);
1007 if (bRead) snew(ir->adress, 1);
1008 gmx_fio_do_int(fio,ir->adress->type);
1009 gmx_fio_do_real(fio,ir->adress->const_wf);
1010 gmx_fio_do_real(fio,ir->adress->ex_width);
1011 gmx_fio_do_real(fio,ir->adress->hy_width);
1012 gmx_fio_do_int(fio,ir->adress->icor);
1013 gmx_fio_do_int(fio,ir->adress->site);
1014 gmx_fio_do_rvec(fio,ir->adress->refs);
1015 gmx_fio_do_int(fio,ir->adress->n_tf_grps);
1016 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1017 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1018 gmx_fio_do_int(fio,ir->adress->do_hybridpairs);
1020 if (bRead)snew(ir->adress->tf_table_index,ir->adress->n_tf_grps);
1021 if (ir->adress->n_tf_grps > 0) {
1022 bDum=gmx_fio_ndo_int(fio,ir->adress->tf_table_index,ir->adress->n_tf_grps);
1024 if (bRead)snew(ir->adress->group_explicit,ir->adress->n_energy_grps);
1025 if (ir->adress->n_energy_grps > 0) {
1026 bDum=gmx_fio_ndo_int(fio, ir->adress->group_explicit,ir->adress->n_energy_grps);
1030 ir->bAdress = FALSE;
1034 if (file_version >= 48) {
1035 gmx_fio_do_int(fio,ir->ePull);
1036 if (ir->ePull != epullNO) {
1039 do_pull(fio, ir->pull,bRead,file_version);
1042 ir->ePull = epullNO;
1045 /* Enforced rotation */
1046 if (file_version >= 74) {
1047 gmx_fio_do_int(fio,ir->bRot);
1048 if (ir->bRot == TRUE) {
1051 do_rot(fio, ir->rot,bRead,file_version);
1058 gmx_fio_do_int(fio,ir->opts.ngtc);
1059 if (file_version >= 69) {
1060 gmx_fio_do_int(fio,ir->opts.nhchainlength);
1062 ir->opts.nhchainlength = 1;
1064 gmx_fio_do_int(fio,ir->opts.ngacc);
1065 gmx_fio_do_int(fio,ir->opts.ngfrz);
1066 gmx_fio_do_int(fio,ir->opts.ngener);
1069 snew(ir->opts.nrdf, ir->opts.ngtc);
1070 snew(ir->opts.ref_t, ir->opts.ngtc);
1071 snew(ir->opts.annealing, ir->opts.ngtc);
1072 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1073 snew(ir->opts.anneal_time, ir->opts.ngtc);
1074 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1075 snew(ir->opts.tau_t, ir->opts.ngtc);
1076 snew(ir->opts.nFreeze,ir->opts.ngfrz);
1077 snew(ir->opts.acc, ir->opts.ngacc);
1078 snew(ir->opts.egp_flags,ir->opts.ngener*ir->opts.ngener);
1080 if (ir->opts.ngtc > 0) {
1081 if (bRead && file_version<13) {
1082 snew(tmp,ir->opts.ngtc);
1083 bDum=gmx_fio_ndo_int(fio,tmp, ir->opts.ngtc);
1084 for(i=0; i<ir->opts.ngtc; i++)
1085 ir->opts.nrdf[i] = tmp[i];
1088 bDum=gmx_fio_ndo_real(fio,ir->opts.nrdf, ir->opts.ngtc);
1090 bDum=gmx_fio_ndo_real(fio,ir->opts.ref_t,ir->opts.ngtc);
1091 bDum=gmx_fio_ndo_real(fio,ir->opts.tau_t,ir->opts.ngtc);
1092 if (file_version<33 && ir->eI==eiBD) {
1093 for(i=0; i<ir->opts.ngtc; i++)
1094 ir->opts.tau_t[i] = bd_temp;
1097 if (ir->opts.ngfrz > 0)
1098 bDum=gmx_fio_ndo_ivec(fio,ir->opts.nFreeze,ir->opts.ngfrz);
1099 if (ir->opts.ngacc > 0)
1100 gmx_fio_ndo_rvec(fio,ir->opts.acc,ir->opts.ngacc);
1101 if (file_version >= 12)
1102 bDum=gmx_fio_ndo_int(fio,ir->opts.egp_flags,
1103 ir->opts.ngener*ir->opts.ngener);
1105 if(bRead && file_version < 26) {
1106 for(i=0;i<ir->opts.ngtc;i++) {
1108 ir->opts.annealing[i] = eannSINGLE;
1109 ir->opts.anneal_npoints[i] = 2;
1110 snew(ir->opts.anneal_time[i],2);
1111 snew(ir->opts.anneal_temp[i],2);
1112 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1113 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1114 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1115 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1116 ir->opts.anneal_time[i][0] = ir->init_t;
1117 ir->opts.anneal_time[i][1] = finish_t;
1118 ir->opts.anneal_temp[i][0] = init_temp;
1119 ir->opts.anneal_temp[i][1] = finish_temp;
1121 ir->opts.annealing[i] = eannNO;
1122 ir->opts.anneal_npoints[i] = 0;
1126 /* file version 26 or later */
1127 /* First read the lists with annealing and npoints for each group */
1128 bDum=gmx_fio_ndo_int(fio,ir->opts.annealing,ir->opts.ngtc);
1129 bDum=gmx_fio_ndo_int(fio,ir->opts.anneal_npoints,ir->opts.ngtc);
1130 for(j=0;j<(ir->opts.ngtc);j++) {
1131 k=ir->opts.anneal_npoints[j];
1133 snew(ir->opts.anneal_time[j],k);
1134 snew(ir->opts.anneal_temp[j],k);
1136 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_time[j],k);
1137 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_temp[j],k);
1141 if (file_version >= 45) {
1142 gmx_fio_do_int(fio,ir->nwall);
1143 gmx_fio_do_int(fio,ir->wall_type);
1144 if (file_version >= 50)
1145 gmx_fio_do_real(fio,ir->wall_r_linpot);
1147 ir->wall_r_linpot = -1;
1148 gmx_fio_do_int(fio,ir->wall_atomtype[0]);
1149 gmx_fio_do_int(fio,ir->wall_atomtype[1]);
1150 gmx_fio_do_real(fio,ir->wall_density[0]);
1151 gmx_fio_do_real(fio,ir->wall_density[1]);
1152 gmx_fio_do_real(fio,ir->wall_ewald_zfac);
1156 ir->wall_atomtype[0] = -1;
1157 ir->wall_atomtype[1] = -1;
1158 ir->wall_density[0] = 0;
1159 ir->wall_density[1] = 0;
1160 ir->wall_ewald_zfac = 3;
1162 /* Cosine stuff for electric fields */
1163 for(j=0; (j<DIM); j++) {
1164 gmx_fio_do_int(fio,ir->ex[j].n);
1165 gmx_fio_do_int(fio,ir->et[j].n);
1167 snew(ir->ex[j].a, ir->ex[j].n);
1168 snew(ir->ex[j].phi,ir->ex[j].n);
1169 snew(ir->et[j].a, ir->et[j].n);
1170 snew(ir->et[j].phi,ir->et[j].n);
1172 bDum=gmx_fio_ndo_real(fio,ir->ex[j].a, ir->ex[j].n);
1173 bDum=gmx_fio_ndo_real(fio,ir->ex[j].phi,ir->ex[j].n);
1174 bDum=gmx_fio_ndo_real(fio,ir->et[j].a, ir->et[j].n);
1175 bDum=gmx_fio_ndo_real(fio,ir->et[j].phi,ir->et[j].n);
1179 if(file_version>=39){
1180 gmx_fio_do_gmx_bool(fio,ir->bQMMM);
1181 gmx_fio_do_int(fio,ir->QMMMscheme);
1182 gmx_fio_do_real(fio,ir->scalefactor);
1183 gmx_fio_do_int(fio,ir->opts.ngQM);
1185 snew(ir->opts.QMmethod, ir->opts.ngQM);
1186 snew(ir->opts.QMbasis, ir->opts.ngQM);
1187 snew(ir->opts.QMcharge, ir->opts.ngQM);
1188 snew(ir->opts.QMmult, ir->opts.ngQM);
1189 snew(ir->opts.bSH, ir->opts.ngQM);
1190 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1191 snew(ir->opts.CASelectrons,ir->opts.ngQM);
1192 snew(ir->opts.SAon, ir->opts.ngQM);
1193 snew(ir->opts.SAoff, ir->opts.ngQM);
1194 snew(ir->opts.SAsteps, ir->opts.ngQM);
1195 snew(ir->opts.bOPT, ir->opts.ngQM);
1196 snew(ir->opts.bTS, ir->opts.ngQM);
1198 if (ir->opts.ngQM > 0) {
1199 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmethod,ir->opts.ngQM);
1200 bDum=gmx_fio_ndo_int(fio,ir->opts.QMbasis,ir->opts.ngQM);
1201 bDum=gmx_fio_ndo_int(fio,ir->opts.QMcharge,ir->opts.ngQM);
1202 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmult,ir->opts.ngQM);
1203 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bSH,ir->opts.ngQM);
1204 bDum=gmx_fio_ndo_int(fio,ir->opts.CASorbitals,ir->opts.ngQM);
1205 bDum=gmx_fio_ndo_int(fio,ir->opts.CASelectrons,ir->opts.ngQM);
1206 bDum=gmx_fio_ndo_real(fio,ir->opts.SAon,ir->opts.ngQM);
1207 bDum=gmx_fio_ndo_real(fio,ir->opts.SAoff,ir->opts.ngQM);
1208 bDum=gmx_fio_ndo_int(fio,ir->opts.SAsteps,ir->opts.ngQM);
1209 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bOPT,ir->opts.ngQM);
1210 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bTS,ir->opts.ngQM);
1212 /* end of QMMM stuff */
1217 static void do_harm(t_fileio *fio, t_iparams *iparams,gmx_bool bRead)
1219 gmx_fio_do_real(fio,iparams->harmonic.rA);
1220 gmx_fio_do_real(fio,iparams->harmonic.krA);
1221 gmx_fio_do_real(fio,iparams->harmonic.rB);
1222 gmx_fio_do_real(fio,iparams->harmonic.krB);
1225 void do_iparams(t_fileio *fio, t_functype ftype,t_iparams *iparams,
1226 gmx_bool bRead, int file_version)
1233 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1241 do_harm(fio, iparams,bRead);
1242 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead) {
1243 /* Correct incorrect storage of parameters */
1244 iparams->pdihs.phiB = iparams->pdihs.phiA;
1245 iparams->pdihs.cpB = iparams->pdihs.cpA;
1248 case F_LINEAR_ANGLES:
1249 gmx_fio_do_real(fio,iparams->linangle.klinA);
1250 gmx_fio_do_real(fio,iparams->linangle.aA);
1251 gmx_fio_do_real(fio,iparams->linangle.klinB);
1252 gmx_fio_do_real(fio,iparams->linangle.aB);
1255 gmx_fio_do_real(fio,iparams->fene.bm);
1256 gmx_fio_do_real(fio,iparams->fene.kb);
1259 gmx_fio_do_real(fio,iparams->restraint.lowA);
1260 gmx_fio_do_real(fio,iparams->restraint.up1A);
1261 gmx_fio_do_real(fio,iparams->restraint.up2A);
1262 gmx_fio_do_real(fio,iparams->restraint.kA);
1263 gmx_fio_do_real(fio,iparams->restraint.lowB);
1264 gmx_fio_do_real(fio,iparams->restraint.up1B);
1265 gmx_fio_do_real(fio,iparams->restraint.up2B);
1266 gmx_fio_do_real(fio,iparams->restraint.kB);
1272 gmx_fio_do_real(fio,iparams->tab.kA);
1273 gmx_fio_do_int(fio,iparams->tab.table);
1274 gmx_fio_do_real(fio,iparams->tab.kB);
1276 case F_CROSS_BOND_BONDS:
1277 gmx_fio_do_real(fio,iparams->cross_bb.r1e);
1278 gmx_fio_do_real(fio,iparams->cross_bb.r2e);
1279 gmx_fio_do_real(fio,iparams->cross_bb.krr);
1281 case F_CROSS_BOND_ANGLES:
1282 gmx_fio_do_real(fio,iparams->cross_ba.r1e);
1283 gmx_fio_do_real(fio,iparams->cross_ba.r2e);
1284 gmx_fio_do_real(fio,iparams->cross_ba.r3e);
1285 gmx_fio_do_real(fio,iparams->cross_ba.krt);
1287 case F_UREY_BRADLEY:
1288 gmx_fio_do_real(fio,iparams->u_b.thetaA);
1289 gmx_fio_do_real(fio,iparams->u_b.kthetaA);
1290 gmx_fio_do_real(fio,iparams->u_b.r13A);
1291 gmx_fio_do_real(fio,iparams->u_b.kUBA);
1292 if (file_version >= 79) {
1293 gmx_fio_do_real(fio,iparams->u_b.thetaB);
1294 gmx_fio_do_real(fio,iparams->u_b.kthetaB);
1295 gmx_fio_do_real(fio,iparams->u_b.r13B);
1296 gmx_fio_do_real(fio,iparams->u_b.kUBB);
1298 iparams->u_b.thetaB=iparams->u_b.thetaA;
1299 iparams->u_b.kthetaB=iparams->u_b.kthetaA;
1300 iparams->u_b.r13B=iparams->u_b.r13A;
1301 iparams->u_b.kUBB=iparams->u_b.kUBA;
1304 case F_QUARTIC_ANGLES:
1305 gmx_fio_do_real(fio,iparams->qangle.theta);
1306 bDum=gmx_fio_ndo_real(fio,iparams->qangle.c,5);
1309 gmx_fio_do_real(fio,iparams->bham.a);
1310 gmx_fio_do_real(fio,iparams->bham.b);
1311 gmx_fio_do_real(fio,iparams->bham.c);
1314 gmx_fio_do_real(fio,iparams->morse.b0A);
1315 gmx_fio_do_real(fio,iparams->morse.cbA);
1316 gmx_fio_do_real(fio,iparams->morse.betaA);
1317 if (file_version >= 79) {
1318 gmx_fio_do_real(fio,iparams->morse.b0B);
1319 gmx_fio_do_real(fio,iparams->morse.cbB);
1320 gmx_fio_do_real(fio,iparams->morse.betaB);
1322 iparams->morse.b0B = iparams->morse.b0A;
1323 iparams->morse.cbB = iparams->morse.cbA;
1324 iparams->morse.betaB = iparams->morse.betaA;
1328 gmx_fio_do_real(fio,iparams->cubic.b0);
1329 gmx_fio_do_real(fio,iparams->cubic.kb);
1330 gmx_fio_do_real(fio,iparams->cubic.kcub);
1334 case F_POLARIZATION:
1335 gmx_fio_do_real(fio,iparams->polarize.alpha);
1338 gmx_fio_do_real(fio,iparams->anharm_polarize.alpha);
1339 gmx_fio_do_real(fio,iparams->anharm_polarize.drcut);
1340 gmx_fio_do_real(fio,iparams->anharm_polarize.khyp);
1343 if (file_version < 31)
1344 gmx_fatal(FARGS,"Old tpr files with water_polarization not supported. Make a new.");
1345 gmx_fio_do_real(fio,iparams->wpol.al_x);
1346 gmx_fio_do_real(fio,iparams->wpol.al_y);
1347 gmx_fio_do_real(fio,iparams->wpol.al_z);
1348 gmx_fio_do_real(fio,iparams->wpol.rOH);
1349 gmx_fio_do_real(fio,iparams->wpol.rHH);
1350 gmx_fio_do_real(fio,iparams->wpol.rOD);
1353 gmx_fio_do_real(fio,iparams->thole.a);
1354 gmx_fio_do_real(fio,iparams->thole.alpha1);
1355 gmx_fio_do_real(fio,iparams->thole.alpha2);
1356 gmx_fio_do_real(fio,iparams->thole.rfac);
1359 gmx_fio_do_real(fio,iparams->lj.c6);
1360 gmx_fio_do_real(fio,iparams->lj.c12);
1363 gmx_fio_do_real(fio,iparams->lj14.c6A);
1364 gmx_fio_do_real(fio,iparams->lj14.c12A);
1365 gmx_fio_do_real(fio,iparams->lj14.c6B);
1366 gmx_fio_do_real(fio,iparams->lj14.c12B);
1369 gmx_fio_do_real(fio,iparams->ljc14.fqq);
1370 gmx_fio_do_real(fio,iparams->ljc14.qi);
1371 gmx_fio_do_real(fio,iparams->ljc14.qj);
1372 gmx_fio_do_real(fio,iparams->ljc14.c6);
1373 gmx_fio_do_real(fio,iparams->ljc14.c12);
1375 case F_LJC_PAIRS_NB:
1376 gmx_fio_do_real(fio,iparams->ljcnb.qi);
1377 gmx_fio_do_real(fio,iparams->ljcnb.qj);
1378 gmx_fio_do_real(fio,iparams->ljcnb.c6);
1379 gmx_fio_do_real(fio,iparams->ljcnb.c12);
1385 gmx_fio_do_real(fio,iparams->pdihs.phiA);
1386 gmx_fio_do_real(fio,iparams->pdihs.cpA);
1387 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42) {
1388 /* Read the incorrectly stored multiplicity */
1389 gmx_fio_do_real(fio,iparams->harmonic.rB);
1390 gmx_fio_do_real(fio,iparams->harmonic.krB);
1391 iparams->pdihs.phiB = iparams->pdihs.phiA;
1392 iparams->pdihs.cpB = iparams->pdihs.cpA;
1394 gmx_fio_do_real(fio,iparams->pdihs.phiB);
1395 gmx_fio_do_real(fio,iparams->pdihs.cpB);
1396 gmx_fio_do_int(fio,iparams->pdihs.mult);
1400 gmx_fio_do_int(fio,iparams->disres.label);
1401 gmx_fio_do_int(fio,iparams->disres.type);
1402 gmx_fio_do_real(fio,iparams->disres.low);
1403 gmx_fio_do_real(fio,iparams->disres.up1);
1404 gmx_fio_do_real(fio,iparams->disres.up2);
1405 gmx_fio_do_real(fio,iparams->disres.kfac);
1408 gmx_fio_do_int(fio,iparams->orires.ex);
1409 gmx_fio_do_int(fio,iparams->orires.label);
1410 gmx_fio_do_int(fio,iparams->orires.power);
1411 gmx_fio_do_real(fio,iparams->orires.c);
1412 gmx_fio_do_real(fio,iparams->orires.obs);
1413 gmx_fio_do_real(fio,iparams->orires.kfac);
1416 gmx_fio_do_real(fio,iparams->dihres.phiA);
1417 gmx_fio_do_real(fio,iparams->dihres.dphiA);
1418 gmx_fio_do_real(fio,iparams->dihres.kfacA);
1419 if (file_version >= 72) {
1420 gmx_fio_do_real(fio,iparams->dihres.phiB);
1421 gmx_fio_do_real(fio,iparams->dihres.dphiB);
1422 gmx_fio_do_real(fio,iparams->dihres.kfacB);
1424 iparams->dihres.phiB=iparams->dihres.phiA;
1425 iparams->dihres.dphiB=iparams->dihres.dphiA;
1426 iparams->dihres.kfacB=iparams->dihres.kfacA;
1430 gmx_fio_do_rvec(fio,iparams->posres.pos0A);
1431 gmx_fio_do_rvec(fio,iparams->posres.fcA);
1432 if (bRead && file_version < 27) {
1433 copy_rvec(iparams->posres.pos0A,iparams->posres.pos0B);
1434 copy_rvec(iparams->posres.fcA,iparams->posres.fcB);
1436 gmx_fio_do_rvec(fio,iparams->posres.pos0B);
1437 gmx_fio_do_rvec(fio,iparams->posres.fcB);
1441 gmx_fio_do_int(fio,iparams->fbposres.geom);
1442 gmx_fio_do_rvec(fio,iparams->fbposres.pos0);
1443 gmx_fio_do_real(fio,iparams->fbposres.r);
1444 gmx_fio_do_real(fio,iparams->fbposres.k);
1447 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA,NR_RBDIHS);
1448 if(file_version>=25)
1449 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB,NR_RBDIHS);
1452 /* Fourier dihedrals are internally represented
1453 * as Ryckaert-Bellemans since those are faster to compute.
1455 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA, NR_RBDIHS);
1456 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB, NR_RBDIHS);
1460 gmx_fio_do_real(fio,iparams->constr.dA);
1461 gmx_fio_do_real(fio,iparams->constr.dB);
1464 gmx_fio_do_real(fio,iparams->settle.doh);
1465 gmx_fio_do_real(fio,iparams->settle.dhh);
1468 gmx_fio_do_real(fio,iparams->vsite.a);
1473 gmx_fio_do_real(fio,iparams->vsite.a);
1474 gmx_fio_do_real(fio,iparams->vsite.b);
1479 gmx_fio_do_real(fio,iparams->vsite.a);
1480 gmx_fio_do_real(fio,iparams->vsite.b);
1481 gmx_fio_do_real(fio,iparams->vsite.c);
1484 gmx_fio_do_int(fio,iparams->vsiten.n);
1485 gmx_fio_do_real(fio,iparams->vsiten.a);
1490 /* We got rid of some parameters in version 68 */
1491 if(bRead && file_version<68)
1493 gmx_fio_do_real(fio,rdum);
1494 gmx_fio_do_real(fio,rdum);
1495 gmx_fio_do_real(fio,rdum);
1496 gmx_fio_do_real(fio,rdum);
1498 gmx_fio_do_real(fio,iparams->gb.sar);
1499 gmx_fio_do_real(fio,iparams->gb.st);
1500 gmx_fio_do_real(fio,iparams->gb.pi);
1501 gmx_fio_do_real(fio,iparams->gb.gbr);
1502 gmx_fio_do_real(fio,iparams->gb.bmlt);
1505 gmx_fio_do_int(fio,iparams->cmap.cmapA);
1506 gmx_fio_do_int(fio,iparams->cmap.cmapB);
1509 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
1510 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
1513 gmx_fio_unset_comment(fio);
1516 static void do_ilist(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,int file_version,
1523 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1525 if (file_version < 44) {
1526 for(i=0; i<MAXNODES; i++)
1527 gmx_fio_do_int(fio,idum);
1529 gmx_fio_do_int(fio,ilist->nr);
1531 snew(ilist->iatoms,ilist->nr);
1532 bDum=gmx_fio_ndo_int(fio,ilist->iatoms,ilist->nr);
1534 gmx_fio_unset_comment(fio);
1537 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1538 gmx_bool bRead, int file_version)
1544 gmx_fio_do_int(fio,ffparams->atnr);
1545 if (file_version < 57) {
1546 gmx_fio_do_int(fio,idum);
1548 gmx_fio_do_int(fio,ffparams->ntypes);
1550 fprintf(debug,"ffparams->atnr = %d, ntypes = %d\n",
1551 ffparams->atnr,ffparams->ntypes);
1553 snew(ffparams->functype,ffparams->ntypes);
1554 snew(ffparams->iparams,ffparams->ntypes);
1556 /* Read/write all the function types */
1557 bDum=gmx_fio_ndo_int(fio,ffparams->functype,ffparams->ntypes);
1559 pr_ivec(debug,0,"functype",ffparams->functype,ffparams->ntypes,TRUE);
1561 if (file_version >= 66) {
1562 gmx_fio_do_double(fio,ffparams->reppow);
1564 ffparams->reppow = 12.0;
1567 if (file_version >= 57) {
1568 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1571 /* Check whether all these function types are supported by the code.
1572 * In practice the code is backwards compatible, which means that the
1573 * numbering may have to be altered from old numbering to new numbering
1575 for (i=0; (i<ffparams->ntypes); i++) {
1577 /* Loop over file versions */
1578 for (k=0; (k<NFTUPD); k++)
1579 /* Compare the read file_version to the update table */
1580 if ((file_version < ftupd[k].fvnr) &&
1581 (ffparams->functype[i] >= ftupd[k].ftype)) {
1582 ffparams->functype[i] += 1;
1584 fprintf(debug,"Incrementing function type %d to %d (due to %s)\n",
1585 i,ffparams->functype[i],
1586 interaction_function[ftupd[k].ftype].longname);
1591 do_iparams(fio, ffparams->functype[i],&ffparams->iparams[i],bRead,
1594 pr_iparams(debug,ffparams->functype[i],&ffparams->iparams[i]);
1598 static void add_settle_atoms(t_ilist *ilist)
1602 /* Settle used to only store the first atom: add the other two */
1603 srenew(ilist->iatoms,2*ilist->nr);
1604 for(i=ilist->nr/2-1; i>=0; i--)
1606 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
1607 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
1608 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
1609 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
1611 ilist->nr = 2*ilist->nr;
1614 static void do_ilists(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,
1617 int i,j,renum[F_NRE];
1618 gmx_bool bDum=TRUE,bClear;
1621 for(j=0; (j<F_NRE); j++) {
1624 for (k=0; k<NFTUPD; k++)
1625 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
1629 ilist[j].iatoms = NULL;
1631 do_ilist(fio, &ilist[j],bRead,file_version,j);
1632 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
1634 add_settle_atoms(&ilist[j]);
1638 if (bRead && gmx_debug_at)
1639 pr_ilist(debug,0,interaction_function[j].longname,
1640 functype,&ilist[j],TRUE);
1645 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams,gmx_moltype_t *molt,
1646 gmx_bool bRead, int file_version)
1648 do_ffparams(fio, ffparams,bRead,file_version);
1650 if (file_version >= 54) {
1651 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1654 do_ilists(fio, molt->ilist,bRead,file_version);
1657 static void do_block(t_fileio *fio, t_block *block,gmx_bool bRead,int file_version)
1659 int i,idum,dum_nra,*dum_a;
1662 if (file_version < 44)
1663 for(i=0; i<MAXNODES; i++)
1664 gmx_fio_do_int(fio,idum);
1665 gmx_fio_do_int(fio,block->nr);
1666 if (file_version < 51)
1667 gmx_fio_do_int(fio,dum_nra);
1669 block->nalloc_index = block->nr+1;
1670 snew(block->index,block->nalloc_index);
1672 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1674 if (file_version < 51 && dum_nra > 0) {
1675 snew(dum_a,dum_nra);
1676 bDum=gmx_fio_ndo_int(fio,dum_a,dum_nra);
1681 static void do_blocka(t_fileio *fio, t_blocka *block,gmx_bool bRead,
1687 if (file_version < 44)
1688 for(i=0; i<MAXNODES; i++)
1689 gmx_fio_do_int(fio,idum);
1690 gmx_fio_do_int(fio,block->nr);
1691 gmx_fio_do_int(fio,block->nra);
1693 block->nalloc_index = block->nr+1;
1694 snew(block->index,block->nalloc_index);
1695 block->nalloc_a = block->nra;
1696 snew(block->a,block->nalloc_a);
1698 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1699 bDum=gmx_fio_ndo_int(fio,block->a,block->nra);
1702 static void do_atom(t_fileio *fio, t_atom *atom,int ngrp,gmx_bool bRead,
1703 int file_version, gmx_groups_t *groups,int atnr)
1707 gmx_fio_do_real(fio,atom->m);
1708 gmx_fio_do_real(fio,atom->q);
1709 gmx_fio_do_real(fio,atom->mB);
1710 gmx_fio_do_real(fio,atom->qB);
1711 gmx_fio_do_ushort(fio, atom->type);
1712 gmx_fio_do_ushort(fio, atom->typeB);
1713 gmx_fio_do_int(fio,atom->ptype);
1714 gmx_fio_do_int(fio,atom->resind);
1715 if (file_version >= 52)
1716 gmx_fio_do_int(fio,atom->atomnumber);
1718 atom->atomnumber = NOTSET;
1719 if (file_version < 23)
1721 else if (file_version < 39)
1726 if (file_version < 57) {
1727 unsigned char uchar[egcNR];
1728 gmx_fio_ndo_uchar(fio,uchar,myngrp);
1729 for(i=myngrp; (i<ngrp); i++) {
1732 /* Copy the old data format to the groups struct */
1733 for(i=0; i<ngrp; i++) {
1734 groups->grpnr[i][atnr] = uchar[i];
1739 static void do_grps(t_fileio *fio, int ngrp,t_grps grps[],gmx_bool bRead,
1745 if (file_version < 23)
1747 else if (file_version < 39)
1752 for(j=0; (j<ngrp); j++) {
1754 gmx_fio_do_int(fio,grps[j].nr);
1756 snew(grps[j].nm_ind,grps[j].nr);
1757 bDum=gmx_fio_ndo_int(fio,grps[j].nm_ind,grps[j].nr);
1761 snew(grps[j].nm_ind,grps[j].nr);
1766 static void do_symstr(t_fileio *fio, char ***nm,gmx_bool bRead,t_symtab *symtab)
1771 gmx_fio_do_int(fio,ls);
1772 *nm = get_symtab_handle(symtab,ls);
1775 ls = lookup_symtab(symtab,*nm);
1776 gmx_fio_do_int(fio,ls);
1780 static void do_strstr(t_fileio *fio, int nstr,char ***nm,gmx_bool bRead,
1785 for (j=0; (j<nstr); j++)
1786 do_symstr(fio, &(nm[j]),bRead,symtab);
1789 static void do_resinfo(t_fileio *fio, int n,t_resinfo *ri,gmx_bool bRead,
1790 t_symtab *symtab, int file_version)
1794 for (j=0; (j<n); j++) {
1795 do_symstr(fio, &(ri[j].name),bRead,symtab);
1796 if (file_version >= 63) {
1797 gmx_fio_do_int(fio,ri[j].nr);
1798 gmx_fio_do_uchar(fio, ri[j].ic);
1806 static void do_atoms(t_fileio *fio, t_atoms *atoms,gmx_bool bRead,t_symtab *symtab,
1808 gmx_groups_t *groups)
1812 gmx_fio_do_int(fio,atoms->nr);
1813 gmx_fio_do_int(fio,atoms->nres);
1814 if (file_version < 57) {
1815 gmx_fio_do_int(fio,groups->ngrpname);
1816 for(i=0; i<egcNR; i++) {
1817 groups->ngrpnr[i] = atoms->nr;
1818 snew(groups->grpnr[i],groups->ngrpnr[i]);
1822 snew(atoms->atom,atoms->nr);
1823 snew(atoms->atomname,atoms->nr);
1824 snew(atoms->atomtype,atoms->nr);
1825 snew(atoms->atomtypeB,atoms->nr);
1826 snew(atoms->resinfo,atoms->nres);
1827 if (file_version < 57) {
1828 snew(groups->grpname,groups->ngrpname);
1830 atoms->pdbinfo = NULL;
1832 for(i=0; (i<atoms->nr); i++) {
1833 do_atom(fio, &atoms->atom[i],egcNR,bRead, file_version,groups,i);
1835 do_strstr(fio, atoms->nr,atoms->atomname,bRead,symtab);
1836 if (bRead && (file_version <= 20)) {
1837 for(i=0; i<atoms->nr; i++) {
1838 atoms->atomtype[i] = put_symtab(symtab,"?");
1839 atoms->atomtypeB[i] = put_symtab(symtab,"?");
1842 do_strstr(fio, atoms->nr,atoms->atomtype,bRead,symtab);
1843 do_strstr(fio, atoms->nr,atoms->atomtypeB,bRead,symtab);
1845 do_resinfo(fio, atoms->nres,atoms->resinfo,bRead,symtab,file_version);
1847 if (file_version < 57) {
1848 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1850 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1854 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
1855 gmx_bool bRead,t_symtab *symtab,
1861 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1862 gmx_fio_do_int(fio,groups->ngrpname);
1864 snew(groups->grpname,groups->ngrpname);
1866 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1867 for(g=0; g<egcNR; g++) {
1868 gmx_fio_do_int(fio,groups->ngrpnr[g]);
1869 if (groups->ngrpnr[g] == 0) {
1871 groups->grpnr[g] = NULL;
1875 snew(groups->grpnr[g],groups->ngrpnr[g]);
1877 bDum=gmx_fio_ndo_uchar(fio, groups->grpnr[g],groups->ngrpnr[g]);
1882 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes,gmx_bool bRead,
1883 t_symtab *symtab,int file_version)
1886 gmx_bool bDum = TRUE;
1888 if (file_version > 25) {
1889 gmx_fio_do_int(fio,atomtypes->nr);
1892 snew(atomtypes->radius,j);
1893 snew(atomtypes->vol,j);
1894 snew(atomtypes->surftens,j);
1895 snew(atomtypes->atomnumber,j);
1896 snew(atomtypes->gb_radius,j);
1897 snew(atomtypes->S_hct,j);
1899 bDum=gmx_fio_ndo_real(fio,atomtypes->radius,j);
1900 bDum=gmx_fio_ndo_real(fio,atomtypes->vol,j);
1901 bDum=gmx_fio_ndo_real(fio,atomtypes->surftens,j);
1902 if(file_version >= 40)
1904 bDum=gmx_fio_ndo_int(fio,atomtypes->atomnumber,j);
1906 if(file_version >= 60)
1908 bDum=gmx_fio_ndo_real(fio,atomtypes->gb_radius,j);
1909 bDum=gmx_fio_ndo_real(fio,atomtypes->S_hct,j);
1912 /* File versions prior to 26 cannot do GBSA,
1913 * so they dont use this structure
1916 atomtypes->radius = NULL;
1917 atomtypes->vol = NULL;
1918 atomtypes->surftens = NULL;
1919 atomtypes->atomnumber = NULL;
1920 atomtypes->gb_radius = NULL;
1921 atomtypes->S_hct = NULL;
1925 static void do_symtab(t_fileio *fio, t_symtab *symtab,gmx_bool bRead)
1931 gmx_fio_do_int(fio,symtab->nr);
1934 snew(symtab->symbuf,1);
1935 symbuf = symtab->symbuf;
1936 symbuf->bufsize = nr;
1937 snew(symbuf->buf,nr);
1938 for (i=0; (i<nr); i++) {
1939 gmx_fio_do_string(fio,buf);
1940 symbuf->buf[i]=strdup(buf);
1944 symbuf = symtab->symbuf;
1945 while (symbuf!=NULL) {
1946 for (i=0; (i<symbuf->bufsize) && (i<nr); i++)
1947 gmx_fio_do_string(fio,symbuf->buf[i]);
1949 symbuf=symbuf->next;
1952 gmx_fatal(FARGS,"nr of symtab strings left: %d",nr);
1956 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
1958 int i,j,ngrid,gs,nelem;
1960 gmx_fio_do_int(fio,cmap_grid->ngrid);
1961 gmx_fio_do_int(fio,cmap_grid->grid_spacing);
1963 ngrid = cmap_grid->ngrid;
1964 gs = cmap_grid->grid_spacing;
1969 snew(cmap_grid->cmapdata,ngrid);
1971 for(i=0;i<cmap_grid->ngrid;i++)
1973 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
1977 for(i=0;i<cmap_grid->ngrid;i++)
1979 for(j=0;j<nelem;j++)
1981 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4]);
1982 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+1]);
1983 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+2]);
1984 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+3]);
1990 void tpx_make_chain_identifiers(t_atoms *atoms,t_block *mols)
1996 /* We always assign a new chain number, but save the chain id characters
1997 * for larger molecules.
1999 #define CHAIN_MIN_ATOMS 15
2003 for(m=0; m<mols->nr; m++)
2006 a1=mols->index[m+1];
2007 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2016 for(a=a0; a<a1; a++)
2018 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2019 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2024 /* Blank out the chain id if there was only one chain */
2027 for(r=0; r<atoms->nres; r++)
2029 atoms->resinfo[r].chainid = ' ';
2034 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt,gmx_bool bRead,
2035 t_symtab *symtab, int file_version,
2036 gmx_groups_t *groups)
2040 if (file_version >= 57) {
2041 do_symstr(fio, &(molt->name),bRead,symtab);
2044 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2046 if (bRead && gmx_debug_at) {
2047 pr_atoms(debug,0,"atoms",&molt->atoms,TRUE);
2050 if (file_version >= 57) {
2051 do_ilists(fio, molt->ilist,bRead,file_version);
2053 do_block(fio, &molt->cgs,bRead,file_version);
2054 if (bRead && gmx_debug_at) {
2055 pr_block(debug,0,"cgs",&molt->cgs,TRUE);
2059 /* This used to be in the atoms struct */
2060 do_blocka(fio, &molt->excls, bRead, file_version);
2063 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb,gmx_bool bRead,
2068 gmx_fio_do_int(fio,molb->type);
2069 gmx_fio_do_int(fio,molb->nmol);
2070 gmx_fio_do_int(fio,molb->natoms_mol);
2071 /* Position restraint coordinates */
2072 gmx_fio_do_int(fio,molb->nposres_xA);
2073 if (molb->nposres_xA > 0) {
2075 snew(molb->posres_xA,molb->nposres_xA);
2077 gmx_fio_ndo_rvec(fio,molb->posres_xA,molb->nposres_xA);
2079 gmx_fio_do_int(fio,molb->nposres_xB);
2080 if (molb->nposres_xB > 0) {
2082 snew(molb->posres_xB,molb->nposres_xB);
2084 gmx_fio_ndo_rvec(fio,molb->posres_xB,molb->nposres_xB);
2089 static t_block mtop_mols(gmx_mtop_t *mtop)
2095 for(mb=0; mb<mtop->nmolblock; mb++) {
2096 mols.nr += mtop->molblock[mb].nmol;
2098 mols.nalloc_index = mols.nr + 1;
2099 snew(mols.index,mols.nalloc_index);
2104 for(mb=0; mb<mtop->nmolblock; mb++) {
2105 for(mol=0; mol<mtop->molblock[mb].nmol; mol++) {
2106 a += mtop->molblock[mb].natoms_mol;
2115 static void add_posres_molblock(gmx_mtop_t *mtop)
2120 gmx_molblock_t *molb;
2123 /* posres reference positions are stored in ip->posres (if present) and
2124 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2125 posres.pos0A are identical to fbposres.pos0. */
2126 il = &mtop->moltype[0].ilist[F_POSRES];
2127 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2128 if (il->nr == 0 && ilfb->nr == 0) {
2133 for(i=0; i<il->nr; i+=2) {
2134 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2135 am = max(am,il->iatoms[i+1]);
2136 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2137 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2138 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ]) {
2142 /* This loop is required if we have only flat-bottomed posres:
2144 - bFE == FALSE (no B-state for flat-bottomed posres) */
2147 for(i=0; i<ilfb->nr; i+=2) {
2148 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2149 am = max(am,ilfb->iatoms[i+1]);
2152 /* Make the posres coordinate block end at a molecule end */
2154 while(am >= mtop->mols.index[mol+1]) {
2157 molb = &mtop->molblock[0];
2158 molb->nposres_xA = mtop->mols.index[mol+1];
2159 snew(molb->posres_xA,molb->nposres_xA);
2161 molb->nposres_xB = molb->nposres_xA;
2162 snew(molb->posres_xB,molb->nposres_xB);
2164 molb->nposres_xB = 0;
2166 for(i=0; i<il->nr; i+=2) {
2167 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2168 a = il->iatoms[i+1];
2169 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2170 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2171 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2173 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2174 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2175 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2180 /* If only flat-bottomed posres are present, take reference pos from them.
2181 Here: bFE == FALSE */
2182 for(i=0; i<ilfb->nr; i+=2)
2184 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2185 a = ilfb->iatoms[i+1];
2186 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2187 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2188 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2193 static void set_disres_npair(gmx_mtop_t *mtop)
2200 ip = mtop->ffparams.iparams;
2202 for(mt=0; mt<mtop->nmoltype; mt++) {
2203 il = &mtop->moltype[mt].ilist[F_DISRES];
2207 for(i=0; i<il->nr; i+=3) {
2209 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label) {
2210 ip[a[i]].disres.npair = npair;
2218 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop,gmx_bool bRead,
2226 do_symtab(fio, &(mtop->symtab),bRead);
2228 pr_symtab(debug,0,"symtab",&mtop->symtab);
2230 do_symstr(fio, &(mtop->name),bRead,&(mtop->symtab));
2232 if (file_version >= 57) {
2233 do_ffparams(fio, &mtop->ffparams,bRead,file_version);
2235 gmx_fio_do_int(fio,mtop->nmoltype);
2240 snew(mtop->moltype,mtop->nmoltype);
2241 if (file_version < 57) {
2242 mtop->moltype[0].name = mtop->name;
2245 for(mt=0; mt<mtop->nmoltype; mt++) {
2246 do_moltype(fio, &mtop->moltype[mt],bRead,&mtop->symtab,file_version,
2250 if (file_version >= 57) {
2251 gmx_fio_do_int(fio,mtop->nmolblock);
2253 mtop->nmolblock = 1;
2256 snew(mtop->molblock,mtop->nmolblock);
2258 if (file_version >= 57) {
2259 for(mb=0; mb<mtop->nmolblock; mb++) {
2260 do_molblock(fio, &mtop->molblock[mb],bRead,file_version);
2262 gmx_fio_do_int(fio,mtop->natoms);
2264 mtop->molblock[0].type = 0;
2265 mtop->molblock[0].nmol = 1;
2266 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
2267 mtop->molblock[0].nposres_xA = 0;
2268 mtop->molblock[0].nposres_xB = 0;
2271 do_atomtypes (fio, &(mtop->atomtypes),bRead,&(mtop->symtab), file_version);
2273 pr_atomtypes(debug,0,"atomtypes",&mtop->atomtypes,TRUE);
2275 if (file_version < 57) {
2276 /* Debug statements are inside do_idef */
2277 do_idef (fio, &mtop->ffparams,&mtop->moltype[0],bRead,file_version);
2278 mtop->natoms = mtop->moltype[0].atoms.nr;
2281 if(file_version >= 65)
2283 do_cmap(fio, &mtop->ffparams.cmap_grid,bRead);
2287 mtop->ffparams.cmap_grid.ngrid = 0;
2288 mtop->ffparams.cmap_grid.grid_spacing = 0;
2289 mtop->ffparams.cmap_grid.cmapdata = NULL;
2292 if (file_version >= 57) {
2293 do_groups(fio, &mtop->groups,bRead,&(mtop->symtab),file_version);
2296 if (file_version < 57) {
2297 do_block(fio, &mtop->moltype[0].cgs,bRead,file_version);
2298 if (bRead && gmx_debug_at) {
2299 pr_block(debug,0,"cgs",&mtop->moltype[0].cgs,TRUE);
2301 do_block(fio, &mtop->mols,bRead,file_version);
2302 /* Add the posres coordinates to the molblock */
2303 add_posres_molblock(mtop);
2306 if (file_version >= 57) {
2307 mtop->mols = mtop_mols(mtop);
2310 pr_block(debug,0,"mols",&mtop->mols,TRUE);
2314 if (file_version < 51) {
2315 /* Here used to be the shake blocks */
2316 do_blocka(fio, &dumb,bRead,file_version);
2324 close_symtab(&(mtop->symtab));
2328 /* If TopOnlyOK is TRUE then we can read even future versions
2329 * of tpx files, provided the file_generation hasn't changed.
2330 * If it is FALSE, we need the inputrecord too, and bail out
2331 * if the file is newer than the program.
2333 * The version and generation if the topology (see top of this file)
2334 * are returned in the two last arguments.
2336 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2338 static void do_tpxheader(t_fileio *fio,gmx_bool bRead,t_tpxheader *tpx,
2339 gmx_bool TopOnlyOK, int *file_version,
2340 int *file_generation)
2343 char file_tag[STRLEN];
2350 gmx_fio_checktype(fio);
2351 gmx_fio_setdebug(fio,bDebugMode());
2353 /* NEW! XDR tpb file */
2354 precision = sizeof(real);
2356 gmx_fio_do_string(fio,buf);
2357 if (strncmp(buf,"VERSION",7))
2358 gmx_fatal(FARGS,"Can not read file %s,\n"
2359 " this file is from a Gromacs version which is older than 2.0\n"
2360 " Make a new one with grompp or use a gro or pdb file, if possible",
2361 gmx_fio_getname(fio));
2362 gmx_fio_do_int(fio,precision);
2363 bDouble = (precision == sizeof(double));
2364 if ((precision != sizeof(float)) && !bDouble)
2365 gmx_fatal(FARGS,"Unknown precision in file %s: real is %d bytes "
2366 "instead of %d or %d",
2367 gmx_fio_getname(fio),precision,sizeof(float),sizeof(double));
2368 gmx_fio_setprecision(fio,bDouble);
2369 fprintf(stderr,"Reading file %s, %s (%s precision)\n",
2370 gmx_fio_getname(fio),buf,bDouble ? "double" : "single");
2373 gmx_fio_write_string(fio,GromacsVersion());
2374 bDouble = (precision == sizeof(double));
2375 gmx_fio_setprecision(fio,bDouble);
2376 gmx_fio_do_int(fio,precision);
2378 sprintf(file_tag,"%s",tpx_tag);
2379 fgen = tpx_generation;
2382 /* Check versions! */
2383 gmx_fio_do_int(fio,fver);
2387 gmx_fio_do_string(fio,file_tag);
2393 /* Versions before 77 don't have the tag, set it to release */
2394 sprintf(file_tag,"%s",TPX_TAG_RELEASE);
2397 if (strcmp(file_tag,tpx_tag) != 0)
2399 fprintf(stderr,"Note: file tpx tag '%s', software tpx tag '%s'\n",
2402 /* We only support reading tpx files with the same tag as the code
2403 * or tpx files with the release tag and with lower version number.
2405 if (!(strcmp(file_tag,TPX_TAG_RELEASE) == 0 && fver < tpx_version))
2407 gmx_fatal(FARGS,"tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2408 gmx_fio_getname(fio),fver,file_tag,
2409 tpx_version,tpx_tag);
2416 gmx_fio_do_int(fio,fgen);
2423 if (file_version != NULL)
2425 *file_version = fver;
2427 if (file_generation != NULL)
2429 *file_generation = fgen;
2433 if ((fver <= tpx_incompatible_version) ||
2434 ((fver > tpx_version) && !TopOnlyOK) ||
2435 (fgen > tpx_generation))
2436 gmx_fatal(FARGS,"reading tpx file (%s) version %d with version %d program",
2437 gmx_fio_getname(fio),fver,tpx_version);
2439 do_section(fio,eitemHEADER,bRead);
2440 gmx_fio_do_int(fio,tpx->natoms);
2442 gmx_fio_do_int(fio,tpx->ngtc);
2446 gmx_fio_do_int(fio,idum);
2447 gmx_fio_do_real(fio,rdum);
2449 /*a better decision will eventually (5.0 or later) need to be made
2450 on how to treat the alchemical state of the system, which can now
2451 vary through a simulation, and cannot be completely described
2452 though a single lambda variable, or even a single state
2453 index. Eventually, should probably be a vector. MRS*/
2456 gmx_fio_do_int(fio,tpx->fep_state);
2458 gmx_fio_do_real(fio,tpx->lambda);
2459 gmx_fio_do_int(fio,tpx->bIr);
2460 gmx_fio_do_int(fio,tpx->bTop);
2461 gmx_fio_do_int(fio,tpx->bX);
2462 gmx_fio_do_int(fio,tpx->bV);
2463 gmx_fio_do_int(fio,tpx->bF);
2464 gmx_fio_do_int(fio,tpx->bBox);
2466 if((fgen > tpx_generation)) {
2467 /* This can only happen if TopOnlyOK=TRUE */
2472 static int do_tpx(t_fileio *fio, gmx_bool bRead,
2473 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop,
2474 gmx_bool bXVallocated)
2479 gmx_bool TopOnlyOK,bDum=TRUE;
2480 int file_version,file_generation;
2484 gmx_bool bPeriodicMols;
2487 tpx.natoms = state->natoms;
2488 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
2489 tpx.fep_state = state->fep_state;
2490 tpx.lambda = state->lambda[efptFEP];
2491 tpx.bIr = (ir != NULL);
2492 tpx.bTop = (mtop != NULL);
2493 tpx.bX = (state->x != NULL);
2494 tpx.bV = (state->v != NULL);
2495 tpx.bF = (f != NULL);
2499 TopOnlyOK = (ir==NULL);
2501 do_tpxheader(fio,bRead,&tpx,TopOnlyOK,&file_version,&file_generation);
2505 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
2506 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
2510 init_state(state,0,tpx.ngtc,0,0,0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
2511 state->natoms = tpx.natoms;
2512 state->nalloc = tpx.natoms;
2516 init_state(state,tpx.natoms,tpx.ngtc,0,0,0); /* nose-hoover chains */
2520 #define do_test(fio,b,p) if (bRead && (p!=NULL) && !b) gmx_fatal(FARGS,"No %s in %s",#p,gmx_fio_getname(fio))
2522 do_test(fio,tpx.bBox,state->box);
2523 do_section(fio,eitemBOX,bRead);
2525 gmx_fio_ndo_rvec(fio,state->box,DIM);
2526 if (file_version >= 51) {
2527 gmx_fio_ndo_rvec(fio,state->box_rel,DIM);
2529 /* We initialize box_rel after reading the inputrec */
2530 clear_mat(state->box_rel);
2532 if (file_version >= 28) {
2533 gmx_fio_ndo_rvec(fio,state->boxv,DIM);
2534 if (file_version < 56) {
2536 gmx_fio_ndo_rvec(fio,mdum,DIM);
2541 if (state->ngtc > 0 && file_version >= 28) {
2543 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
2544 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
2545 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
2546 snew(dumv,state->ngtc);
2547 if (file_version < 69) {
2548 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2550 /* These used to be the Berendsen tcoupl_lambda's */
2551 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2555 /* Prior to tpx version 26, the inputrec was here.
2556 * I moved it to enable partial forward-compatibility
2557 * for analysis/viewer programs.
2559 if(file_version<26) {
2560 do_test(fio,tpx.bIr,ir);
2561 do_section(fio,eitemIR,bRead);
2564 do_inputrec(fio, ir,bRead,file_version,
2565 mtop ? &mtop->ffparams.fudgeQQ : NULL);
2567 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2570 do_inputrec(fio, &dum_ir,bRead,file_version,
2571 mtop ? &mtop->ffparams.fudgeQQ :NULL);
2573 pr_inputrec(debug,0,"inputrec",&dum_ir,FALSE);
2574 done_inputrec(&dum_ir);
2580 do_test(fio,tpx.bTop,mtop);
2581 do_section(fio,eitemTOP,bRead);
2584 do_mtop(fio,mtop,bRead, file_version);
2586 do_mtop(fio,&dum_top,bRead,file_version);
2587 done_mtop(&dum_top,TRUE);
2590 do_test(fio,tpx.bX,state->x);
2591 do_section(fio,eitemX,bRead);
2594 state->flags |= (1<<estX);
2596 gmx_fio_ndo_rvec(fio,state->x,state->natoms);
2599 do_test(fio,tpx.bV,state->v);
2600 do_section(fio,eitemV,bRead);
2603 state->flags |= (1<<estV);
2605 gmx_fio_ndo_rvec(fio,state->v,state->natoms);
2608 do_test(fio,tpx.bF,f);
2609 do_section(fio,eitemF,bRead);
2610 if (tpx.bF) gmx_fio_ndo_rvec(fio,f,state->natoms);
2612 /* Starting with tpx version 26, we have the inputrec
2613 * at the end of the file, so we can ignore it
2614 * if the file is never than the software (but still the
2615 * same generation - see comments at the top of this file.
2620 bPeriodicMols = FALSE;
2621 if (file_version >= 26) {
2622 do_test(fio,tpx.bIr,ir);
2623 do_section(fio,eitemIR,bRead);
2625 if (file_version >= 53) {
2626 /* Removed the pbc info from do_inputrec, since we always want it */
2629 bPeriodicMols = ir->bPeriodicMols;
2631 gmx_fio_do_int(fio,ePBC);
2632 gmx_fio_do_gmx_bool(fio,bPeriodicMols);
2634 if (file_generation <= tpx_generation && ir) {
2635 do_inputrec(fio, ir,bRead,file_version,mtop ? &mtop->ffparams.fudgeQQ : NULL);
2637 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2638 if (file_version < 51)
2639 set_box_rel(ir,state);
2640 if (file_version < 53) {
2642 bPeriodicMols = ir->bPeriodicMols;
2645 if (bRead && ir && file_version >= 53) {
2646 /* We need to do this after do_inputrec, since that initializes ir */
2648 ir->bPeriodicMols = bPeriodicMols;
2657 if (state->ngtc == 0)
2659 /* Reading old version without tcoupl state data: set it */
2660 init_gtc_state(state,ir->opts.ngtc,0,ir->opts.nhchainlength);
2662 if (tpx.bTop && mtop)
2664 if (file_version < 57)
2666 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
2668 ir->eDisre = edrSimple;
2672 ir->eDisre = edrNone;
2675 set_disres_npair(mtop);
2679 if (tpx.bTop && mtop)
2681 gmx_mtop_finalize(mtop);
2684 if (file_version >= 57)
2688 env = getenv("GMX_NOCHARGEGROUPS");
2691 sscanf(env,"%d",&ienv);
2692 fprintf(stderr,"\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
2697 "Will make single atomic charge groups in non-solvent%s\n",
2698 ienv > 1 ? " and solvent" : "");
2699 gmx_mtop_make_atomic_charge_groups(mtop,ienv==1);
2701 fprintf(stderr,"\n");
2709 /************************************************************
2711 * The following routines are the exported ones
2713 ************************************************************/
2715 t_fileio *open_tpx(const char *fn,const char *mode)
2717 return gmx_fio_open(fn,mode);
2720 void close_tpx(t_fileio *fio)
2725 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
2726 int *file_version, int *file_generation)
2730 fio = open_tpx(fn,"r");
2731 do_tpxheader(fio,TRUE,tpx,TopOnlyOK,file_version,file_generation);
2735 void write_tpx_state(const char *fn,
2736 t_inputrec *ir,t_state *state,gmx_mtop_t *mtop)
2740 fio = open_tpx(fn,"w");
2741 do_tpx(fio,FALSE,ir,state,NULL,mtop,FALSE);
2745 void read_tpx_state(const char *fn,
2746 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop)
2750 fio = open_tpx(fn,"r");
2751 do_tpx(fio,TRUE,ir,state,f,mtop,FALSE);
2755 int read_tpx(const char *fn,
2756 t_inputrec *ir, matrix box,int *natoms,
2757 rvec *x,rvec *v,rvec *f,gmx_mtop_t *mtop)
2765 fio = open_tpx(fn,"r");
2766 ePBC = do_tpx(fio,TRUE,ir,&state,f,mtop,TRUE);
2768 *natoms = state.natoms;
2770 copy_mat(state.box,box);
2778 int read_tpx_top(const char *fn,
2779 t_inputrec *ir, matrix box,int *natoms,
2780 rvec *x,rvec *v,rvec *f,t_topology *top)
2786 ePBC = read_tpx(fn,ir,box,natoms,x,v,f,&mtop);
2788 *top = gmx_mtop_t_to_t_topology(&mtop);
2793 gmx_bool fn2bTPX(const char *file)
2795 switch (fn2ftp(file)) {
2805 gmx_bool read_tps_conf(const char *infile,char *title,t_topology *top,int *ePBC,
2806 rvec **x,rvec **v,matrix box,gmx_bool bMass)
2809 int natoms,i,version,generation;
2810 gmx_bool bTop,bXNULL=FALSE;
2812 t_topology *topconv;
2815 bTop = fn2bTPX(infile);
2818 read_tpxheader(infile,&header,TRUE,&version,&generation);
2820 snew(*x,header.natoms);
2822 snew(*v,header.natoms);
2824 *ePBC = read_tpx(infile,NULL,box,&natoms,
2825 (x==NULL) ? NULL : *x,(v==NULL) ? NULL : *v,NULL,mtop);
2826 *top = gmx_mtop_t_to_t_topology(mtop);
2828 strcpy(title,*top->name);
2829 tpx_make_chain_identifiers(&top->atoms,&top->mols);
2832 get_stx_coordnum(infile,&natoms);
2833 init_t_atoms(&top->atoms,natoms,(fn2ftp(infile) == efPDB));
2842 read_stx_conf(infile,title,&top->atoms,*x,(v==NULL) ? NULL : *v,ePBC,box);
2849 aps = gmx_atomprop_init();
2850 for(i=0; (i<natoms); i++)
2851 if (!gmx_atomprop_query(aps,epropMass,
2852 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2853 *top->atoms.atomname[i],
2854 &(top->atoms.atom[i].m))) {
2856 fprintf(debug,"Can not find mass for atom %s %d %s, setting to 1\n",
2857 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2858 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
2859 *top->atoms.atomname[i]);
2861 gmx_atomprop_destroy(aps);
2863 top->idef.ntypes=-1;