2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * Copyright (c) 2013, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
41 /* This file is completely threadsafe - keep it that way! */
47 #include "gmx_fatal.h"
61 #include "mtop_util.h"
63 #define TPX_TAG_RELEASE "release"
65 /* This is the tag string which is stored in the tpx file.
66 * Change this if you want to change the tpx format in a feature branch.
67 * This ensures that there will not be different tpx formats around which
68 * can not be distinguished.
70 static const char *tpx_tag = TPX_TAG_RELEASE;
72 /* This number should be increased whenever the file format changes! */
73 static const int tpx_version = 92;
75 /* This number should only be increased when you edit the TOPOLOGY section
76 * or the HEADER of the tpx format.
77 * This way we can maintain forward compatibility too for all analysis tools
78 * and/or external programs that only need to know the atom/residue names,
79 * charges, and bond connectivity.
81 * It first appeared in tpx version 26, when I also moved the inputrecord
82 * to the end of the tpx file, so we can just skip it if we only
85 static const int tpx_generation = 25;
87 /* This number should be the most recent backwards incompatible version
88 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
90 static const int tpx_incompatible_version = 9;
94 /* Struct used to maintain tpx compatibility when function types are added */
96 int fvnr; /* file version number in which the function type first appeared */
97 int ftype; /* function type */
101 * The entries should be ordered in:
102 * 1. ascending file version number
103 * 2. ascending function type number
105 /*static const t_ftupd ftupd[] = {
106 { 20, F_CUBICBONDS },
110 { 22, F_DISRESVIOL },
116 { 26, F_DIHRESVIOL },
117 { 30, F_CROSS_BOND_BONDS },
118 { 30, F_CROSS_BOND_ANGLES },
119 { 30, F_UREY_BRADLEY },
120 { 30, F_POLARIZATION },
124 * The entries should be ordered in:
125 * 1. ascending function type number
126 * 2. ascending file version number
128 /* question; what is the purpose of the commented code above? */
129 static const t_ftupd ftupd[] = {
130 { 20, F_CUBICBONDS },
135 { 43, F_TABBONDSNC },
136 { 70, F_RESTRBONDS },
137 { 76, F_LINEAR_ANGLES },
138 { 30, F_CROSS_BOND_BONDS },
139 { 30, F_CROSS_BOND_ANGLES },
140 { 30, F_UREY_BRADLEY },
141 { 34, F_QUARTIC_ANGLES },
151 { 72, F_NPSOLVATION },
153 { 41, F_LJC_PAIRS_NB },
156 { 32, F_COUL_RECIP },
158 { 30, F_POLARIZATION },
161 { 22, F_DISRESVIOL },
165 { 26, F_DIHRESVIOL },
170 { 46, F_ECONSERVED },
171 { 69, F_VTEMP_NOLONGERUSED},
173 { 54, F_DVDL_CONSTR },
174 { 76, F_ANHARM_POL },
177 { 79, F_DVDL_BONDED, },
178 { 79, F_DVDL_RESTRAINT },
179 { 79, F_DVDL_TEMPERATURE },
181 #define NFTUPD asize(ftupd)
183 /* Needed for backward compatibility */
186 static void _do_section(t_fileio *fio, int key, gmx_bool bRead, const char *src,
192 if (gmx_fio_getftp(fio) == efTPA)
196 gmx_fio_write_string(fio, itemstr[key]);
197 bDbg = gmx_fio_getdebug(fio);
198 gmx_fio_setdebug(fio, FALSE);
199 gmx_fio_write_string(fio, comment_str[key]);
200 gmx_fio_setdebug(fio, bDbg);
204 if (gmx_fio_getdebug(fio))
206 fprintf(stderr, "Looking for section %s (%s, %d)",
207 itemstr[key], src, line);
212 gmx_fio_do_string(fio, buf);
214 while ((gmx_strcasecmp(buf, itemstr[key]) != 0));
216 if (gmx_strcasecmp(buf, itemstr[key]) != 0)
218 gmx_fatal(FARGS, "\nCould not find section heading %s", itemstr[key]);
220 else if (gmx_fio_getdebug(fio))
222 fprintf(stderr, " and found it\n");
228 #define do_section(fio, key, bRead) _do_section(fio, key, bRead, __FILE__, __LINE__)
230 /**************************************************************
232 * Now the higer level routines that do io of the structures and arrays
234 **************************************************************/
235 static void do_pullgrp(t_fileio *fio, t_pullgrp *pgrp, gmx_bool bRead,
240 gmx_fio_do_int(fio, pgrp->nat);
243 snew(pgrp->ind, pgrp->nat);
245 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
246 gmx_fio_do_int(fio, pgrp->nweight);
249 snew(pgrp->weight, pgrp->nweight);
251 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
252 gmx_fio_do_int(fio, pgrp->pbcatom);
253 gmx_fio_do_rvec(fio, pgrp->vec);
254 gmx_fio_do_rvec(fio, pgrp->init);
255 gmx_fio_do_real(fio, pgrp->rate);
256 gmx_fio_do_real(fio, pgrp->k);
257 if (file_version >= 56)
259 gmx_fio_do_real(fio, pgrp->kB);
267 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
269 /* i is used in the ndo_double macro*/
273 int n_lambda = fepvals->n_lambda;
275 /* reset the lambda calculation window */
276 fepvals->lambda_start_n = 0;
277 fepvals->lambda_stop_n = n_lambda;
278 if (file_version >= 79)
284 snew(expand->init_lambda_weights, n_lambda);
286 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
287 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
290 gmx_fio_do_int(fio, expand->nstexpanded);
291 gmx_fio_do_int(fio, expand->elmcmove);
292 gmx_fio_do_int(fio, expand->elamstats);
293 gmx_fio_do_int(fio, expand->lmc_repeats);
294 gmx_fio_do_int(fio, expand->gibbsdeltalam);
295 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
296 gmx_fio_do_int(fio, expand->lmc_seed);
297 gmx_fio_do_real(fio, expand->mc_temp);
298 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
299 gmx_fio_do_int(fio, expand->nstTij);
300 gmx_fio_do_int(fio, expand->minvarmin);
301 gmx_fio_do_int(fio, expand->c_range);
302 gmx_fio_do_real(fio, expand->wl_scale);
303 gmx_fio_do_real(fio, expand->wl_ratio);
304 gmx_fio_do_real(fio, expand->init_wl_delta);
305 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
306 gmx_fio_do_int(fio, expand->elmceq);
307 gmx_fio_do_int(fio, expand->equil_steps);
308 gmx_fio_do_int(fio, expand->equil_samples);
309 gmx_fio_do_int(fio, expand->equil_n_at_lam);
310 gmx_fio_do_real(fio, expand->equil_wl_delta);
311 gmx_fio_do_real(fio, expand->equil_ratio);
315 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
318 if (file_version >= 79)
320 gmx_fio_do_int(fio, simtemp->eSimTempScale);
321 gmx_fio_do_real(fio, simtemp->simtemp_high);
322 gmx_fio_do_real(fio, simtemp->simtemp_low);
327 snew(simtemp->temperatures, n_lambda);
329 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
334 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
336 /* i is defined in the ndo_double macro; use g to iterate. */
341 /* free energy values */
343 if (file_version >= 79)
345 gmx_fio_do_int(fio, fepvals->init_fep_state);
346 gmx_fio_do_double(fio, fepvals->init_lambda);
347 gmx_fio_do_double(fio, fepvals->delta_lambda);
349 else if (file_version >= 59)
351 gmx_fio_do_double(fio, fepvals->init_lambda);
352 gmx_fio_do_double(fio, fepvals->delta_lambda);
356 gmx_fio_do_real(fio, rdum);
357 fepvals->init_lambda = rdum;
358 gmx_fio_do_real(fio, rdum);
359 fepvals->delta_lambda = rdum;
361 if (file_version >= 79)
363 gmx_fio_do_int(fio, fepvals->n_lambda);
366 snew(fepvals->all_lambda, efptNR);
368 for (g = 0; g < efptNR; g++)
370 if (fepvals->n_lambda > 0)
374 snew(fepvals->all_lambda[g], fepvals->n_lambda);
376 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
377 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
379 else if (fepvals->init_lambda >= 0)
381 fepvals->separate_dvdl[efptFEP] = TRUE;
385 else if (file_version >= 64)
387 gmx_fio_do_int(fio, fepvals->n_lambda);
392 snew(fepvals->all_lambda, efptNR);
393 /* still allocate the all_lambda array's contents. */
394 for (g = 0; g < efptNR; g++)
396 if (fepvals->n_lambda > 0)
398 snew(fepvals->all_lambda[g], fepvals->n_lambda);
402 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
404 if (fepvals->init_lambda >= 0)
408 fepvals->separate_dvdl[efptFEP] = TRUE;
412 /* copy the contents of the efptFEP lambda component to all
413 the other components */
414 for (g = 0; g < efptNR; g++)
416 for (h = 0; h < fepvals->n_lambda; h++)
420 fepvals->all_lambda[g][h] =
421 fepvals->all_lambda[efptFEP][h];
430 fepvals->n_lambda = 0;
431 fepvals->all_lambda = NULL;
432 if (fepvals->init_lambda >= 0)
434 fepvals->separate_dvdl[efptFEP] = TRUE;
437 if (file_version >= 13)
439 gmx_fio_do_real(fio, fepvals->sc_alpha);
443 fepvals->sc_alpha = 0;
445 if (file_version >= 38)
447 gmx_fio_do_int(fio, fepvals->sc_power);
451 fepvals->sc_power = 2;
453 if (file_version >= 79)
455 gmx_fio_do_real(fio, fepvals->sc_r_power);
459 fepvals->sc_r_power = 6.0;
461 if (file_version >= 15)
463 gmx_fio_do_real(fio, fepvals->sc_sigma);
467 fepvals->sc_sigma = 0.3;
471 if (file_version >= 71)
473 fepvals->sc_sigma_min = fepvals->sc_sigma;
477 fepvals->sc_sigma_min = 0;
480 if (file_version >= 79)
482 gmx_fio_do_int(fio, fepvals->bScCoul);
486 fepvals->bScCoul = TRUE;
488 if (file_version >= 64)
490 gmx_fio_do_int(fio, fepvals->nstdhdl);
494 fepvals->nstdhdl = 1;
497 if (file_version >= 73)
499 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
500 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
504 fepvals->separate_dhdl_file = esepdhdlfileYES;
505 fepvals->dhdl_derivatives = edhdlderivativesYES;
507 if (file_version >= 71)
509 gmx_fio_do_int(fio, fepvals->dh_hist_size);
510 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
514 fepvals->dh_hist_size = 0;
515 fepvals->dh_hist_spacing = 0.1;
517 if (file_version >= 79)
519 gmx_fio_do_int(fio, fepvals->bPrintEnergy);
523 fepvals->bPrintEnergy = FALSE;
526 /* handle lambda_neighbors */
527 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
529 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
530 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
531 (fepvals->init_lambda < 0) )
533 fepvals->lambda_start_n = (fepvals->init_fep_state -
534 fepvals->lambda_neighbors);
535 fepvals->lambda_stop_n = (fepvals->init_fep_state +
536 fepvals->lambda_neighbors + 1);
537 if (fepvals->lambda_start_n < 0)
539 fepvals->lambda_start_n = 0;;
541 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
543 fepvals->lambda_stop_n = fepvals->n_lambda;
548 fepvals->lambda_start_n = 0;
549 fepvals->lambda_stop_n = fepvals->n_lambda;
554 fepvals->lambda_start_n = 0;
555 fepvals->lambda_stop_n = fepvals->n_lambda;
559 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
563 gmx_fio_do_int(fio, pull->ngrp);
564 gmx_fio_do_int(fio, pull->eGeom);
565 gmx_fio_do_ivec(fio, pull->dim);
566 gmx_fio_do_real(fio, pull->cyl_r1);
567 gmx_fio_do_real(fio, pull->cyl_r0);
568 gmx_fio_do_real(fio, pull->constr_tol);
569 gmx_fio_do_int(fio, pull->nstxout);
570 gmx_fio_do_int(fio, pull->nstfout);
573 snew(pull->grp, pull->ngrp+1);
575 for (g = 0; g < pull->ngrp+1; g++)
577 do_pullgrp(fio, &pull->grp[g], bRead, file_version);
582 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
586 gmx_fio_do_int(fio, rotg->eType);
587 gmx_fio_do_int(fio, rotg->bMassW);
588 gmx_fio_do_int(fio, rotg->nat);
591 snew(rotg->ind, rotg->nat);
593 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
596 snew(rotg->x_ref, rotg->nat);
598 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
599 gmx_fio_do_rvec(fio, rotg->vec);
600 gmx_fio_do_rvec(fio, rotg->pivot);
601 gmx_fio_do_real(fio, rotg->rate);
602 gmx_fio_do_real(fio, rotg->k);
603 gmx_fio_do_real(fio, rotg->slab_dist);
604 gmx_fio_do_real(fio, rotg->min_gaussian);
605 gmx_fio_do_real(fio, rotg->eps);
606 gmx_fio_do_int(fio, rotg->eFittype);
607 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
608 gmx_fio_do_real(fio, rotg->PotAngle_step);
611 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
615 gmx_fio_do_int(fio, rot->ngrp);
616 gmx_fio_do_int(fio, rot->nstrout);
617 gmx_fio_do_int(fio, rot->nstsout);
620 snew(rot->grp, rot->ngrp);
622 for (g = 0; g < rot->ngrp; g++)
624 do_rotgrp(fio, &rot->grp[g], bRead);
629 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
630 int file_version, real *fudgeQQ)
632 int i, j, k, *tmp, idum = 0;
636 real zerotemptime, finish_t, init_temp, finish_temp;
638 if (file_version != tpx_version)
640 /* Give a warning about features that are not accessible */
641 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
642 file_version, tpx_version);
650 if (file_version == 0)
655 /* Basic inputrec stuff */
656 gmx_fio_do_int(fio, ir->eI);
657 if (file_version >= 62)
659 gmx_fio_do_gmx_large_int(fio, ir->nsteps);
663 gmx_fio_do_int(fio, idum);
666 if (file_version > 25)
668 if (file_version >= 62)
670 gmx_fio_do_gmx_large_int(fio, ir->init_step);
674 gmx_fio_do_int(fio, idum);
675 ir->init_step = idum;
683 if (file_version >= 58)
685 gmx_fio_do_int(fio, ir->simulation_part);
689 ir->simulation_part = 1;
692 if (file_version >= 67)
694 gmx_fio_do_int(fio, ir->nstcalcenergy);
698 ir->nstcalcenergy = 1;
700 if (file_version < 53)
702 /* The pbc info has been moved out of do_inputrec,
703 * since we always want it, also without reading the inputrec.
705 gmx_fio_do_int(fio, ir->ePBC);
706 if ((file_version <= 15) && (ir->ePBC == 2))
710 if (file_version >= 45)
712 gmx_fio_do_int(fio, ir->bPeriodicMols);
719 ir->bPeriodicMols = TRUE;
723 ir->bPeriodicMols = FALSE;
727 if (file_version >= 81)
729 gmx_fio_do_int(fio, ir->cutoff_scheme);
733 ir->cutoff_scheme = ecutsGROUP;
735 gmx_fio_do_int(fio, ir->ns_type);
736 gmx_fio_do_int(fio, ir->nstlist);
737 gmx_fio_do_int(fio, ir->ndelta);
738 if (file_version < 41)
740 gmx_fio_do_int(fio, idum);
741 gmx_fio_do_int(fio, idum);
743 if (file_version >= 45)
745 gmx_fio_do_real(fio, ir->rtpi);
751 gmx_fio_do_int(fio, ir->nstcomm);
752 if (file_version > 34)
754 gmx_fio_do_int(fio, ir->comm_mode);
756 else if (ir->nstcomm < 0)
758 ir->comm_mode = ecmANGULAR;
762 ir->comm_mode = ecmLINEAR;
764 ir->nstcomm = abs(ir->nstcomm);
766 if (file_version > 25)
768 gmx_fio_do_int(fio, ir->nstcheckpoint);
772 ir->nstcheckpoint = 0;
775 gmx_fio_do_int(fio, ir->nstcgsteep);
777 if (file_version >= 30)
779 gmx_fio_do_int(fio, ir->nbfgscorr);
786 gmx_fio_do_int(fio, ir->nstlog);
787 gmx_fio_do_int(fio, ir->nstxout);
788 gmx_fio_do_int(fio, ir->nstvout);
789 gmx_fio_do_int(fio, ir->nstfout);
790 gmx_fio_do_int(fio, ir->nstenergy);
791 gmx_fio_do_int(fio, ir->nstxtcout);
792 if (file_version >= 59)
794 gmx_fio_do_double(fio, ir->init_t);
795 gmx_fio_do_double(fio, ir->delta_t);
799 gmx_fio_do_real(fio, rdum);
801 gmx_fio_do_real(fio, rdum);
804 gmx_fio_do_real(fio, ir->xtcprec);
805 if (file_version < 19)
807 gmx_fio_do_int(fio, idum);
808 gmx_fio_do_int(fio, idum);
810 if (file_version < 18)
812 gmx_fio_do_int(fio, idum);
814 if (file_version >= 81)
816 gmx_fio_do_real(fio, ir->verletbuf_drift);
820 ir->verletbuf_drift = 0;
822 gmx_fio_do_real(fio, ir->rlist);
823 if (file_version >= 67)
825 gmx_fio_do_real(fio, ir->rlistlong);
827 if (file_version >= 82 && file_version != 90)
829 gmx_fio_do_int(fio, ir->nstcalclr);
833 /* Calculate at NS steps */
834 ir->nstcalclr = ir->nstlist;
836 gmx_fio_do_int(fio, ir->coulombtype);
837 if (file_version < 32 && ir->coulombtype == eelRF)
839 ir->coulombtype = eelRF_NEC;
841 if (file_version >= 81)
843 gmx_fio_do_int(fio, ir->coulomb_modifier);
847 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
849 gmx_fio_do_real(fio, ir->rcoulomb_switch);
850 gmx_fio_do_real(fio, ir->rcoulomb);
851 gmx_fio_do_int(fio, ir->vdwtype);
852 if (file_version >= 81)
854 gmx_fio_do_int(fio, ir->vdw_modifier);
858 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
860 gmx_fio_do_real(fio, ir->rvdw_switch);
861 gmx_fio_do_real(fio, ir->rvdw);
862 if (file_version < 67)
864 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
866 gmx_fio_do_int(fio, ir->eDispCorr);
867 gmx_fio_do_real(fio, ir->epsilon_r);
868 if (file_version >= 37)
870 gmx_fio_do_real(fio, ir->epsilon_rf);
874 if (EEL_RF(ir->coulombtype))
876 ir->epsilon_rf = ir->epsilon_r;
881 ir->epsilon_rf = 1.0;
884 if (file_version >= 29)
886 gmx_fio_do_real(fio, ir->tabext);
893 if (file_version > 25)
895 gmx_fio_do_int(fio, ir->gb_algorithm);
896 gmx_fio_do_int(fio, ir->nstgbradii);
897 gmx_fio_do_real(fio, ir->rgbradii);
898 gmx_fio_do_real(fio, ir->gb_saltconc);
899 gmx_fio_do_int(fio, ir->implicit_solvent);
903 ir->gb_algorithm = egbSTILL;
907 ir->implicit_solvent = eisNO;
909 if (file_version >= 55)
911 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
912 gmx_fio_do_real(fio, ir->gb_obc_alpha);
913 gmx_fio_do_real(fio, ir->gb_obc_beta);
914 gmx_fio_do_real(fio, ir->gb_obc_gamma);
915 if (file_version >= 60)
917 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
918 gmx_fio_do_int(fio, ir->sa_algorithm);
922 ir->gb_dielectric_offset = 0.009;
923 ir->sa_algorithm = esaAPPROX;
925 gmx_fio_do_real(fio, ir->sa_surface_tension);
927 /* Override sa_surface_tension if it is not changed in the mpd-file */
928 if (ir->sa_surface_tension < 0)
930 if (ir->gb_algorithm == egbSTILL)
932 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
934 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
936 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
943 /* Better use sensible values than insane (0.0) ones... */
944 ir->gb_epsilon_solvent = 80;
945 ir->gb_obc_alpha = 1.0;
946 ir->gb_obc_beta = 0.8;
947 ir->gb_obc_gamma = 4.85;
948 ir->sa_surface_tension = 2.092;
952 if (file_version >= 81)
954 gmx_fio_do_real(fio, ir->fourier_spacing);
958 ir->fourier_spacing = 0.0;
960 gmx_fio_do_int(fio, ir->nkx);
961 gmx_fio_do_int(fio, ir->nky);
962 gmx_fio_do_int(fio, ir->nkz);
963 gmx_fio_do_int(fio, ir->pme_order);
964 gmx_fio_do_real(fio, ir->ewald_rtol);
966 if (file_version >= 24)
968 gmx_fio_do_int(fio, ir->ewald_geometry);
972 ir->ewald_geometry = eewg3D;
975 if (file_version <= 17)
977 ir->epsilon_surface = 0;
978 if (file_version == 17)
980 gmx_fio_do_int(fio, idum);
985 gmx_fio_do_real(fio, ir->epsilon_surface);
988 gmx_fio_do_gmx_bool(fio, ir->bOptFFT);
990 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
991 gmx_fio_do_int(fio, ir->etc);
992 /* before version 18, ir->etc was a gmx_bool (ir->btc),
993 * but the values 0 and 1 still mean no and
994 * berendsen temperature coupling, respectively.
996 if (file_version >= 79)
998 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1000 if (file_version >= 71)
1002 gmx_fio_do_int(fio, ir->nsttcouple);
1006 ir->nsttcouple = ir->nstcalcenergy;
1008 if (file_version <= 15)
1010 gmx_fio_do_int(fio, idum);
1012 if (file_version <= 17)
1014 gmx_fio_do_int(fio, ir->epct);
1015 if (file_version <= 15)
1019 ir->epct = epctSURFACETENSION;
1021 gmx_fio_do_int(fio, idum);
1024 /* we have removed the NO alternative at the beginning */
1028 ir->epct = epctISOTROPIC;
1032 ir->epc = epcBERENDSEN;
1037 gmx_fio_do_int(fio, ir->epc);
1038 gmx_fio_do_int(fio, ir->epct);
1040 if (file_version >= 71)
1042 gmx_fio_do_int(fio, ir->nstpcouple);
1046 ir->nstpcouple = ir->nstcalcenergy;
1048 gmx_fio_do_real(fio, ir->tau_p);
1049 if (file_version <= 15)
1051 gmx_fio_do_rvec(fio, vdum);
1052 clear_mat(ir->ref_p);
1053 for (i = 0; i < DIM; i++)
1055 ir->ref_p[i][i] = vdum[i];
1060 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1061 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1062 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1064 if (file_version <= 15)
1066 gmx_fio_do_rvec(fio, vdum);
1067 clear_mat(ir->compress);
1068 for (i = 0; i < DIM; i++)
1070 ir->compress[i][i] = vdum[i];
1075 gmx_fio_do_rvec(fio, ir->compress[XX]);
1076 gmx_fio_do_rvec(fio, ir->compress[YY]);
1077 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1079 if (file_version >= 47)
1081 gmx_fio_do_int(fio, ir->refcoord_scaling);
1082 gmx_fio_do_rvec(fio, ir->posres_com);
1083 gmx_fio_do_rvec(fio, ir->posres_comB);
1087 ir->refcoord_scaling = erscNO;
1088 clear_rvec(ir->posres_com);
1089 clear_rvec(ir->posres_comB);
1091 if ((file_version > 25) && (file_version < 79))
1093 gmx_fio_do_int(fio, ir->andersen_seed);
1097 ir->andersen_seed = 0;
1099 if (file_version < 26)
1101 gmx_fio_do_gmx_bool(fio, bSimAnn);
1102 gmx_fio_do_real(fio, zerotemptime);
1105 if (file_version < 37)
1107 gmx_fio_do_real(fio, rdum);
1110 gmx_fio_do_real(fio, ir->shake_tol);
1111 if (file_version < 54)
1113 gmx_fio_do_real(fio, *fudgeQQ);
1116 gmx_fio_do_int(fio, ir->efep);
1117 if (file_version <= 14 && ir->efep != efepNO)
1121 do_fepvals(fio, ir->fepvals, bRead, file_version);
1123 if (file_version >= 79)
1125 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1128 ir->bSimTemp = TRUE;
1133 ir->bSimTemp = FALSE;
1137 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1140 if (file_version >= 79)
1142 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1145 ir->bExpanded = TRUE;
1149 ir->bExpanded = FALSE;
1154 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1156 if (file_version >= 57)
1158 gmx_fio_do_int(fio, ir->eDisre);
1160 gmx_fio_do_int(fio, ir->eDisreWeighting);
1161 if (file_version < 22)
1163 if (ir->eDisreWeighting == 0)
1165 ir->eDisreWeighting = edrwEqual;
1169 ir->eDisreWeighting = edrwConservative;
1172 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1173 gmx_fio_do_real(fio, ir->dr_fc);
1174 gmx_fio_do_real(fio, ir->dr_tau);
1175 gmx_fio_do_int(fio, ir->nstdisreout);
1176 if (file_version >= 22)
1178 gmx_fio_do_real(fio, ir->orires_fc);
1179 gmx_fio_do_real(fio, ir->orires_tau);
1180 gmx_fio_do_int(fio, ir->nstorireout);
1186 ir->nstorireout = 0;
1188 if (file_version >= 26 && file_version < 79)
1190 gmx_fio_do_real(fio, ir->dihre_fc);
1191 if (file_version < 56)
1193 gmx_fio_do_real(fio, rdum);
1194 gmx_fio_do_int(fio, idum);
1202 gmx_fio_do_real(fio, ir->em_stepsize);
1203 gmx_fio_do_real(fio, ir->em_tol);
1204 if (file_version >= 22)
1206 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1210 ir->bShakeSOR = TRUE;
1212 if (file_version >= 11)
1214 gmx_fio_do_int(fio, ir->niter);
1219 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1222 if (file_version >= 21)
1224 gmx_fio_do_real(fio, ir->fc_stepsize);
1228 ir->fc_stepsize = 0;
1230 gmx_fio_do_int(fio, ir->eConstrAlg);
1231 gmx_fio_do_int(fio, ir->nProjOrder);
1232 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1233 if (file_version <= 14)
1235 gmx_fio_do_int(fio, idum);
1237 if (file_version >= 26)
1239 gmx_fio_do_int(fio, ir->nLincsIter);
1244 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1247 if (file_version < 33)
1249 gmx_fio_do_real(fio, bd_temp);
1251 gmx_fio_do_real(fio, ir->bd_fric);
1252 gmx_fio_do_int(fio, ir->ld_seed);
1253 if (file_version >= 33)
1255 for (i = 0; i < DIM; i++)
1257 gmx_fio_do_rvec(fio, ir->deform[i]);
1262 for (i = 0; i < DIM; i++)
1264 clear_rvec(ir->deform[i]);
1267 if (file_version >= 14)
1269 gmx_fio_do_real(fio, ir->cos_accel);
1275 gmx_fio_do_int(fio, ir->userint1);
1276 gmx_fio_do_int(fio, ir->userint2);
1277 gmx_fio_do_int(fio, ir->userint3);
1278 gmx_fio_do_int(fio, ir->userint4);
1279 gmx_fio_do_real(fio, ir->userreal1);
1280 gmx_fio_do_real(fio, ir->userreal2);
1281 gmx_fio_do_real(fio, ir->userreal3);
1282 gmx_fio_do_real(fio, ir->userreal4);
1285 if (file_version >= 77)
1287 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1292 snew(ir->adress, 1);
1294 gmx_fio_do_int(fio, ir->adress->type);
1295 gmx_fio_do_real(fio, ir->adress->const_wf);
1296 gmx_fio_do_real(fio, ir->adress->ex_width);
1297 gmx_fio_do_real(fio, ir->adress->hy_width);
1298 gmx_fio_do_int(fio, ir->adress->icor);
1299 gmx_fio_do_int(fio, ir->adress->site);
1300 gmx_fio_do_rvec(fio, ir->adress->refs);
1301 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1302 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1303 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1304 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1308 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1310 if (ir->adress->n_tf_grps > 0)
1312 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1316 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1318 if (ir->adress->n_energy_grps > 0)
1320 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1326 ir->bAdress = FALSE;
1330 if (file_version >= 48)
1332 gmx_fio_do_int(fio, ir->ePull);
1333 if (ir->ePull != epullNO)
1339 do_pull(fio, ir->pull, bRead, file_version);
1344 ir->ePull = epullNO;
1347 /* Enforced rotation */
1348 if (file_version >= 74)
1350 gmx_fio_do_int(fio, ir->bRot);
1351 if (ir->bRot == TRUE)
1357 do_rot(fio, ir->rot, bRead);
1366 gmx_fio_do_int(fio, ir->opts.ngtc);
1367 if (file_version >= 69)
1369 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1373 ir->opts.nhchainlength = 1;
1375 gmx_fio_do_int(fio, ir->opts.ngacc);
1376 gmx_fio_do_int(fio, ir->opts.ngfrz);
1377 gmx_fio_do_int(fio, ir->opts.ngener);
1381 snew(ir->opts.nrdf, ir->opts.ngtc);
1382 snew(ir->opts.ref_t, ir->opts.ngtc);
1383 snew(ir->opts.annealing, ir->opts.ngtc);
1384 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1385 snew(ir->opts.anneal_time, ir->opts.ngtc);
1386 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1387 snew(ir->opts.tau_t, ir->opts.ngtc);
1388 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1389 snew(ir->opts.acc, ir->opts.ngacc);
1390 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1392 if (ir->opts.ngtc > 0)
1394 if (bRead && file_version < 13)
1396 snew(tmp, ir->opts.ngtc);
1397 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1398 for (i = 0; i < ir->opts.ngtc; i++)
1400 ir->opts.nrdf[i] = tmp[i];
1406 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1408 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1409 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1410 if (file_version < 33 && ir->eI == eiBD)
1412 for (i = 0; i < ir->opts.ngtc; i++)
1414 ir->opts.tau_t[i] = bd_temp;
1418 if (ir->opts.ngfrz > 0)
1420 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1422 if (ir->opts.ngacc > 0)
1424 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1426 if (file_version >= 12)
1428 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1429 ir->opts.ngener*ir->opts.ngener);
1432 if (bRead && file_version < 26)
1434 for (i = 0; i < ir->opts.ngtc; i++)
1438 ir->opts.annealing[i] = eannSINGLE;
1439 ir->opts.anneal_npoints[i] = 2;
1440 snew(ir->opts.anneal_time[i], 2);
1441 snew(ir->opts.anneal_temp[i], 2);
1442 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1443 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1444 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1445 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1446 ir->opts.anneal_time[i][0] = ir->init_t;
1447 ir->opts.anneal_time[i][1] = finish_t;
1448 ir->opts.anneal_temp[i][0] = init_temp;
1449 ir->opts.anneal_temp[i][1] = finish_temp;
1453 ir->opts.annealing[i] = eannNO;
1454 ir->opts.anneal_npoints[i] = 0;
1460 /* file version 26 or later */
1461 /* First read the lists with annealing and npoints for each group */
1462 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1463 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1464 for (j = 0; j < (ir->opts.ngtc); j++)
1466 k = ir->opts.anneal_npoints[j];
1469 snew(ir->opts.anneal_time[j], k);
1470 snew(ir->opts.anneal_temp[j], k);
1472 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1473 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1477 if (file_version >= 45)
1479 gmx_fio_do_int(fio, ir->nwall);
1480 gmx_fio_do_int(fio, ir->wall_type);
1481 if (file_version >= 50)
1483 gmx_fio_do_real(fio, ir->wall_r_linpot);
1487 ir->wall_r_linpot = -1;
1489 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1490 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1491 gmx_fio_do_real(fio, ir->wall_density[0]);
1492 gmx_fio_do_real(fio, ir->wall_density[1]);
1493 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1499 ir->wall_atomtype[0] = -1;
1500 ir->wall_atomtype[1] = -1;
1501 ir->wall_density[0] = 0;
1502 ir->wall_density[1] = 0;
1503 ir->wall_ewald_zfac = 3;
1505 /* Cosine stuff for electric fields */
1506 for (j = 0; (j < DIM); j++)
1508 gmx_fio_do_int(fio, ir->ex[j].n);
1509 gmx_fio_do_int(fio, ir->et[j].n);
1512 snew(ir->ex[j].a, ir->ex[j].n);
1513 snew(ir->ex[j].phi, ir->ex[j].n);
1514 snew(ir->et[j].a, ir->et[j].n);
1515 snew(ir->et[j].phi, ir->et[j].n);
1517 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1518 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1519 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1520 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1524 if (file_version >= 39)
1526 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1527 gmx_fio_do_int(fio, ir->QMMMscheme);
1528 gmx_fio_do_real(fio, ir->scalefactor);
1529 gmx_fio_do_int(fio, ir->opts.ngQM);
1532 snew(ir->opts.QMmethod, ir->opts.ngQM);
1533 snew(ir->opts.QMbasis, ir->opts.ngQM);
1534 snew(ir->opts.QMcharge, ir->opts.ngQM);
1535 snew(ir->opts.QMmult, ir->opts.ngQM);
1536 snew(ir->opts.bSH, ir->opts.ngQM);
1537 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1538 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1539 snew(ir->opts.SAon, ir->opts.ngQM);
1540 snew(ir->opts.SAoff, ir->opts.ngQM);
1541 snew(ir->opts.SAsteps, ir->opts.ngQM);
1542 snew(ir->opts.bOPT, ir->opts.ngQM);
1543 snew(ir->opts.bTS, ir->opts.ngQM);
1545 if (ir->opts.ngQM > 0)
1547 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1548 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1549 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1550 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1551 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1552 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1553 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1554 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1555 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1556 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1557 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1558 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1560 /* end of QMMM stuff */
1565 static void do_harm(t_fileio *fio, t_iparams *iparams)
1567 gmx_fio_do_real(fio, iparams->harmonic.rA);
1568 gmx_fio_do_real(fio, iparams->harmonic.krA);
1569 gmx_fio_do_real(fio, iparams->harmonic.rB);
1570 gmx_fio_do_real(fio, iparams->harmonic.krB);
1573 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1574 gmx_bool bRead, int file_version)
1581 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1591 do_harm(fio, iparams);
1592 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1594 /* Correct incorrect storage of parameters */
1595 iparams->pdihs.phiB = iparams->pdihs.phiA;
1596 iparams->pdihs.cpB = iparams->pdihs.cpA;
1599 case F_LINEAR_ANGLES:
1600 gmx_fio_do_real(fio, iparams->linangle.klinA);
1601 gmx_fio_do_real(fio, iparams->linangle.aA);
1602 gmx_fio_do_real(fio, iparams->linangle.klinB);
1603 gmx_fio_do_real(fio, iparams->linangle.aB);
1606 gmx_fio_do_real(fio, iparams->fene.bm);
1607 gmx_fio_do_real(fio, iparams->fene.kb);
1610 gmx_fio_do_real(fio, iparams->restraint.lowA);
1611 gmx_fio_do_real(fio, iparams->restraint.up1A);
1612 gmx_fio_do_real(fio, iparams->restraint.up2A);
1613 gmx_fio_do_real(fio, iparams->restraint.kA);
1614 gmx_fio_do_real(fio, iparams->restraint.lowB);
1615 gmx_fio_do_real(fio, iparams->restraint.up1B);
1616 gmx_fio_do_real(fio, iparams->restraint.up2B);
1617 gmx_fio_do_real(fio, iparams->restraint.kB);
1623 gmx_fio_do_real(fio, iparams->tab.kA);
1624 gmx_fio_do_int(fio, iparams->tab.table);
1625 gmx_fio_do_real(fio, iparams->tab.kB);
1627 case F_CROSS_BOND_BONDS:
1628 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1629 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1630 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1632 case F_CROSS_BOND_ANGLES:
1633 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1634 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1635 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1636 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1638 case F_UREY_BRADLEY:
1639 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1640 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1641 gmx_fio_do_real(fio, iparams->u_b.r13A);
1642 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1643 if (file_version >= 79)
1645 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1646 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1647 gmx_fio_do_real(fio, iparams->u_b.r13B);
1648 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1652 iparams->u_b.thetaB = iparams->u_b.thetaA;
1653 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1654 iparams->u_b.r13B = iparams->u_b.r13A;
1655 iparams->u_b.kUBB = iparams->u_b.kUBA;
1658 case F_QUARTIC_ANGLES:
1659 gmx_fio_do_real(fio, iparams->qangle.theta);
1660 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1663 gmx_fio_do_real(fio, iparams->bham.a);
1664 gmx_fio_do_real(fio, iparams->bham.b);
1665 gmx_fio_do_real(fio, iparams->bham.c);
1668 gmx_fio_do_real(fio, iparams->morse.b0A);
1669 gmx_fio_do_real(fio, iparams->morse.cbA);
1670 gmx_fio_do_real(fio, iparams->morse.betaA);
1671 if (file_version >= 79)
1673 gmx_fio_do_real(fio, iparams->morse.b0B);
1674 gmx_fio_do_real(fio, iparams->morse.cbB);
1675 gmx_fio_do_real(fio, iparams->morse.betaB);
1679 iparams->morse.b0B = iparams->morse.b0A;
1680 iparams->morse.cbB = iparams->morse.cbA;
1681 iparams->morse.betaB = iparams->morse.betaA;
1685 gmx_fio_do_real(fio, iparams->cubic.b0);
1686 gmx_fio_do_real(fio, iparams->cubic.kb);
1687 gmx_fio_do_real(fio, iparams->cubic.kcub);
1691 case F_POLARIZATION:
1692 gmx_fio_do_real(fio, iparams->polarize.alpha);
1695 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1696 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1697 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1700 if (file_version < 31)
1702 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1704 gmx_fio_do_real(fio, iparams->wpol.al_x);
1705 gmx_fio_do_real(fio, iparams->wpol.al_y);
1706 gmx_fio_do_real(fio, iparams->wpol.al_z);
1707 gmx_fio_do_real(fio, iparams->wpol.rOH);
1708 gmx_fio_do_real(fio, iparams->wpol.rHH);
1709 gmx_fio_do_real(fio, iparams->wpol.rOD);
1712 gmx_fio_do_real(fio, iparams->thole.a);
1713 gmx_fio_do_real(fio, iparams->thole.alpha1);
1714 gmx_fio_do_real(fio, iparams->thole.alpha2);
1715 gmx_fio_do_real(fio, iparams->thole.rfac);
1718 gmx_fio_do_real(fio, iparams->lj.c6);
1719 gmx_fio_do_real(fio, iparams->lj.c12);
1722 gmx_fio_do_real(fio, iparams->lj14.c6A);
1723 gmx_fio_do_real(fio, iparams->lj14.c12A);
1724 gmx_fio_do_real(fio, iparams->lj14.c6B);
1725 gmx_fio_do_real(fio, iparams->lj14.c12B);
1728 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1729 gmx_fio_do_real(fio, iparams->ljc14.qi);
1730 gmx_fio_do_real(fio, iparams->ljc14.qj);
1731 gmx_fio_do_real(fio, iparams->ljc14.c6);
1732 gmx_fio_do_real(fio, iparams->ljc14.c12);
1734 case F_LJC_PAIRS_NB:
1735 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1736 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1737 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1738 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1744 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1745 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1746 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1748 /* Read the incorrectly stored multiplicity */
1749 gmx_fio_do_real(fio, iparams->harmonic.rB);
1750 gmx_fio_do_real(fio, iparams->harmonic.krB);
1751 iparams->pdihs.phiB = iparams->pdihs.phiA;
1752 iparams->pdihs.cpB = iparams->pdihs.cpA;
1756 gmx_fio_do_real(fio, iparams->pdihs.phiB);
1757 gmx_fio_do_real(fio, iparams->pdihs.cpB);
1758 gmx_fio_do_int(fio, iparams->pdihs.mult);
1762 gmx_fio_do_int(fio, iparams->disres.label);
1763 gmx_fio_do_int(fio, iparams->disres.type);
1764 gmx_fio_do_real(fio, iparams->disres.low);
1765 gmx_fio_do_real(fio, iparams->disres.up1);
1766 gmx_fio_do_real(fio, iparams->disres.up2);
1767 gmx_fio_do_real(fio, iparams->disres.kfac);
1770 gmx_fio_do_int(fio, iparams->orires.ex);
1771 gmx_fio_do_int(fio, iparams->orires.label);
1772 gmx_fio_do_int(fio, iparams->orires.power);
1773 gmx_fio_do_real(fio, iparams->orires.c);
1774 gmx_fio_do_real(fio, iparams->orires.obs);
1775 gmx_fio_do_real(fio, iparams->orires.kfac);
1778 if (file_version < 82)
1780 gmx_fio_do_int(fio, idum);
1781 gmx_fio_do_int(fio, idum);
1783 gmx_fio_do_real(fio, iparams->dihres.phiA);
1784 gmx_fio_do_real(fio, iparams->dihres.dphiA);
1785 gmx_fio_do_real(fio, iparams->dihres.kfacA);
1786 if (file_version >= 82)
1788 gmx_fio_do_real(fio, iparams->dihres.phiB);
1789 gmx_fio_do_real(fio, iparams->dihres.dphiB);
1790 gmx_fio_do_real(fio, iparams->dihres.kfacB);
1794 iparams->dihres.phiB = iparams->dihres.phiA;
1795 iparams->dihres.dphiB = iparams->dihres.dphiA;
1796 iparams->dihres.kfacB = iparams->dihres.kfacA;
1800 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
1801 gmx_fio_do_rvec(fio, iparams->posres.fcA);
1802 if (bRead && file_version < 27)
1804 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
1805 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
1809 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
1810 gmx_fio_do_rvec(fio, iparams->posres.fcB);
1814 gmx_fio_do_int(fio, iparams->fbposres.geom);
1815 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
1816 gmx_fio_do_real(fio, iparams->fbposres.r);
1817 gmx_fio_do_real(fio, iparams->fbposres.k);
1820 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1821 if (file_version >= 25)
1823 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1827 /* Fourier dihedrals are internally represented
1828 * as Ryckaert-Bellemans since those are faster to compute.
1830 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1831 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1835 gmx_fio_do_real(fio, iparams->constr.dA);
1836 gmx_fio_do_real(fio, iparams->constr.dB);
1839 gmx_fio_do_real(fio, iparams->settle.doh);
1840 gmx_fio_do_real(fio, iparams->settle.dhh);
1843 gmx_fio_do_real(fio, iparams->vsite.a);
1848 gmx_fio_do_real(fio, iparams->vsite.a);
1849 gmx_fio_do_real(fio, iparams->vsite.b);
1854 gmx_fio_do_real(fio, iparams->vsite.a);
1855 gmx_fio_do_real(fio, iparams->vsite.b);
1856 gmx_fio_do_real(fio, iparams->vsite.c);
1859 gmx_fio_do_int(fio, iparams->vsiten.n);
1860 gmx_fio_do_real(fio, iparams->vsiten.a);
1865 /* We got rid of some parameters in version 68 */
1866 if (bRead && file_version < 68)
1868 gmx_fio_do_real(fio, rdum);
1869 gmx_fio_do_real(fio, rdum);
1870 gmx_fio_do_real(fio, rdum);
1871 gmx_fio_do_real(fio, rdum);
1873 gmx_fio_do_real(fio, iparams->gb.sar);
1874 gmx_fio_do_real(fio, iparams->gb.st);
1875 gmx_fio_do_real(fio, iparams->gb.pi);
1876 gmx_fio_do_real(fio, iparams->gb.gbr);
1877 gmx_fio_do_real(fio, iparams->gb.bmlt);
1880 gmx_fio_do_int(fio, iparams->cmap.cmapA);
1881 gmx_fio_do_int(fio, iparams->cmap.cmapB);
1884 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1885 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1889 gmx_fio_unset_comment(fio);
1893 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
1900 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1902 if (file_version < 44)
1904 for (i = 0; i < MAXNODES; i++)
1906 gmx_fio_do_int(fio, idum);
1909 gmx_fio_do_int(fio, ilist->nr);
1912 snew(ilist->iatoms, ilist->nr);
1914 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
1917 gmx_fio_unset_comment(fio);
1921 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1922 gmx_bool bRead, int file_version)
1927 gmx_fio_do_int(fio, ffparams->atnr);
1928 if (file_version < 57)
1930 gmx_fio_do_int(fio, idum);
1932 gmx_fio_do_int(fio, ffparams->ntypes);
1935 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
1936 ffparams->atnr, ffparams->ntypes);
1940 snew(ffparams->functype, ffparams->ntypes);
1941 snew(ffparams->iparams, ffparams->ntypes);
1943 /* Read/write all the function types */
1944 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
1947 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
1950 if (file_version >= 66)
1952 gmx_fio_do_double(fio, ffparams->reppow);
1956 ffparams->reppow = 12.0;
1959 if (file_version >= 57)
1961 gmx_fio_do_real(fio, ffparams->fudgeQQ);
1964 /* Check whether all these function types are supported by the code.
1965 * In practice the code is backwards compatible, which means that the
1966 * numbering may have to be altered from old numbering to new numbering
1968 for (i = 0; (i < ffparams->ntypes); i++)
1972 /* Loop over file versions */
1973 for (k = 0; (k < NFTUPD); k++)
1975 /* Compare the read file_version to the update table */
1976 if ((file_version < ftupd[k].fvnr) &&
1977 (ffparams->functype[i] >= ftupd[k].ftype))
1979 ffparams->functype[i] += 1;
1982 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
1983 i, ffparams->functype[i],
1984 interaction_function[ftupd[k].ftype].longname);
1991 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
1995 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2000 static void add_settle_atoms(t_ilist *ilist)
2004 /* Settle used to only store the first atom: add the other two */
2005 srenew(ilist->iatoms, 2*ilist->nr);
2006 for (i = ilist->nr/2-1; i >= 0; i--)
2008 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2009 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2010 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2011 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2013 ilist->nr = 2*ilist->nr;
2016 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2019 int i, j, renum[F_NRE];
2023 for (j = 0; (j < F_NRE); j++)
2028 for (k = 0; k < NFTUPD; k++)
2030 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2039 ilist[j].iatoms = NULL;
2043 do_ilist(fio, &ilist[j], bRead, file_version, j);
2044 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2046 add_settle_atoms(&ilist[j]);
2050 if (bRead && gmx_debug_at)
2051 pr_ilist(debug,0,interaction_function[j].longname,
2052 functype,&ilist[j],TRUE);
2057 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2058 gmx_bool bRead, int file_version)
2060 do_ffparams(fio, ffparams, bRead, file_version);
2062 if (file_version >= 54)
2064 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2067 do_ilists(fio, molt->ilist, bRead, file_version);
2070 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2072 int i, idum, dum_nra, *dum_a;
2074 if (file_version < 44)
2076 for (i = 0; i < MAXNODES; i++)
2078 gmx_fio_do_int(fio, idum);
2081 gmx_fio_do_int(fio, block->nr);
2082 if (file_version < 51)
2084 gmx_fio_do_int(fio, dum_nra);
2088 if ((block->nalloc_index > 0) && (NULL != block->index))
2090 sfree(block->index);
2092 block->nalloc_index = block->nr+1;
2093 snew(block->index, block->nalloc_index);
2095 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2097 if (file_version < 51 && dum_nra > 0)
2099 snew(dum_a, dum_nra);
2100 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2105 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2110 if (file_version < 44)
2112 for (i = 0; i < MAXNODES; i++)
2114 gmx_fio_do_int(fio, idum);
2117 gmx_fio_do_int(fio, block->nr);
2118 gmx_fio_do_int(fio, block->nra);
2121 block->nalloc_index = block->nr+1;
2122 snew(block->index, block->nalloc_index);
2123 block->nalloc_a = block->nra;
2124 snew(block->a, block->nalloc_a);
2126 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2127 gmx_fio_ndo_int(fio, block->a, block->nra);
2130 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2131 int file_version, gmx_groups_t *groups, int atnr)
2135 gmx_fio_do_real(fio, atom->m);
2136 gmx_fio_do_real(fio, atom->q);
2137 gmx_fio_do_real(fio, atom->mB);
2138 gmx_fio_do_real(fio, atom->qB);
2139 gmx_fio_do_ushort(fio, atom->type);
2140 gmx_fio_do_ushort(fio, atom->typeB);
2141 gmx_fio_do_int(fio, atom->ptype);
2142 gmx_fio_do_int(fio, atom->resind);
2143 if (file_version >= 52)
2145 gmx_fio_do_int(fio, atom->atomnumber);
2149 atom->atomnumber = NOTSET;
2151 if (file_version < 23)
2155 else if (file_version < 39)
2164 if (file_version < 57)
2166 unsigned char uchar[egcNR];
2167 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2168 for (i = myngrp; (i < ngrp); i++)
2172 /* Copy the old data format to the groups struct */
2173 for (i = 0; i < ngrp; i++)
2175 groups->grpnr[i][atnr] = uchar[i];
2180 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2185 if (file_version < 23)
2189 else if (file_version < 39)
2198 for (j = 0; (j < ngrp); j++)
2202 gmx_fio_do_int(fio, grps[j].nr);
2205 snew(grps[j].nm_ind, grps[j].nr);
2207 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2212 snew(grps[j].nm_ind, grps[j].nr);
2217 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2223 gmx_fio_do_int(fio, ls);
2224 *nm = get_symtab_handle(symtab, ls);
2228 ls = lookup_symtab(symtab, *nm);
2229 gmx_fio_do_int(fio, ls);
2233 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2238 for (j = 0; (j < nstr); j++)
2240 do_symstr(fio, &(nm[j]), bRead, symtab);
2244 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2245 t_symtab *symtab, int file_version)
2249 for (j = 0; (j < n); j++)
2251 do_symstr(fio, &(ri[j].name), bRead, symtab);
2252 if (file_version >= 63)
2254 gmx_fio_do_int(fio, ri[j].nr);
2255 gmx_fio_do_uchar(fio, ri[j].ic);
2265 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2267 gmx_groups_t *groups)
2271 gmx_fio_do_int(fio, atoms->nr);
2272 gmx_fio_do_int(fio, atoms->nres);
2273 if (file_version < 57)
2275 gmx_fio_do_int(fio, groups->ngrpname);
2276 for (i = 0; i < egcNR; i++)
2278 groups->ngrpnr[i] = atoms->nr;
2279 snew(groups->grpnr[i], groups->ngrpnr[i]);
2284 snew(atoms->atom, atoms->nr);
2285 snew(atoms->atomname, atoms->nr);
2286 snew(atoms->atomtype, atoms->nr);
2287 snew(atoms->atomtypeB, atoms->nr);
2288 snew(atoms->resinfo, atoms->nres);
2289 if (file_version < 57)
2291 snew(groups->grpname, groups->ngrpname);
2293 atoms->pdbinfo = NULL;
2295 for (i = 0; (i < atoms->nr); i++)
2297 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2299 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2300 if (bRead && (file_version <= 20))
2302 for (i = 0; i < atoms->nr; i++)
2304 atoms->atomtype[i] = put_symtab(symtab, "?");
2305 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2310 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2311 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2313 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2315 if (file_version < 57)
2317 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2319 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2323 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2324 gmx_bool bRead, t_symtab *symtab,
2329 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2330 gmx_fio_do_int(fio, groups->ngrpname);
2333 snew(groups->grpname, groups->ngrpname);
2335 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2336 for (g = 0; g < egcNR; g++)
2338 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2339 if (groups->ngrpnr[g] == 0)
2343 groups->grpnr[g] = NULL;
2350 snew(groups->grpnr[g], groups->ngrpnr[g]);
2352 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2357 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2362 if (file_version > 25)
2364 gmx_fio_do_int(fio, atomtypes->nr);
2368 snew(atomtypes->radius, j);
2369 snew(atomtypes->vol, j);
2370 snew(atomtypes->surftens, j);
2371 snew(atomtypes->atomnumber, j);
2372 snew(atomtypes->gb_radius, j);
2373 snew(atomtypes->S_hct, j);
2375 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2376 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2377 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2378 if (file_version >= 40)
2380 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2382 if (file_version >= 60)
2384 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2385 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2390 /* File versions prior to 26 cannot do GBSA,
2391 * so they dont use this structure
2394 atomtypes->radius = NULL;
2395 atomtypes->vol = NULL;
2396 atomtypes->surftens = NULL;
2397 atomtypes->atomnumber = NULL;
2398 atomtypes->gb_radius = NULL;
2399 atomtypes->S_hct = NULL;
2403 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2409 gmx_fio_do_int(fio, symtab->nr);
2413 snew(symtab->symbuf, 1);
2414 symbuf = symtab->symbuf;
2415 symbuf->bufsize = nr;
2416 snew(symbuf->buf, nr);
2417 for (i = 0; (i < nr); i++)
2419 gmx_fio_do_string(fio, buf);
2420 symbuf->buf[i] = strdup(buf);
2425 symbuf = symtab->symbuf;
2426 while (symbuf != NULL)
2428 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2430 gmx_fio_do_string(fio, symbuf->buf[i]);
2433 symbuf = symbuf->next;
2437 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2442 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2444 int i, j, ngrid, gs, nelem;
2446 gmx_fio_do_int(fio, cmap_grid->ngrid);
2447 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2449 ngrid = cmap_grid->ngrid;
2450 gs = cmap_grid->grid_spacing;
2455 snew(cmap_grid->cmapdata, ngrid);
2457 for (i = 0; i < cmap_grid->ngrid; i++)
2459 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2463 for (i = 0; i < cmap_grid->ngrid; i++)
2465 for (j = 0; j < nelem; j++)
2467 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2468 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2469 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2470 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2476 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2478 int m, a, a0, a1, r;
2482 /* We always assign a new chain number, but save the chain id characters
2483 * for larger molecules.
2485 #define CHAIN_MIN_ATOMS 15
2489 for (m = 0; m < mols->nr; m++)
2491 a0 = mols->index[m];
2492 a1 = mols->index[m+1];
2493 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2502 for (a = a0; a < a1; a++)
2504 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2505 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2510 /* Blank out the chain id if there was only one chain */
2513 for (r = 0; r < atoms->nres; r++)
2515 atoms->resinfo[r].chainid = ' ';
2520 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2521 t_symtab *symtab, int file_version,
2522 gmx_groups_t *groups)
2526 if (file_version >= 57)
2528 do_symstr(fio, &(molt->name), bRead, symtab);
2531 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2533 if (bRead && gmx_debug_at)
2535 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2538 if (file_version >= 57)
2540 do_ilists(fio, molt->ilist, bRead, file_version);
2542 do_block(fio, &molt->cgs, bRead, file_version);
2543 if (bRead && gmx_debug_at)
2545 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2549 /* This used to be in the atoms struct */
2550 do_blocka(fio, &molt->excls, bRead, file_version);
2553 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2557 gmx_fio_do_int(fio, molb->type);
2558 gmx_fio_do_int(fio, molb->nmol);
2559 gmx_fio_do_int(fio, molb->natoms_mol);
2560 /* Position restraint coordinates */
2561 gmx_fio_do_int(fio, molb->nposres_xA);
2562 if (molb->nposres_xA > 0)
2566 snew(molb->posres_xA, molb->nposres_xA);
2568 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2570 gmx_fio_do_int(fio, molb->nposres_xB);
2571 if (molb->nposres_xB > 0)
2575 snew(molb->posres_xB, molb->nposres_xB);
2577 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2582 static t_block mtop_mols(gmx_mtop_t *mtop)
2588 for (mb = 0; mb < mtop->nmolblock; mb++)
2590 mols.nr += mtop->molblock[mb].nmol;
2592 mols.nalloc_index = mols.nr + 1;
2593 snew(mols.index, mols.nalloc_index);
2598 for (mb = 0; mb < mtop->nmolblock; mb++)
2600 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2602 a += mtop->molblock[mb].natoms_mol;
2611 static void add_posres_molblock(gmx_mtop_t *mtop)
2616 gmx_molblock_t *molb;
2619 /* posres reference positions are stored in ip->posres (if present) and
2620 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2621 posres.pos0A are identical to fbposres.pos0. */
2622 il = &mtop->moltype[0].ilist[F_POSRES];
2623 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2624 if (il->nr == 0 && ilfb->nr == 0)
2630 for (i = 0; i < il->nr; i += 2)
2632 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2633 am = max(am, il->iatoms[i+1]);
2634 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2635 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2636 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2641 /* This loop is required if we have only flat-bottomed posres:
2643 - bFE == FALSE (no B-state for flat-bottomed posres) */
2646 for (i = 0; i < ilfb->nr; i += 2)
2648 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2649 am = max(am, ilfb->iatoms[i+1]);
2652 /* Make the posres coordinate block end at a molecule end */
2654 while (am >= mtop->mols.index[mol+1])
2658 molb = &mtop->molblock[0];
2659 molb->nposres_xA = mtop->mols.index[mol+1];
2660 snew(molb->posres_xA, molb->nposres_xA);
2663 molb->nposres_xB = molb->nposres_xA;
2664 snew(molb->posres_xB, molb->nposres_xB);
2668 molb->nposres_xB = 0;
2670 for (i = 0; i < il->nr; i += 2)
2672 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2673 a = il->iatoms[i+1];
2674 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2675 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2676 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2679 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2680 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2681 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2686 /* If only flat-bottomed posres are present, take reference pos from them.
2687 Here: bFE == FALSE */
2688 for (i = 0; i < ilfb->nr; i += 2)
2690 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2691 a = ilfb->iatoms[i+1];
2692 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2693 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2694 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2699 static void set_disres_npair(gmx_mtop_t *mtop)
2706 ip = mtop->ffparams.iparams;
2708 for (mt = 0; mt < mtop->nmoltype; mt++)
2710 il = &mtop->moltype[mt].ilist[F_DISRES];
2715 for (i = 0; i < il->nr; i += 3)
2718 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2720 ip[a[i]].disres.npair = npair;
2728 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2738 do_symtab(fio, &(mtop->symtab), bRead);
2741 pr_symtab(debug, 0, "symtab", &mtop->symtab);
2744 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2746 if (file_version >= 57)
2748 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
2750 gmx_fio_do_int(fio, mtop->nmoltype);
2758 snew(mtop->moltype, mtop->nmoltype);
2759 if (file_version < 57)
2761 mtop->moltype[0].name = mtop->name;
2764 for (mt = 0; mt < mtop->nmoltype; mt++)
2766 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
2770 if (file_version >= 57)
2772 gmx_fio_do_int(fio, mtop->nmolblock);
2776 mtop->nmolblock = 1;
2780 snew(mtop->molblock, mtop->nmolblock);
2782 if (file_version >= 57)
2784 for (mb = 0; mb < mtop->nmolblock; mb++)
2786 do_molblock(fio, &mtop->molblock[mb], bRead);
2788 gmx_fio_do_int(fio, mtop->natoms);
2792 mtop->molblock[0].type = 0;
2793 mtop->molblock[0].nmol = 1;
2794 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
2795 mtop->molblock[0].nposres_xA = 0;
2796 mtop->molblock[0].nposres_xB = 0;
2799 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
2802 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
2805 if (file_version < 57)
2807 /* Debug statements are inside do_idef */
2808 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
2809 mtop->natoms = mtop->moltype[0].atoms.nr;
2812 if (file_version >= 65)
2814 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
2818 mtop->ffparams.cmap_grid.ngrid = 0;
2819 mtop->ffparams.cmap_grid.grid_spacing = 0;
2820 mtop->ffparams.cmap_grid.cmapdata = NULL;
2823 if (file_version >= 57)
2825 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
2828 if (file_version < 57)
2830 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
2831 if (bRead && gmx_debug_at)
2833 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
2835 do_block(fio, &mtop->mols, bRead, file_version);
2836 /* Add the posres coordinates to the molblock */
2837 add_posres_molblock(mtop);
2841 if (file_version >= 57)
2843 done_block(&mtop->mols);
2844 mtop->mols = mtop_mols(mtop);
2848 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
2852 if (file_version < 51)
2854 /* Here used to be the shake blocks */
2855 do_blocka(fio, &dumb, bRead, file_version);
2868 close_symtab(&(mtop->symtab));
2872 /* If TopOnlyOK is TRUE then we can read even future versions
2873 * of tpx files, provided the file_generation hasn't changed.
2874 * If it is FALSE, we need the inputrecord too, and bail out
2875 * if the file is newer than the program.
2877 * The version and generation if the topology (see top of this file)
2878 * are returned in the two last arguments.
2880 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2882 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
2883 gmx_bool TopOnlyOK, int *file_version,
2884 int *file_generation)
2887 char file_tag[STRLEN];
2894 gmx_fio_checktype(fio);
2895 gmx_fio_setdebug(fio, bDebugMode());
2897 /* NEW! XDR tpb file */
2898 precision = sizeof(real);
2901 gmx_fio_do_string(fio, buf);
2902 if (strncmp(buf, "VERSION", 7))
2904 gmx_fatal(FARGS, "Can not read file %s,\n"
2905 " this file is from a Gromacs version which is older than 2.0\n"
2906 " Make a new one with grompp or use a gro or pdb file, if possible",
2907 gmx_fio_getname(fio));
2909 gmx_fio_do_int(fio, precision);
2910 bDouble = (precision == sizeof(double));
2911 if ((precision != sizeof(float)) && !bDouble)
2913 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
2914 "instead of %d or %d",
2915 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
2917 gmx_fio_setprecision(fio, bDouble);
2918 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
2919 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
2923 gmx_fio_write_string(fio, GromacsVersion());
2924 bDouble = (precision == sizeof(double));
2925 gmx_fio_setprecision(fio, bDouble);
2926 gmx_fio_do_int(fio, precision);
2928 sprintf(file_tag, "%s", tpx_tag);
2929 fgen = tpx_generation;
2932 /* Check versions! */
2933 gmx_fio_do_int(fio, fver);
2935 /* This is for backward compatibility with development versions 77-79
2936 * where the tag was, mistakenly, placed before the generation,
2937 * which would cause a segv instead of a proper error message
2938 * when reading the topology only from tpx with <77 code.
2940 if (fver >= 77 && fver <= 79)
2942 gmx_fio_do_string(fio, file_tag);
2947 gmx_fio_do_int(fio, fgen);
2956 gmx_fio_do_string(fio, file_tag);
2962 /* Versions before 77 don't have the tag, set it to release */
2963 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
2966 if (strcmp(file_tag, tpx_tag) != 0)
2968 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
2971 /* We only support reading tpx files with the same tag as the code
2972 * or tpx files with the release tag and with lower version number.
2974 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
2976 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2977 gmx_fio_getname(fio), fver, file_tag,
2978 tpx_version, tpx_tag);
2983 if (file_version != NULL)
2985 *file_version = fver;
2987 if (file_generation != NULL)
2989 *file_generation = fgen;
2993 if ((fver <= tpx_incompatible_version) ||
2994 ((fver > tpx_version) && !TopOnlyOK) ||
2995 (fgen > tpx_generation) ||
2996 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
2998 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
2999 gmx_fio_getname(fio), fver, tpx_version);
3002 do_section(fio, eitemHEADER, bRead);
3003 gmx_fio_do_int(fio, tpx->natoms);
3006 gmx_fio_do_int(fio, tpx->ngtc);
3014 gmx_fio_do_int(fio, idum);
3015 gmx_fio_do_real(fio, rdum);
3017 /*a better decision will eventually (5.0 or later) need to be made
3018 on how to treat the alchemical state of the system, which can now
3019 vary through a simulation, and cannot be completely described
3020 though a single lambda variable, or even a single state
3021 index. Eventually, should probably be a vector. MRS*/
3024 gmx_fio_do_int(fio, tpx->fep_state);
3026 gmx_fio_do_real(fio, tpx->lambda);
3027 gmx_fio_do_int(fio, tpx->bIr);
3028 gmx_fio_do_int(fio, tpx->bTop);
3029 gmx_fio_do_int(fio, tpx->bX);
3030 gmx_fio_do_int(fio, tpx->bV);
3031 gmx_fio_do_int(fio, tpx->bF);
3032 gmx_fio_do_int(fio, tpx->bBox);
3034 if ((fgen > tpx_generation))
3036 /* This can only happen if TopOnlyOK=TRUE */
3041 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3042 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3043 gmx_bool bXVallocated)
3049 int file_version, file_generation;
3053 gmx_bool bPeriodicMols;
3057 tpx.natoms = state->natoms;
3058 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3059 tpx.fep_state = state->fep_state;
3060 tpx.lambda = state->lambda[efptFEP];
3061 tpx.bIr = (ir != NULL);
3062 tpx.bTop = (mtop != NULL);
3063 tpx.bX = (state->x != NULL);
3064 tpx.bV = (state->v != NULL);
3065 tpx.bF = (f != NULL);
3069 TopOnlyOK = (ir == NULL);
3071 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3076 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3077 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3082 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3083 state->natoms = tpx.natoms;
3084 state->nalloc = tpx.natoms;
3090 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3094 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3096 do_test(fio, tpx.bBox, state->box);
3097 do_section(fio, eitemBOX, bRead);
3100 gmx_fio_ndo_rvec(fio, state->box, DIM);
3101 if (file_version >= 51)
3103 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3107 /* We initialize box_rel after reading the inputrec */
3108 clear_mat(state->box_rel);
3110 if (file_version >= 28)
3112 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3113 if (file_version < 56)
3116 gmx_fio_ndo_rvec(fio, mdum, DIM);
3121 if (state->ngtc > 0 && file_version >= 28)
3124 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3125 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3126 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3127 snew(dumv, state->ngtc);
3128 if (file_version < 69)
3130 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3132 /* These used to be the Berendsen tcoupl_lambda's */
3133 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3137 /* Prior to tpx version 26, the inputrec was here.
3138 * I moved it to enable partial forward-compatibility
3139 * for analysis/viewer programs.
3141 if (file_version < 26)
3143 do_test(fio, tpx.bIr, ir);
3144 do_section(fio, eitemIR, bRead);
3149 do_inputrec(fio, ir, bRead, file_version,
3150 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3153 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3158 do_inputrec(fio, &dum_ir, bRead, file_version,
3159 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3162 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3164 done_inputrec(&dum_ir);
3170 do_test(fio, tpx.bTop, mtop);
3171 do_section(fio, eitemTOP, bRead);
3176 do_mtop(fio, mtop, bRead, file_version);
3180 do_mtop(fio, &dum_top, bRead, file_version);
3181 done_mtop(&dum_top, TRUE);
3184 do_test(fio, tpx.bX, state->x);
3185 do_section(fio, eitemX, bRead);
3190 state->flags |= (1<<estX);
3192 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3195 do_test(fio, tpx.bV, state->v);
3196 do_section(fio, eitemV, bRead);
3201 state->flags |= (1<<estV);
3203 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3206 do_test(fio, tpx.bF, f);
3207 do_section(fio, eitemF, bRead);
3210 gmx_fio_ndo_rvec(fio, f, state->natoms);
3213 /* Starting with tpx version 26, we have the inputrec
3214 * at the end of the file, so we can ignore it
3215 * if the file is never than the software (but still the
3216 * same generation - see comments at the top of this file.
3221 bPeriodicMols = FALSE;
3222 if (file_version >= 26)
3224 do_test(fio, tpx.bIr, ir);
3225 do_section(fio, eitemIR, bRead);
3228 if (file_version >= 53)
3230 /* Removed the pbc info from do_inputrec, since we always want it */
3234 bPeriodicMols = ir->bPeriodicMols;
3236 gmx_fio_do_int(fio, ePBC);
3237 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3239 if (file_generation <= tpx_generation && ir)
3241 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3244 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3246 if (file_version < 51)
3248 set_box_rel(ir, state);
3250 if (file_version < 53)
3253 bPeriodicMols = ir->bPeriodicMols;
3256 if (bRead && ir && file_version >= 53)
3258 /* We need to do this after do_inputrec, since that initializes ir */
3260 ir->bPeriodicMols = bPeriodicMols;
3269 if (state->ngtc == 0)
3271 /* Reading old version without tcoupl state data: set it */
3272 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3274 if (tpx.bTop && mtop)
3276 if (file_version < 57)
3278 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3280 ir->eDisre = edrSimple;
3284 ir->eDisre = edrNone;
3287 set_disres_npair(mtop);
3291 if (tpx.bTop && mtop)
3293 gmx_mtop_finalize(mtop);
3296 if (file_version >= 57)
3300 env = getenv("GMX_NOCHARGEGROUPS");
3303 sscanf(env, "%d", &ienv);
3304 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3309 "Will make single atomic charge groups in non-solvent%s\n",
3310 ienv > 1 ? " and solvent" : "");
3311 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3313 fprintf(stderr, "\n");
3321 /************************************************************
3323 * The following routines are the exported ones
3325 ************************************************************/
3327 t_fileio *open_tpx(const char *fn, const char *mode)
3329 return gmx_fio_open(fn, mode);
3332 void close_tpx(t_fileio *fio)
3337 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3338 int *file_version, int *file_generation)
3342 fio = open_tpx(fn, "r");
3343 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3347 void write_tpx_state(const char *fn,
3348 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3352 fio = open_tpx(fn, "w");
3353 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3357 void read_tpx_state(const char *fn,
3358 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3362 fio = open_tpx(fn, "r");
3363 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3367 int read_tpx(const char *fn,
3368 t_inputrec *ir, matrix box, int *natoms,
3369 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3377 fio = open_tpx(fn, "r");
3378 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3380 *natoms = state.natoms;
3383 copy_mat(state.box, box);
3392 int read_tpx_top(const char *fn,
3393 t_inputrec *ir, matrix box, int *natoms,
3394 rvec *x, rvec *v, rvec *f, t_topology *top)
3400 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3402 *top = gmx_mtop_t_to_t_topology(&mtop);
3407 gmx_bool fn2bTPX(const char *file)
3409 switch (fn2ftp(file))
3420 static void done_gmx_groups_t(gmx_groups_t *g)
3424 for (i = 0; (i < egcNR); i++)
3426 if (NULL != g->grps[i].nm_ind)
3428 sfree(g->grps[i].nm_ind);
3429 g->grps[i].nm_ind = NULL;
3431 if (NULL != g->grpnr[i])
3437 /* The contents of this array is in symtab, don't free it here */
3441 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3442 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3445 int natoms, i, version, generation;
3446 gmx_bool bTop, bXNULL = FALSE;
3448 t_topology *topconv;
3451 bTop = fn2bTPX(infile);
3455 read_tpxheader(infile, &header, TRUE, &version, &generation);
3458 snew(*x, header.natoms);
3462 snew(*v, header.natoms);
3465 *ePBC = read_tpx(infile, NULL, box, &natoms,
3466 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3467 *top = gmx_mtop_t_to_t_topology(mtop);
3468 /* In this case we need to throw away the group data too */
3469 done_gmx_groups_t(&mtop->groups);
3471 strcpy(title, *top->name);
3472 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3476 get_stx_coordnum(infile, &natoms);
3477 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3488 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3496 aps = gmx_atomprop_init();
3497 for (i = 0; (i < natoms); i++)
3499 if (!gmx_atomprop_query(aps, epropMass,
3500 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3501 *top->atoms.atomname[i],
3502 &(top->atoms.atom[i].m)))
3506 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3507 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3508 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3509 *top->atoms.atomname[i]);
3513 gmx_atomprop_destroy(aps);
3515 top->idef.ntypes = -1;