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! */
43 #include <thread_mpi.h>
51 #include "gmx_fatal.h"
65 #include "mtop_util.h"
67 #define TPX_TAG_RELEASE "release"
69 /* This is the tag string which is stored in the tpx file.
70 * Change this if you want to change the tpx format in a feature branch.
71 * This ensures that there will not be different tpx formats around which
72 * can not be distinguished.
74 static const char *tpx_tag = TPX_TAG_RELEASE;
76 /* This number should be increased whenever the file format changes! */
77 static const int tpx_version = 92;
79 /* This number should only be increased when you edit the TOPOLOGY section
80 * or the HEADER of the tpx format.
81 * This way we can maintain forward compatibility too for all analysis tools
82 * and/or external programs that only need to know the atom/residue names,
83 * charges, and bond connectivity.
85 * It first appeared in tpx version 26, when I also moved the inputrecord
86 * to the end of the tpx file, so we can just skip it if we only
89 static const int tpx_generation = 25;
91 /* This number should be the most recent backwards incompatible version
92 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
94 static const int tpx_incompatible_version = 9;
98 /* Struct used to maintain tpx compatibility when function types are added */
100 int fvnr; /* file version number in which the function type first appeared */
101 int ftype; /* function type */
105 * The entries should be ordered in:
106 * 1. ascending file version number
107 * 2. ascending function type number
109 /*static const t_ftupd ftupd[] = {
110 { 20, F_CUBICBONDS },
114 { 22, F_DISRESVIOL },
120 { 26, F_DIHRESVIOL },
121 { 30, F_CROSS_BOND_BONDS },
122 { 30, F_CROSS_BOND_ANGLES },
123 { 30, F_UREY_BRADLEY },
124 { 30, F_POLARIZATION },
128 * The entries should be ordered in:
129 * 1. ascending function type number
130 * 2. ascending file version number
132 /* question; what is the purpose of the commented code above? */
133 static const t_ftupd ftupd[] = {
134 { 20, F_CUBICBONDS },
139 { 43, F_TABBONDSNC },
140 { 70, F_RESTRBONDS },
141 { 76, F_LINEAR_ANGLES },
142 { 30, F_CROSS_BOND_BONDS },
143 { 30, F_CROSS_BOND_ANGLES },
144 { 30, F_UREY_BRADLEY },
145 { 34, F_QUARTIC_ANGLES },
155 { 72, F_NPSOLVATION },
157 { 41, F_LJC_PAIRS_NB },
160 { 32, F_COUL_RECIP },
162 { 30, F_POLARIZATION },
165 { 22, F_DISRESVIOL },
169 { 26, F_DIHRESVIOL },
174 { 46, F_ECONSERVED },
175 { 69, F_VTEMP_NOLONGERUSED},
177 { 54, F_DVDL_CONSTR },
178 { 76, F_ANHARM_POL },
181 { 79, F_DVDL_BONDED, },
182 { 79, F_DVDL_RESTRAINT },
183 { 79, F_DVDL_TEMPERATURE },
185 #define NFTUPD asize(ftupd)
187 /* Needed for backward compatibility */
190 static void _do_section(t_fileio *fio, int key, gmx_bool bRead, const char *src,
196 if (gmx_fio_getftp(fio) == efTPA)
200 gmx_fio_write_string(fio, itemstr[key]);
201 bDbg = gmx_fio_getdebug(fio);
202 gmx_fio_setdebug(fio, FALSE);
203 gmx_fio_write_string(fio, comment_str[key]);
204 gmx_fio_setdebug(fio, bDbg);
208 if (gmx_fio_getdebug(fio))
210 fprintf(stderr, "Looking for section %s (%s, %d)",
211 itemstr[key], src, line);
216 gmx_fio_do_string(fio, buf);
218 while ((gmx_strcasecmp(buf, itemstr[key]) != 0));
220 if (gmx_strcasecmp(buf, itemstr[key]) != 0)
222 gmx_fatal(FARGS, "\nCould not find section heading %s", itemstr[key]);
224 else if (gmx_fio_getdebug(fio))
226 fprintf(stderr, " and found it\n");
232 #define do_section(fio, key, bRead) _do_section(fio, key, bRead, __FILE__, __LINE__)
234 /**************************************************************
236 * Now the higer level routines that do io of the structures and arrays
238 **************************************************************/
239 static void do_pullgrp(t_fileio *fio, t_pullgrp *pgrp, gmx_bool bRead,
244 gmx_fio_do_int(fio, pgrp->nat);
247 snew(pgrp->ind, pgrp->nat);
249 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
250 gmx_fio_do_int(fio, pgrp->nweight);
253 snew(pgrp->weight, pgrp->nweight);
255 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
256 gmx_fio_do_int(fio, pgrp->pbcatom);
257 gmx_fio_do_rvec(fio, pgrp->vec);
258 gmx_fio_do_rvec(fio, pgrp->init);
259 gmx_fio_do_real(fio, pgrp->rate);
260 gmx_fio_do_real(fio, pgrp->k);
261 if (file_version >= 56)
263 gmx_fio_do_real(fio, pgrp->kB);
271 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
273 /* i is used in the ndo_double macro*/
277 int n_lambda = fepvals->n_lambda;
279 /* reset the lambda calculation window */
280 fepvals->lambda_start_n = 0;
281 fepvals->lambda_stop_n = n_lambda;
282 if (file_version >= 79)
288 snew(expand->init_lambda_weights, n_lambda);
290 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
291 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
294 gmx_fio_do_int(fio, expand->nstexpanded);
295 gmx_fio_do_int(fio, expand->elmcmove);
296 gmx_fio_do_int(fio, expand->elamstats);
297 gmx_fio_do_int(fio, expand->lmc_repeats);
298 gmx_fio_do_int(fio, expand->gibbsdeltalam);
299 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
300 gmx_fio_do_int(fio, expand->lmc_seed);
301 gmx_fio_do_real(fio, expand->mc_temp);
302 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
303 gmx_fio_do_int(fio, expand->nstTij);
304 gmx_fio_do_int(fio, expand->minvarmin);
305 gmx_fio_do_int(fio, expand->c_range);
306 gmx_fio_do_real(fio, expand->wl_scale);
307 gmx_fio_do_real(fio, expand->wl_ratio);
308 gmx_fio_do_real(fio, expand->init_wl_delta);
309 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
310 gmx_fio_do_int(fio, expand->elmceq);
311 gmx_fio_do_int(fio, expand->equil_steps);
312 gmx_fio_do_int(fio, expand->equil_samples);
313 gmx_fio_do_int(fio, expand->equil_n_at_lam);
314 gmx_fio_do_real(fio, expand->equil_wl_delta);
315 gmx_fio_do_real(fio, expand->equil_ratio);
319 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
322 if (file_version >= 79)
324 gmx_fio_do_int(fio, simtemp->eSimTempScale);
325 gmx_fio_do_real(fio, simtemp->simtemp_high);
326 gmx_fio_do_real(fio, simtemp->simtemp_low);
331 snew(simtemp->temperatures, n_lambda);
333 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
338 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
340 /* i is defined in the ndo_double macro; use g to iterate. */
345 /* free energy values */
347 if (file_version >= 79)
349 gmx_fio_do_int(fio, fepvals->init_fep_state);
350 gmx_fio_do_double(fio, fepvals->init_lambda);
351 gmx_fio_do_double(fio, fepvals->delta_lambda);
353 else if (file_version >= 59)
355 gmx_fio_do_double(fio, fepvals->init_lambda);
356 gmx_fio_do_double(fio, fepvals->delta_lambda);
360 gmx_fio_do_real(fio, rdum);
361 fepvals->init_lambda = rdum;
362 gmx_fio_do_real(fio, rdum);
363 fepvals->delta_lambda = rdum;
365 if (file_version >= 79)
367 gmx_fio_do_int(fio, fepvals->n_lambda);
370 snew(fepvals->all_lambda, efptNR);
372 for (g = 0; g < efptNR; g++)
374 if (fepvals->n_lambda > 0)
378 snew(fepvals->all_lambda[g], fepvals->n_lambda);
380 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
381 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
383 else if (fepvals->init_lambda >= 0)
385 fepvals->separate_dvdl[efptFEP] = TRUE;
389 else if (file_version >= 64)
391 gmx_fio_do_int(fio, fepvals->n_lambda);
396 snew(fepvals->all_lambda, efptNR);
397 /* still allocate the all_lambda array's contents. */
398 for (g = 0; g < efptNR; g++)
400 if (fepvals->n_lambda > 0)
402 snew(fepvals->all_lambda[g], fepvals->n_lambda);
406 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
408 if (fepvals->init_lambda >= 0)
412 fepvals->separate_dvdl[efptFEP] = TRUE;
416 /* copy the contents of the efptFEP lambda component to all
417 the other components */
418 for (g = 0; g < efptNR; g++)
420 for (h = 0; h < fepvals->n_lambda; h++)
424 fepvals->all_lambda[g][h] =
425 fepvals->all_lambda[efptFEP][h];
434 fepvals->n_lambda = 0;
435 fepvals->all_lambda = NULL;
436 if (fepvals->init_lambda >= 0)
438 fepvals->separate_dvdl[efptFEP] = TRUE;
441 if (file_version >= 13)
443 gmx_fio_do_real(fio, fepvals->sc_alpha);
447 fepvals->sc_alpha = 0;
449 if (file_version >= 38)
451 gmx_fio_do_int(fio, fepvals->sc_power);
455 fepvals->sc_power = 2;
457 if (file_version >= 79)
459 gmx_fio_do_real(fio, fepvals->sc_r_power);
463 fepvals->sc_r_power = 6.0;
465 if (file_version >= 15)
467 gmx_fio_do_real(fio, fepvals->sc_sigma);
471 fepvals->sc_sigma = 0.3;
475 if (file_version >= 71)
477 fepvals->sc_sigma_min = fepvals->sc_sigma;
481 fepvals->sc_sigma_min = 0;
484 if (file_version >= 79)
486 gmx_fio_do_int(fio, fepvals->bScCoul);
490 fepvals->bScCoul = TRUE;
492 if (file_version >= 64)
494 gmx_fio_do_int(fio, fepvals->nstdhdl);
498 fepvals->nstdhdl = 1;
501 if (file_version >= 73)
503 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
504 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
508 fepvals->separate_dhdl_file = esepdhdlfileYES;
509 fepvals->dhdl_derivatives = edhdlderivativesYES;
511 if (file_version >= 71)
513 gmx_fio_do_int(fio, fepvals->dh_hist_size);
514 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
518 fepvals->dh_hist_size = 0;
519 fepvals->dh_hist_spacing = 0.1;
521 if (file_version >= 79)
523 gmx_fio_do_int(fio, fepvals->bPrintEnergy);
527 fepvals->bPrintEnergy = FALSE;
530 /* handle lambda_neighbors */
531 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
533 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
534 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
535 (fepvals->init_lambda < 0) )
537 fepvals->lambda_start_n = (fepvals->init_fep_state -
538 fepvals->lambda_neighbors);
539 fepvals->lambda_stop_n = (fepvals->init_fep_state +
540 fepvals->lambda_neighbors + 1);
541 if (fepvals->lambda_start_n < 0)
543 fepvals->lambda_start_n = 0;;
545 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
547 fepvals->lambda_stop_n = fepvals->n_lambda;
552 fepvals->lambda_start_n = 0;
553 fepvals->lambda_stop_n = fepvals->n_lambda;
558 fepvals->lambda_start_n = 0;
559 fepvals->lambda_stop_n = fepvals->n_lambda;
563 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
567 gmx_fio_do_int(fio, pull->ngrp);
568 gmx_fio_do_int(fio, pull->eGeom);
569 gmx_fio_do_ivec(fio, pull->dim);
570 gmx_fio_do_real(fio, pull->cyl_r1);
571 gmx_fio_do_real(fio, pull->cyl_r0);
572 gmx_fio_do_real(fio, pull->constr_tol);
573 gmx_fio_do_int(fio, pull->nstxout);
574 gmx_fio_do_int(fio, pull->nstfout);
577 snew(pull->grp, pull->ngrp+1);
579 for (g = 0; g < pull->ngrp+1; g++)
581 do_pullgrp(fio, &pull->grp[g], bRead, file_version);
586 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
590 gmx_fio_do_int(fio, rotg->eType);
591 gmx_fio_do_int(fio, rotg->bMassW);
592 gmx_fio_do_int(fio, rotg->nat);
595 snew(rotg->ind, rotg->nat);
597 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
600 snew(rotg->x_ref, rotg->nat);
602 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
603 gmx_fio_do_rvec(fio, rotg->vec);
604 gmx_fio_do_rvec(fio, rotg->pivot);
605 gmx_fio_do_real(fio, rotg->rate);
606 gmx_fio_do_real(fio, rotg->k);
607 gmx_fio_do_real(fio, rotg->slab_dist);
608 gmx_fio_do_real(fio, rotg->min_gaussian);
609 gmx_fio_do_real(fio, rotg->eps);
610 gmx_fio_do_int(fio, rotg->eFittype);
611 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
612 gmx_fio_do_real(fio, rotg->PotAngle_step);
615 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
619 gmx_fio_do_int(fio, rot->ngrp);
620 gmx_fio_do_int(fio, rot->nstrout);
621 gmx_fio_do_int(fio, rot->nstsout);
624 snew(rot->grp, rot->ngrp);
626 for (g = 0; g < rot->ngrp; g++)
628 do_rotgrp(fio, &rot->grp[g], bRead);
633 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
634 int file_version, real *fudgeQQ)
636 int i, j, k, *tmp, idum = 0;
640 real zerotemptime, finish_t, init_temp, finish_temp;
642 if (file_version != tpx_version)
644 /* Give a warning about features that are not accessible */
645 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
646 file_version, tpx_version);
654 if (file_version == 0)
659 /* Basic inputrec stuff */
660 gmx_fio_do_int(fio, ir->eI);
661 if (file_version >= 62)
663 gmx_fio_do_gmx_large_int(fio, ir->nsteps);
667 gmx_fio_do_int(fio, idum);
670 if (file_version > 25)
672 if (file_version >= 62)
674 gmx_fio_do_gmx_large_int(fio, ir->init_step);
678 gmx_fio_do_int(fio, idum);
679 ir->init_step = idum;
687 if (file_version >= 58)
689 gmx_fio_do_int(fio, ir->simulation_part);
693 ir->simulation_part = 1;
696 if (file_version >= 67)
698 gmx_fio_do_int(fio, ir->nstcalcenergy);
702 ir->nstcalcenergy = 1;
704 if (file_version < 53)
706 /* The pbc info has been moved out of do_inputrec,
707 * since we always want it, also without reading the inputrec.
709 gmx_fio_do_int(fio, ir->ePBC);
710 if ((file_version <= 15) && (ir->ePBC == 2))
714 if (file_version >= 45)
716 gmx_fio_do_int(fio, ir->bPeriodicMols);
723 ir->bPeriodicMols = TRUE;
727 ir->bPeriodicMols = FALSE;
731 if (file_version >= 81)
733 gmx_fio_do_int(fio, ir->cutoff_scheme);
737 ir->cutoff_scheme = ecutsGROUP;
739 gmx_fio_do_int(fio, ir->ns_type);
740 gmx_fio_do_int(fio, ir->nstlist);
741 gmx_fio_do_int(fio, ir->ndelta);
742 if (file_version < 41)
744 gmx_fio_do_int(fio, idum);
745 gmx_fio_do_int(fio, idum);
747 if (file_version >= 45)
749 gmx_fio_do_real(fio, ir->rtpi);
755 gmx_fio_do_int(fio, ir->nstcomm);
756 if (file_version > 34)
758 gmx_fio_do_int(fio, ir->comm_mode);
760 else if (ir->nstcomm < 0)
762 ir->comm_mode = ecmANGULAR;
766 ir->comm_mode = ecmLINEAR;
768 ir->nstcomm = abs(ir->nstcomm);
770 if (file_version > 25)
772 gmx_fio_do_int(fio, ir->nstcheckpoint);
776 ir->nstcheckpoint = 0;
779 gmx_fio_do_int(fio, ir->nstcgsteep);
781 if (file_version >= 30)
783 gmx_fio_do_int(fio, ir->nbfgscorr);
790 gmx_fio_do_int(fio, ir->nstlog);
791 gmx_fio_do_int(fio, ir->nstxout);
792 gmx_fio_do_int(fio, ir->nstvout);
793 gmx_fio_do_int(fio, ir->nstfout);
794 gmx_fio_do_int(fio, ir->nstenergy);
795 gmx_fio_do_int(fio, ir->nstxtcout);
796 if (file_version >= 59)
798 gmx_fio_do_double(fio, ir->init_t);
799 gmx_fio_do_double(fio, ir->delta_t);
803 gmx_fio_do_real(fio, rdum);
805 gmx_fio_do_real(fio, rdum);
808 gmx_fio_do_real(fio, ir->xtcprec);
809 if (file_version < 19)
811 gmx_fio_do_int(fio, idum);
812 gmx_fio_do_int(fio, idum);
814 if (file_version < 18)
816 gmx_fio_do_int(fio, idum);
818 if (file_version >= 81)
820 gmx_fio_do_real(fio, ir->verletbuf_drift);
824 ir->verletbuf_drift = 0;
826 gmx_fio_do_real(fio, ir->rlist);
827 if (file_version >= 67)
829 gmx_fio_do_real(fio, ir->rlistlong);
831 if (file_version >= 82 && file_version != 90)
833 gmx_fio_do_int(fio, ir->nstcalclr);
837 /* Calculate at NS steps */
838 ir->nstcalclr = ir->nstlist;
840 gmx_fio_do_int(fio, ir->coulombtype);
841 if (file_version < 32 && ir->coulombtype == eelRF)
843 ir->coulombtype = eelRF_NEC;
845 if (file_version >= 81)
847 gmx_fio_do_int(fio, ir->coulomb_modifier);
851 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
853 gmx_fio_do_real(fio, ir->rcoulomb_switch);
854 gmx_fio_do_real(fio, ir->rcoulomb);
855 gmx_fio_do_int(fio, ir->vdwtype);
856 if (file_version >= 81)
858 gmx_fio_do_int(fio, ir->vdw_modifier);
862 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
864 gmx_fio_do_real(fio, ir->rvdw_switch);
865 gmx_fio_do_real(fio, ir->rvdw);
866 if (file_version < 67)
868 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
870 gmx_fio_do_int(fio, ir->eDispCorr);
871 gmx_fio_do_real(fio, ir->epsilon_r);
872 if (file_version >= 37)
874 gmx_fio_do_real(fio, ir->epsilon_rf);
878 if (EEL_RF(ir->coulombtype))
880 ir->epsilon_rf = ir->epsilon_r;
885 ir->epsilon_rf = 1.0;
888 if (file_version >= 29)
890 gmx_fio_do_real(fio, ir->tabext);
897 if (file_version > 25)
899 gmx_fio_do_int(fio, ir->gb_algorithm);
900 gmx_fio_do_int(fio, ir->nstgbradii);
901 gmx_fio_do_real(fio, ir->rgbradii);
902 gmx_fio_do_real(fio, ir->gb_saltconc);
903 gmx_fio_do_int(fio, ir->implicit_solvent);
907 ir->gb_algorithm = egbSTILL;
911 ir->implicit_solvent = eisNO;
913 if (file_version >= 55)
915 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
916 gmx_fio_do_real(fio, ir->gb_obc_alpha);
917 gmx_fio_do_real(fio, ir->gb_obc_beta);
918 gmx_fio_do_real(fio, ir->gb_obc_gamma);
919 if (file_version >= 60)
921 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
922 gmx_fio_do_int(fio, ir->sa_algorithm);
926 ir->gb_dielectric_offset = 0.009;
927 ir->sa_algorithm = esaAPPROX;
929 gmx_fio_do_real(fio, ir->sa_surface_tension);
931 /* Override sa_surface_tension if it is not changed in the mpd-file */
932 if (ir->sa_surface_tension < 0)
934 if (ir->gb_algorithm == egbSTILL)
936 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
938 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
940 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
947 /* Better use sensible values than insane (0.0) ones... */
948 ir->gb_epsilon_solvent = 80;
949 ir->gb_obc_alpha = 1.0;
950 ir->gb_obc_beta = 0.8;
951 ir->gb_obc_gamma = 4.85;
952 ir->sa_surface_tension = 2.092;
956 if (file_version >= 81)
958 gmx_fio_do_real(fio, ir->fourier_spacing);
962 ir->fourier_spacing = 0.0;
964 gmx_fio_do_int(fio, ir->nkx);
965 gmx_fio_do_int(fio, ir->nky);
966 gmx_fio_do_int(fio, ir->nkz);
967 gmx_fio_do_int(fio, ir->pme_order);
968 gmx_fio_do_real(fio, ir->ewald_rtol);
970 if (file_version >= 24)
972 gmx_fio_do_int(fio, ir->ewald_geometry);
976 ir->ewald_geometry = eewg3D;
979 if (file_version <= 17)
981 ir->epsilon_surface = 0;
982 if (file_version == 17)
984 gmx_fio_do_int(fio, idum);
989 gmx_fio_do_real(fio, ir->epsilon_surface);
992 gmx_fio_do_gmx_bool(fio, ir->bOptFFT);
994 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
995 gmx_fio_do_int(fio, ir->etc);
996 /* before version 18, ir->etc was a gmx_bool (ir->btc),
997 * but the values 0 and 1 still mean no and
998 * berendsen temperature coupling, respectively.
1000 if (file_version >= 79)
1002 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1004 if (file_version >= 71)
1006 gmx_fio_do_int(fio, ir->nsttcouple);
1010 ir->nsttcouple = ir->nstcalcenergy;
1012 if (file_version <= 15)
1014 gmx_fio_do_int(fio, idum);
1016 if (file_version <= 17)
1018 gmx_fio_do_int(fio, ir->epct);
1019 if (file_version <= 15)
1023 ir->epct = epctSURFACETENSION;
1025 gmx_fio_do_int(fio, idum);
1028 /* we have removed the NO alternative at the beginning */
1032 ir->epct = epctISOTROPIC;
1036 ir->epc = epcBERENDSEN;
1041 gmx_fio_do_int(fio, ir->epc);
1042 gmx_fio_do_int(fio, ir->epct);
1044 if (file_version >= 71)
1046 gmx_fio_do_int(fio, ir->nstpcouple);
1050 ir->nstpcouple = ir->nstcalcenergy;
1052 gmx_fio_do_real(fio, ir->tau_p);
1053 if (file_version <= 15)
1055 gmx_fio_do_rvec(fio, vdum);
1056 clear_mat(ir->ref_p);
1057 for (i = 0; i < DIM; i++)
1059 ir->ref_p[i][i] = vdum[i];
1064 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1065 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1066 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1068 if (file_version <= 15)
1070 gmx_fio_do_rvec(fio, vdum);
1071 clear_mat(ir->compress);
1072 for (i = 0; i < DIM; i++)
1074 ir->compress[i][i] = vdum[i];
1079 gmx_fio_do_rvec(fio, ir->compress[XX]);
1080 gmx_fio_do_rvec(fio, ir->compress[YY]);
1081 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1083 if (file_version >= 47)
1085 gmx_fio_do_int(fio, ir->refcoord_scaling);
1086 gmx_fio_do_rvec(fio, ir->posres_com);
1087 gmx_fio_do_rvec(fio, ir->posres_comB);
1091 ir->refcoord_scaling = erscNO;
1092 clear_rvec(ir->posres_com);
1093 clear_rvec(ir->posres_comB);
1095 if ((file_version > 25) && (file_version < 79))
1097 gmx_fio_do_int(fio, ir->andersen_seed);
1101 ir->andersen_seed = 0;
1103 if (file_version < 26)
1105 gmx_fio_do_gmx_bool(fio, bSimAnn);
1106 gmx_fio_do_real(fio, zerotemptime);
1109 if (file_version < 37)
1111 gmx_fio_do_real(fio, rdum);
1114 gmx_fio_do_real(fio, ir->shake_tol);
1115 if (file_version < 54)
1117 gmx_fio_do_real(fio, *fudgeQQ);
1120 gmx_fio_do_int(fio, ir->efep);
1121 if (file_version <= 14 && ir->efep != efepNO)
1125 do_fepvals(fio, ir->fepvals, bRead, file_version);
1127 if (file_version >= 79)
1129 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1132 ir->bSimTemp = TRUE;
1137 ir->bSimTemp = FALSE;
1141 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1144 if (file_version >= 79)
1146 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1149 ir->bExpanded = TRUE;
1153 ir->bExpanded = FALSE;
1158 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1160 if (file_version >= 57)
1162 gmx_fio_do_int(fio, ir->eDisre);
1164 gmx_fio_do_int(fio, ir->eDisreWeighting);
1165 if (file_version < 22)
1167 if (ir->eDisreWeighting == 0)
1169 ir->eDisreWeighting = edrwEqual;
1173 ir->eDisreWeighting = edrwConservative;
1176 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1177 gmx_fio_do_real(fio, ir->dr_fc);
1178 gmx_fio_do_real(fio, ir->dr_tau);
1179 gmx_fio_do_int(fio, ir->nstdisreout);
1180 if (file_version >= 22)
1182 gmx_fio_do_real(fio, ir->orires_fc);
1183 gmx_fio_do_real(fio, ir->orires_tau);
1184 gmx_fio_do_int(fio, ir->nstorireout);
1190 ir->nstorireout = 0;
1192 if (file_version >= 26 && file_version < 79)
1194 gmx_fio_do_real(fio, ir->dihre_fc);
1195 if (file_version < 56)
1197 gmx_fio_do_real(fio, rdum);
1198 gmx_fio_do_int(fio, idum);
1206 gmx_fio_do_real(fio, ir->em_stepsize);
1207 gmx_fio_do_real(fio, ir->em_tol);
1208 if (file_version >= 22)
1210 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1214 ir->bShakeSOR = TRUE;
1216 if (file_version >= 11)
1218 gmx_fio_do_int(fio, ir->niter);
1223 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1226 if (file_version >= 21)
1228 gmx_fio_do_real(fio, ir->fc_stepsize);
1232 ir->fc_stepsize = 0;
1234 gmx_fio_do_int(fio, ir->eConstrAlg);
1235 gmx_fio_do_int(fio, ir->nProjOrder);
1236 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1237 if (file_version <= 14)
1239 gmx_fio_do_int(fio, idum);
1241 if (file_version >= 26)
1243 gmx_fio_do_int(fio, ir->nLincsIter);
1248 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1251 if (file_version < 33)
1253 gmx_fio_do_real(fio, bd_temp);
1255 gmx_fio_do_real(fio, ir->bd_fric);
1256 gmx_fio_do_int(fio, ir->ld_seed);
1257 if (file_version >= 33)
1259 for (i = 0; i < DIM; i++)
1261 gmx_fio_do_rvec(fio, ir->deform[i]);
1266 for (i = 0; i < DIM; i++)
1268 clear_rvec(ir->deform[i]);
1271 if (file_version >= 14)
1273 gmx_fio_do_real(fio, ir->cos_accel);
1279 gmx_fio_do_int(fio, ir->userint1);
1280 gmx_fio_do_int(fio, ir->userint2);
1281 gmx_fio_do_int(fio, ir->userint3);
1282 gmx_fio_do_int(fio, ir->userint4);
1283 gmx_fio_do_real(fio, ir->userreal1);
1284 gmx_fio_do_real(fio, ir->userreal2);
1285 gmx_fio_do_real(fio, ir->userreal3);
1286 gmx_fio_do_real(fio, ir->userreal4);
1289 if (file_version >= 77)
1291 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1296 snew(ir->adress, 1);
1298 gmx_fio_do_int(fio, ir->adress->type);
1299 gmx_fio_do_real(fio, ir->adress->const_wf);
1300 gmx_fio_do_real(fio, ir->adress->ex_width);
1301 gmx_fio_do_real(fio, ir->adress->hy_width);
1302 gmx_fio_do_int(fio, ir->adress->icor);
1303 gmx_fio_do_int(fio, ir->adress->site);
1304 gmx_fio_do_rvec(fio, ir->adress->refs);
1305 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1306 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1307 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1308 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1312 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1314 if (ir->adress->n_tf_grps > 0)
1316 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1320 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1322 if (ir->adress->n_energy_grps > 0)
1324 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1330 ir->bAdress = FALSE;
1334 if (file_version >= 48)
1336 gmx_fio_do_int(fio, ir->ePull);
1337 if (ir->ePull != epullNO)
1343 do_pull(fio, ir->pull, bRead, file_version);
1348 ir->ePull = epullNO;
1351 /* Enforced rotation */
1352 if (file_version >= 74)
1354 gmx_fio_do_int(fio, ir->bRot);
1355 if (ir->bRot == TRUE)
1361 do_rot(fio, ir->rot, bRead);
1370 gmx_fio_do_int(fio, ir->opts.ngtc);
1371 if (file_version >= 69)
1373 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1377 ir->opts.nhchainlength = 1;
1379 gmx_fio_do_int(fio, ir->opts.ngacc);
1380 gmx_fio_do_int(fio, ir->opts.ngfrz);
1381 gmx_fio_do_int(fio, ir->opts.ngener);
1385 snew(ir->opts.nrdf, ir->opts.ngtc);
1386 snew(ir->opts.ref_t, ir->opts.ngtc);
1387 snew(ir->opts.annealing, ir->opts.ngtc);
1388 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1389 snew(ir->opts.anneal_time, ir->opts.ngtc);
1390 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1391 snew(ir->opts.tau_t, ir->opts.ngtc);
1392 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1393 snew(ir->opts.acc, ir->opts.ngacc);
1394 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1396 if (ir->opts.ngtc > 0)
1398 if (bRead && file_version < 13)
1400 snew(tmp, ir->opts.ngtc);
1401 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1402 for (i = 0; i < ir->opts.ngtc; i++)
1404 ir->opts.nrdf[i] = tmp[i];
1410 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1412 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1413 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1414 if (file_version < 33 && ir->eI == eiBD)
1416 for (i = 0; i < ir->opts.ngtc; i++)
1418 ir->opts.tau_t[i] = bd_temp;
1422 if (ir->opts.ngfrz > 0)
1424 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1426 if (ir->opts.ngacc > 0)
1428 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1430 if (file_version >= 12)
1432 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1433 ir->opts.ngener*ir->opts.ngener);
1436 if (bRead && file_version < 26)
1438 for (i = 0; i < ir->opts.ngtc; i++)
1442 ir->opts.annealing[i] = eannSINGLE;
1443 ir->opts.anneal_npoints[i] = 2;
1444 snew(ir->opts.anneal_time[i], 2);
1445 snew(ir->opts.anneal_temp[i], 2);
1446 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1447 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1448 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1449 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1450 ir->opts.anneal_time[i][0] = ir->init_t;
1451 ir->opts.anneal_time[i][1] = finish_t;
1452 ir->opts.anneal_temp[i][0] = init_temp;
1453 ir->opts.anneal_temp[i][1] = finish_temp;
1457 ir->opts.annealing[i] = eannNO;
1458 ir->opts.anneal_npoints[i] = 0;
1464 /* file version 26 or later */
1465 /* First read the lists with annealing and npoints for each group */
1466 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1467 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1468 for (j = 0; j < (ir->opts.ngtc); j++)
1470 k = ir->opts.anneal_npoints[j];
1473 snew(ir->opts.anneal_time[j], k);
1474 snew(ir->opts.anneal_temp[j], k);
1476 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1477 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1481 if (file_version >= 45)
1483 gmx_fio_do_int(fio, ir->nwall);
1484 gmx_fio_do_int(fio, ir->wall_type);
1485 if (file_version >= 50)
1487 gmx_fio_do_real(fio, ir->wall_r_linpot);
1491 ir->wall_r_linpot = -1;
1493 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1494 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1495 gmx_fio_do_real(fio, ir->wall_density[0]);
1496 gmx_fio_do_real(fio, ir->wall_density[1]);
1497 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1503 ir->wall_atomtype[0] = -1;
1504 ir->wall_atomtype[1] = -1;
1505 ir->wall_density[0] = 0;
1506 ir->wall_density[1] = 0;
1507 ir->wall_ewald_zfac = 3;
1509 /* Cosine stuff for electric fields */
1510 for (j = 0; (j < DIM); j++)
1512 gmx_fio_do_int(fio, ir->ex[j].n);
1513 gmx_fio_do_int(fio, ir->et[j].n);
1516 snew(ir->ex[j].a, ir->ex[j].n);
1517 snew(ir->ex[j].phi, ir->ex[j].n);
1518 snew(ir->et[j].a, ir->et[j].n);
1519 snew(ir->et[j].phi, ir->et[j].n);
1521 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1522 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1523 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1524 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1528 if (file_version >= 39)
1530 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1531 gmx_fio_do_int(fio, ir->QMMMscheme);
1532 gmx_fio_do_real(fio, ir->scalefactor);
1533 gmx_fio_do_int(fio, ir->opts.ngQM);
1536 snew(ir->opts.QMmethod, ir->opts.ngQM);
1537 snew(ir->opts.QMbasis, ir->opts.ngQM);
1538 snew(ir->opts.QMcharge, ir->opts.ngQM);
1539 snew(ir->opts.QMmult, ir->opts.ngQM);
1540 snew(ir->opts.bSH, ir->opts.ngQM);
1541 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1542 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1543 snew(ir->opts.SAon, ir->opts.ngQM);
1544 snew(ir->opts.SAoff, ir->opts.ngQM);
1545 snew(ir->opts.SAsteps, ir->opts.ngQM);
1546 snew(ir->opts.bOPT, ir->opts.ngQM);
1547 snew(ir->opts.bTS, ir->opts.ngQM);
1549 if (ir->opts.ngQM > 0)
1551 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1552 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1553 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1554 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1555 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1556 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1557 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1558 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1559 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1560 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1561 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1562 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1564 /* end of QMMM stuff */
1569 static void do_harm(t_fileio *fio, t_iparams *iparams)
1571 gmx_fio_do_real(fio, iparams->harmonic.rA);
1572 gmx_fio_do_real(fio, iparams->harmonic.krA);
1573 gmx_fio_do_real(fio, iparams->harmonic.rB);
1574 gmx_fio_do_real(fio, iparams->harmonic.krB);
1577 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1578 gmx_bool bRead, int file_version)
1585 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1595 do_harm(fio, iparams);
1596 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1598 /* Correct incorrect storage of parameters */
1599 iparams->pdihs.phiB = iparams->pdihs.phiA;
1600 iparams->pdihs.cpB = iparams->pdihs.cpA;
1603 case F_LINEAR_ANGLES:
1604 gmx_fio_do_real(fio, iparams->linangle.klinA);
1605 gmx_fio_do_real(fio, iparams->linangle.aA);
1606 gmx_fio_do_real(fio, iparams->linangle.klinB);
1607 gmx_fio_do_real(fio, iparams->linangle.aB);
1610 gmx_fio_do_real(fio, iparams->fene.bm);
1611 gmx_fio_do_real(fio, iparams->fene.kb);
1614 gmx_fio_do_real(fio, iparams->restraint.lowA);
1615 gmx_fio_do_real(fio, iparams->restraint.up1A);
1616 gmx_fio_do_real(fio, iparams->restraint.up2A);
1617 gmx_fio_do_real(fio, iparams->restraint.kA);
1618 gmx_fio_do_real(fio, iparams->restraint.lowB);
1619 gmx_fio_do_real(fio, iparams->restraint.up1B);
1620 gmx_fio_do_real(fio, iparams->restraint.up2B);
1621 gmx_fio_do_real(fio, iparams->restraint.kB);
1627 gmx_fio_do_real(fio, iparams->tab.kA);
1628 gmx_fio_do_int(fio, iparams->tab.table);
1629 gmx_fio_do_real(fio, iparams->tab.kB);
1631 case F_CROSS_BOND_BONDS:
1632 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1633 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1634 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1636 case F_CROSS_BOND_ANGLES:
1637 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1638 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1639 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1640 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1642 case F_UREY_BRADLEY:
1643 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1644 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1645 gmx_fio_do_real(fio, iparams->u_b.r13A);
1646 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1647 if (file_version >= 79)
1649 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1650 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1651 gmx_fio_do_real(fio, iparams->u_b.r13B);
1652 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1656 iparams->u_b.thetaB = iparams->u_b.thetaA;
1657 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1658 iparams->u_b.r13B = iparams->u_b.r13A;
1659 iparams->u_b.kUBB = iparams->u_b.kUBA;
1662 case F_QUARTIC_ANGLES:
1663 gmx_fio_do_real(fio, iparams->qangle.theta);
1664 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1667 gmx_fio_do_real(fio, iparams->bham.a);
1668 gmx_fio_do_real(fio, iparams->bham.b);
1669 gmx_fio_do_real(fio, iparams->bham.c);
1672 gmx_fio_do_real(fio, iparams->morse.b0A);
1673 gmx_fio_do_real(fio, iparams->morse.cbA);
1674 gmx_fio_do_real(fio, iparams->morse.betaA);
1675 if (file_version >= 79)
1677 gmx_fio_do_real(fio, iparams->morse.b0B);
1678 gmx_fio_do_real(fio, iparams->morse.cbB);
1679 gmx_fio_do_real(fio, iparams->morse.betaB);
1683 iparams->morse.b0B = iparams->morse.b0A;
1684 iparams->morse.cbB = iparams->morse.cbA;
1685 iparams->morse.betaB = iparams->morse.betaA;
1689 gmx_fio_do_real(fio, iparams->cubic.b0);
1690 gmx_fio_do_real(fio, iparams->cubic.kb);
1691 gmx_fio_do_real(fio, iparams->cubic.kcub);
1695 case F_POLARIZATION:
1696 gmx_fio_do_real(fio, iparams->polarize.alpha);
1699 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1700 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1701 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1704 if (file_version < 31)
1706 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1708 gmx_fio_do_real(fio, iparams->wpol.al_x);
1709 gmx_fio_do_real(fio, iparams->wpol.al_y);
1710 gmx_fio_do_real(fio, iparams->wpol.al_z);
1711 gmx_fio_do_real(fio, iparams->wpol.rOH);
1712 gmx_fio_do_real(fio, iparams->wpol.rHH);
1713 gmx_fio_do_real(fio, iparams->wpol.rOD);
1716 gmx_fio_do_real(fio, iparams->thole.a);
1717 gmx_fio_do_real(fio, iparams->thole.alpha1);
1718 gmx_fio_do_real(fio, iparams->thole.alpha2);
1719 gmx_fio_do_real(fio, iparams->thole.rfac);
1722 gmx_fio_do_real(fio, iparams->lj.c6);
1723 gmx_fio_do_real(fio, iparams->lj.c12);
1726 gmx_fio_do_real(fio, iparams->lj14.c6A);
1727 gmx_fio_do_real(fio, iparams->lj14.c12A);
1728 gmx_fio_do_real(fio, iparams->lj14.c6B);
1729 gmx_fio_do_real(fio, iparams->lj14.c12B);
1732 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1733 gmx_fio_do_real(fio, iparams->ljc14.qi);
1734 gmx_fio_do_real(fio, iparams->ljc14.qj);
1735 gmx_fio_do_real(fio, iparams->ljc14.c6);
1736 gmx_fio_do_real(fio, iparams->ljc14.c12);
1738 case F_LJC_PAIRS_NB:
1739 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1740 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1741 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1742 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1748 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1749 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1750 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1752 /* Read the incorrectly stored multiplicity */
1753 gmx_fio_do_real(fio, iparams->harmonic.rB);
1754 gmx_fio_do_real(fio, iparams->harmonic.krB);
1755 iparams->pdihs.phiB = iparams->pdihs.phiA;
1756 iparams->pdihs.cpB = iparams->pdihs.cpA;
1760 gmx_fio_do_real(fio, iparams->pdihs.phiB);
1761 gmx_fio_do_real(fio, iparams->pdihs.cpB);
1762 gmx_fio_do_int(fio, iparams->pdihs.mult);
1766 gmx_fio_do_int(fio, iparams->disres.label);
1767 gmx_fio_do_int(fio, iparams->disres.type);
1768 gmx_fio_do_real(fio, iparams->disres.low);
1769 gmx_fio_do_real(fio, iparams->disres.up1);
1770 gmx_fio_do_real(fio, iparams->disres.up2);
1771 gmx_fio_do_real(fio, iparams->disres.kfac);
1774 gmx_fio_do_int(fio, iparams->orires.ex);
1775 gmx_fio_do_int(fio, iparams->orires.label);
1776 gmx_fio_do_int(fio, iparams->orires.power);
1777 gmx_fio_do_real(fio, iparams->orires.c);
1778 gmx_fio_do_real(fio, iparams->orires.obs);
1779 gmx_fio_do_real(fio, iparams->orires.kfac);
1782 if (file_version < 82)
1784 gmx_fio_do_int(fio, idum);
1785 gmx_fio_do_int(fio, idum);
1787 gmx_fio_do_real(fio, iparams->dihres.phiA);
1788 gmx_fio_do_real(fio, iparams->dihres.dphiA);
1789 gmx_fio_do_real(fio, iparams->dihres.kfacA);
1790 if (file_version >= 82)
1792 gmx_fio_do_real(fio, iparams->dihres.phiB);
1793 gmx_fio_do_real(fio, iparams->dihres.dphiB);
1794 gmx_fio_do_real(fio, iparams->dihres.kfacB);
1798 iparams->dihres.phiB = iparams->dihres.phiA;
1799 iparams->dihres.dphiB = iparams->dihres.dphiA;
1800 iparams->dihres.kfacB = iparams->dihres.kfacA;
1804 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
1805 gmx_fio_do_rvec(fio, iparams->posres.fcA);
1806 if (bRead && file_version < 27)
1808 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
1809 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
1813 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
1814 gmx_fio_do_rvec(fio, iparams->posres.fcB);
1818 gmx_fio_do_int(fio, iparams->fbposres.geom);
1819 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
1820 gmx_fio_do_real(fio, iparams->fbposres.r);
1821 gmx_fio_do_real(fio, iparams->fbposres.k);
1824 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1825 if (file_version >= 25)
1827 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1831 /* Fourier dihedrals are internally represented
1832 * as Ryckaert-Bellemans since those are faster to compute.
1834 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1835 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1839 gmx_fio_do_real(fio, iparams->constr.dA);
1840 gmx_fio_do_real(fio, iparams->constr.dB);
1843 gmx_fio_do_real(fio, iparams->settle.doh);
1844 gmx_fio_do_real(fio, iparams->settle.dhh);
1847 gmx_fio_do_real(fio, iparams->vsite.a);
1852 gmx_fio_do_real(fio, iparams->vsite.a);
1853 gmx_fio_do_real(fio, iparams->vsite.b);
1858 gmx_fio_do_real(fio, iparams->vsite.a);
1859 gmx_fio_do_real(fio, iparams->vsite.b);
1860 gmx_fio_do_real(fio, iparams->vsite.c);
1863 gmx_fio_do_int(fio, iparams->vsiten.n);
1864 gmx_fio_do_real(fio, iparams->vsiten.a);
1869 /* We got rid of some parameters in version 68 */
1870 if (bRead && file_version < 68)
1872 gmx_fio_do_real(fio, rdum);
1873 gmx_fio_do_real(fio, rdum);
1874 gmx_fio_do_real(fio, rdum);
1875 gmx_fio_do_real(fio, rdum);
1877 gmx_fio_do_real(fio, iparams->gb.sar);
1878 gmx_fio_do_real(fio, iparams->gb.st);
1879 gmx_fio_do_real(fio, iparams->gb.pi);
1880 gmx_fio_do_real(fio, iparams->gb.gbr);
1881 gmx_fio_do_real(fio, iparams->gb.bmlt);
1884 gmx_fio_do_int(fio, iparams->cmap.cmapA);
1885 gmx_fio_do_int(fio, iparams->cmap.cmapB);
1888 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1889 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1893 gmx_fio_unset_comment(fio);
1897 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
1904 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1906 if (file_version < 44)
1908 for (i = 0; i < MAXNODES; i++)
1910 gmx_fio_do_int(fio, idum);
1913 gmx_fio_do_int(fio, ilist->nr);
1916 snew(ilist->iatoms, ilist->nr);
1918 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
1921 gmx_fio_unset_comment(fio);
1925 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1926 gmx_bool bRead, int file_version)
1931 gmx_fio_do_int(fio, ffparams->atnr);
1932 if (file_version < 57)
1934 gmx_fio_do_int(fio, idum);
1936 gmx_fio_do_int(fio, ffparams->ntypes);
1939 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
1940 ffparams->atnr, ffparams->ntypes);
1944 snew(ffparams->functype, ffparams->ntypes);
1945 snew(ffparams->iparams, ffparams->ntypes);
1947 /* Read/write all the function types */
1948 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
1951 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
1954 if (file_version >= 66)
1956 gmx_fio_do_double(fio, ffparams->reppow);
1960 ffparams->reppow = 12.0;
1963 if (file_version >= 57)
1965 gmx_fio_do_real(fio, ffparams->fudgeQQ);
1968 /* Check whether all these function types are supported by the code.
1969 * In practice the code is backwards compatible, which means that the
1970 * numbering may have to be altered from old numbering to new numbering
1972 for (i = 0; (i < ffparams->ntypes); i++)
1976 /* Loop over file versions */
1977 for (k = 0; (k < NFTUPD); k++)
1979 /* Compare the read file_version to the update table */
1980 if ((file_version < ftupd[k].fvnr) &&
1981 (ffparams->functype[i] >= ftupd[k].ftype))
1983 ffparams->functype[i] += 1;
1986 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
1987 i, ffparams->functype[i],
1988 interaction_function[ftupd[k].ftype].longname);
1995 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
1999 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2004 static void add_settle_atoms(t_ilist *ilist)
2008 /* Settle used to only store the first atom: add the other two */
2009 srenew(ilist->iatoms, 2*ilist->nr);
2010 for (i = ilist->nr/2-1; i >= 0; i--)
2012 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2013 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2014 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2015 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2017 ilist->nr = 2*ilist->nr;
2020 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2023 int i, j, renum[F_NRE];
2027 for (j = 0; (j < F_NRE); j++)
2032 for (k = 0; k < NFTUPD; k++)
2034 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2043 ilist[j].iatoms = NULL;
2047 do_ilist(fio, &ilist[j], bRead, file_version, j);
2048 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2050 add_settle_atoms(&ilist[j]);
2054 if (bRead && gmx_debug_at)
2055 pr_ilist(debug,0,interaction_function[j].longname,
2056 functype,&ilist[j],TRUE);
2061 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2062 gmx_bool bRead, int file_version)
2064 do_ffparams(fio, ffparams, bRead, file_version);
2066 if (file_version >= 54)
2068 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2071 do_ilists(fio, molt->ilist, bRead, file_version);
2074 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2076 int i, idum, dum_nra, *dum_a;
2078 if (file_version < 44)
2080 for (i = 0; i < MAXNODES; i++)
2082 gmx_fio_do_int(fio, idum);
2085 gmx_fio_do_int(fio, block->nr);
2086 if (file_version < 51)
2088 gmx_fio_do_int(fio, dum_nra);
2092 if ((block->nalloc_index > 0) && (NULL != block->index))
2094 sfree(block->index);
2096 block->nalloc_index = block->nr+1;
2097 snew(block->index, block->nalloc_index);
2099 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2101 if (file_version < 51 && dum_nra > 0)
2103 snew(dum_a, dum_nra);
2104 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2109 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2114 if (file_version < 44)
2116 for (i = 0; i < MAXNODES; i++)
2118 gmx_fio_do_int(fio, idum);
2121 gmx_fio_do_int(fio, block->nr);
2122 gmx_fio_do_int(fio, block->nra);
2125 block->nalloc_index = block->nr+1;
2126 snew(block->index, block->nalloc_index);
2127 block->nalloc_a = block->nra;
2128 snew(block->a, block->nalloc_a);
2130 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2131 gmx_fio_ndo_int(fio, block->a, block->nra);
2134 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2135 int file_version, gmx_groups_t *groups, int atnr)
2139 gmx_fio_do_real(fio, atom->m);
2140 gmx_fio_do_real(fio, atom->q);
2141 gmx_fio_do_real(fio, atom->mB);
2142 gmx_fio_do_real(fio, atom->qB);
2143 gmx_fio_do_ushort(fio, atom->type);
2144 gmx_fio_do_ushort(fio, atom->typeB);
2145 gmx_fio_do_int(fio, atom->ptype);
2146 gmx_fio_do_int(fio, atom->resind);
2147 if (file_version >= 52)
2149 gmx_fio_do_int(fio, atom->atomnumber);
2153 atom->atomnumber = NOTSET;
2155 if (file_version < 23)
2159 else if (file_version < 39)
2168 if (file_version < 57)
2170 unsigned char uchar[egcNR];
2171 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2172 for (i = myngrp; (i < ngrp); i++)
2176 /* Copy the old data format to the groups struct */
2177 for (i = 0; i < ngrp; i++)
2179 groups->grpnr[i][atnr] = uchar[i];
2184 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2189 if (file_version < 23)
2193 else if (file_version < 39)
2202 for (j = 0; (j < ngrp); j++)
2206 gmx_fio_do_int(fio, grps[j].nr);
2209 snew(grps[j].nm_ind, grps[j].nr);
2211 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2216 snew(grps[j].nm_ind, grps[j].nr);
2221 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2227 gmx_fio_do_int(fio, ls);
2228 *nm = get_symtab_handle(symtab, ls);
2232 ls = lookup_symtab(symtab, *nm);
2233 gmx_fio_do_int(fio, ls);
2237 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2242 for (j = 0; (j < nstr); j++)
2244 do_symstr(fio, &(nm[j]), bRead, symtab);
2248 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2249 t_symtab *symtab, int file_version)
2253 for (j = 0; (j < n); j++)
2255 do_symstr(fio, &(ri[j].name), bRead, symtab);
2256 if (file_version >= 63)
2258 gmx_fio_do_int(fio, ri[j].nr);
2259 gmx_fio_do_uchar(fio, ri[j].ic);
2269 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2271 gmx_groups_t *groups)
2275 gmx_fio_do_int(fio, atoms->nr);
2276 gmx_fio_do_int(fio, atoms->nres);
2277 if (file_version < 57)
2279 gmx_fio_do_int(fio, groups->ngrpname);
2280 for (i = 0; i < egcNR; i++)
2282 groups->ngrpnr[i] = atoms->nr;
2283 snew(groups->grpnr[i], groups->ngrpnr[i]);
2288 snew(atoms->atom, atoms->nr);
2289 snew(atoms->atomname, atoms->nr);
2290 snew(atoms->atomtype, atoms->nr);
2291 snew(atoms->atomtypeB, atoms->nr);
2292 snew(atoms->resinfo, atoms->nres);
2293 if (file_version < 57)
2295 snew(groups->grpname, groups->ngrpname);
2297 atoms->pdbinfo = NULL;
2299 for (i = 0; (i < atoms->nr); i++)
2301 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2303 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2304 if (bRead && (file_version <= 20))
2306 for (i = 0; i < atoms->nr; i++)
2308 atoms->atomtype[i] = put_symtab(symtab, "?");
2309 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2314 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2315 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2317 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2319 if (file_version < 57)
2321 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2323 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2327 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2328 gmx_bool bRead, t_symtab *symtab,
2333 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2334 gmx_fio_do_int(fio, groups->ngrpname);
2337 snew(groups->grpname, groups->ngrpname);
2339 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2340 for (g = 0; g < egcNR; g++)
2342 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2343 if (groups->ngrpnr[g] == 0)
2347 groups->grpnr[g] = NULL;
2354 snew(groups->grpnr[g], groups->ngrpnr[g]);
2356 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2361 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2366 if (file_version > 25)
2368 gmx_fio_do_int(fio, atomtypes->nr);
2372 snew(atomtypes->radius, j);
2373 snew(atomtypes->vol, j);
2374 snew(atomtypes->surftens, j);
2375 snew(atomtypes->atomnumber, j);
2376 snew(atomtypes->gb_radius, j);
2377 snew(atomtypes->S_hct, j);
2379 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2380 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2381 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2382 if (file_version >= 40)
2384 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2386 if (file_version >= 60)
2388 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2389 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2394 /* File versions prior to 26 cannot do GBSA,
2395 * so they dont use this structure
2398 atomtypes->radius = NULL;
2399 atomtypes->vol = NULL;
2400 atomtypes->surftens = NULL;
2401 atomtypes->atomnumber = NULL;
2402 atomtypes->gb_radius = NULL;
2403 atomtypes->S_hct = NULL;
2407 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2413 gmx_fio_do_int(fio, symtab->nr);
2417 snew(symtab->symbuf, 1);
2418 symbuf = symtab->symbuf;
2419 symbuf->bufsize = nr;
2420 snew(symbuf->buf, nr);
2421 for (i = 0; (i < nr); i++)
2423 gmx_fio_do_string(fio, buf);
2424 symbuf->buf[i] = strdup(buf);
2429 symbuf = symtab->symbuf;
2430 while (symbuf != NULL)
2432 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2434 gmx_fio_do_string(fio, symbuf->buf[i]);
2437 symbuf = symbuf->next;
2441 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2446 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2448 int i, j, ngrid, gs, nelem;
2450 gmx_fio_do_int(fio, cmap_grid->ngrid);
2451 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2453 ngrid = cmap_grid->ngrid;
2454 gs = cmap_grid->grid_spacing;
2459 snew(cmap_grid->cmapdata, ngrid);
2461 for (i = 0; i < cmap_grid->ngrid; i++)
2463 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2467 for (i = 0; i < cmap_grid->ngrid; i++)
2469 for (j = 0; j < nelem; j++)
2471 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2472 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2473 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2474 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2480 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2482 int m, a, a0, a1, r;
2486 /* We always assign a new chain number, but save the chain id characters
2487 * for larger molecules.
2489 #define CHAIN_MIN_ATOMS 15
2493 for (m = 0; m < mols->nr; m++)
2495 a0 = mols->index[m];
2496 a1 = mols->index[m+1];
2497 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2506 for (a = a0; a < a1; a++)
2508 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2509 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2514 /* Blank out the chain id if there was only one chain */
2517 for (r = 0; r < atoms->nres; r++)
2519 atoms->resinfo[r].chainid = ' ';
2524 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2525 t_symtab *symtab, int file_version,
2526 gmx_groups_t *groups)
2530 if (file_version >= 57)
2532 do_symstr(fio, &(molt->name), bRead, symtab);
2535 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2537 if (bRead && gmx_debug_at)
2539 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2542 if (file_version >= 57)
2544 do_ilists(fio, molt->ilist, bRead, file_version);
2546 do_block(fio, &molt->cgs, bRead, file_version);
2547 if (bRead && gmx_debug_at)
2549 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2553 /* This used to be in the atoms struct */
2554 do_blocka(fio, &molt->excls, bRead, file_version);
2557 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2561 gmx_fio_do_int(fio, molb->type);
2562 gmx_fio_do_int(fio, molb->nmol);
2563 gmx_fio_do_int(fio, molb->natoms_mol);
2564 /* Position restraint coordinates */
2565 gmx_fio_do_int(fio, molb->nposres_xA);
2566 if (molb->nposres_xA > 0)
2570 snew(molb->posres_xA, molb->nposres_xA);
2572 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2574 gmx_fio_do_int(fio, molb->nposres_xB);
2575 if (molb->nposres_xB > 0)
2579 snew(molb->posres_xB, molb->nposres_xB);
2581 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2586 static t_block mtop_mols(gmx_mtop_t *mtop)
2592 for (mb = 0; mb < mtop->nmolblock; mb++)
2594 mols.nr += mtop->molblock[mb].nmol;
2596 mols.nalloc_index = mols.nr + 1;
2597 snew(mols.index, mols.nalloc_index);
2602 for (mb = 0; mb < mtop->nmolblock; mb++)
2604 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2606 a += mtop->molblock[mb].natoms_mol;
2615 static void add_posres_molblock(gmx_mtop_t *mtop)
2620 gmx_molblock_t *molb;
2623 /* posres reference positions are stored in ip->posres (if present) and
2624 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2625 posres.pos0A are identical to fbposres.pos0. */
2626 il = &mtop->moltype[0].ilist[F_POSRES];
2627 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2628 if (il->nr == 0 && ilfb->nr == 0)
2634 for (i = 0; i < il->nr; i += 2)
2636 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2637 am = max(am, il->iatoms[i+1]);
2638 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2639 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2640 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2645 /* This loop is required if we have only flat-bottomed posres:
2647 - bFE == FALSE (no B-state for flat-bottomed posres) */
2650 for (i = 0; i < ilfb->nr; i += 2)
2652 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2653 am = max(am, ilfb->iatoms[i+1]);
2656 /* Make the posres coordinate block end at a molecule end */
2658 while (am >= mtop->mols.index[mol+1])
2662 molb = &mtop->molblock[0];
2663 molb->nposres_xA = mtop->mols.index[mol+1];
2664 snew(molb->posres_xA, molb->nposres_xA);
2667 molb->nposres_xB = molb->nposres_xA;
2668 snew(molb->posres_xB, molb->nposres_xB);
2672 molb->nposres_xB = 0;
2674 for (i = 0; i < il->nr; i += 2)
2676 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2677 a = il->iatoms[i+1];
2678 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2679 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2680 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2683 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2684 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2685 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2690 /* If only flat-bottomed posres are present, take reference pos from them.
2691 Here: bFE == FALSE */
2692 for (i = 0; i < ilfb->nr; i += 2)
2694 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2695 a = ilfb->iatoms[i+1];
2696 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2697 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2698 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2703 static void set_disres_npair(gmx_mtop_t *mtop)
2710 ip = mtop->ffparams.iparams;
2712 for (mt = 0; mt < mtop->nmoltype; mt++)
2714 il = &mtop->moltype[mt].ilist[F_DISRES];
2719 for (i = 0; i < il->nr; i += 3)
2722 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2724 ip[a[i]].disres.npair = npair;
2732 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2742 do_symtab(fio, &(mtop->symtab), bRead);
2745 pr_symtab(debug, 0, "symtab", &mtop->symtab);
2748 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2750 if (file_version >= 57)
2752 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
2754 gmx_fio_do_int(fio, mtop->nmoltype);
2762 snew(mtop->moltype, mtop->nmoltype);
2763 if (file_version < 57)
2765 mtop->moltype[0].name = mtop->name;
2768 for (mt = 0; mt < mtop->nmoltype; mt++)
2770 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
2774 if (file_version >= 57)
2776 gmx_fio_do_int(fio, mtop->nmolblock);
2780 mtop->nmolblock = 1;
2784 snew(mtop->molblock, mtop->nmolblock);
2786 if (file_version >= 57)
2788 for (mb = 0; mb < mtop->nmolblock; mb++)
2790 do_molblock(fio, &mtop->molblock[mb], bRead);
2792 gmx_fio_do_int(fio, mtop->natoms);
2796 mtop->molblock[0].type = 0;
2797 mtop->molblock[0].nmol = 1;
2798 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
2799 mtop->molblock[0].nposres_xA = 0;
2800 mtop->molblock[0].nposres_xB = 0;
2803 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
2806 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
2809 if (file_version < 57)
2811 /* Debug statements are inside do_idef */
2812 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
2813 mtop->natoms = mtop->moltype[0].atoms.nr;
2816 if (file_version >= 65)
2818 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
2822 mtop->ffparams.cmap_grid.ngrid = 0;
2823 mtop->ffparams.cmap_grid.grid_spacing = 0;
2824 mtop->ffparams.cmap_grid.cmapdata = NULL;
2827 if (file_version >= 57)
2829 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
2832 if (file_version < 57)
2834 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
2835 if (bRead && gmx_debug_at)
2837 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
2839 do_block(fio, &mtop->mols, bRead, file_version);
2840 /* Add the posres coordinates to the molblock */
2841 add_posres_molblock(mtop);
2845 if (file_version >= 57)
2847 done_block(&mtop->mols);
2848 mtop->mols = mtop_mols(mtop);
2852 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
2856 if (file_version < 51)
2858 /* Here used to be the shake blocks */
2859 do_blocka(fio, &dumb, bRead, file_version);
2872 close_symtab(&(mtop->symtab));
2876 /* If TopOnlyOK is TRUE then we can read even future versions
2877 * of tpx files, provided the file_generation hasn't changed.
2878 * If it is FALSE, we need the inputrecord too, and bail out
2879 * if the file is newer than the program.
2881 * The version and generation if the topology (see top of this file)
2882 * are returned in the two last arguments.
2884 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2886 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
2887 gmx_bool TopOnlyOK, int *file_version,
2888 int *file_generation)
2891 char file_tag[STRLEN];
2898 gmx_fio_checktype(fio);
2899 gmx_fio_setdebug(fio, bDebugMode());
2901 /* NEW! XDR tpb file */
2902 precision = sizeof(real);
2905 gmx_fio_do_string(fio, buf);
2906 if (strncmp(buf, "VERSION", 7))
2908 gmx_fatal(FARGS, "Can not read file %s,\n"
2909 " this file is from a Gromacs version which is older than 2.0\n"
2910 " Make a new one with grompp or use a gro or pdb file, if possible",
2911 gmx_fio_getname(fio));
2913 gmx_fio_do_int(fio, precision);
2914 bDouble = (precision == sizeof(double));
2915 if ((precision != sizeof(float)) && !bDouble)
2917 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
2918 "instead of %d or %d",
2919 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
2921 gmx_fio_setprecision(fio, bDouble);
2922 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
2923 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
2927 gmx_fio_write_string(fio, GromacsVersion());
2928 bDouble = (precision == sizeof(double));
2929 gmx_fio_setprecision(fio, bDouble);
2930 gmx_fio_do_int(fio, precision);
2932 sprintf(file_tag, "%s", tpx_tag);
2933 fgen = tpx_generation;
2936 /* Check versions! */
2937 gmx_fio_do_int(fio, fver);
2939 /* This is for backward compatibility with development versions 77-79
2940 * where the tag was, mistakenly, placed before the generation,
2941 * which would cause a segv instead of a proper error message
2942 * when reading the topology only from tpx with <77 code.
2944 if (fver >= 77 && fver <= 79)
2946 gmx_fio_do_string(fio, file_tag);
2951 gmx_fio_do_int(fio, fgen);
2960 gmx_fio_do_string(fio, file_tag);
2966 /* Versions before 77 don't have the tag, set it to release */
2967 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
2970 if (strcmp(file_tag, tpx_tag) != 0)
2972 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
2975 /* We only support reading tpx files with the same tag as the code
2976 * or tpx files with the release tag and with lower version number.
2978 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
2980 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2981 gmx_fio_getname(fio), fver, file_tag,
2982 tpx_version, tpx_tag);
2987 if (file_version != NULL)
2989 *file_version = fver;
2991 if (file_generation != NULL)
2993 *file_generation = fgen;
2997 if ((fver <= tpx_incompatible_version) ||
2998 ((fver > tpx_version) && !TopOnlyOK) ||
2999 (fgen > tpx_generation) ||
3000 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3002 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3003 gmx_fio_getname(fio), fver, tpx_version);
3006 do_section(fio, eitemHEADER, bRead);
3007 gmx_fio_do_int(fio, tpx->natoms);
3010 gmx_fio_do_int(fio, tpx->ngtc);
3018 gmx_fio_do_int(fio, idum);
3019 gmx_fio_do_real(fio, rdum);
3021 /*a better decision will eventually (5.0 or later) need to be made
3022 on how to treat the alchemical state of the system, which can now
3023 vary through a simulation, and cannot be completely described
3024 though a single lambda variable, or even a single state
3025 index. Eventually, should probably be a vector. MRS*/
3028 gmx_fio_do_int(fio, tpx->fep_state);
3030 gmx_fio_do_real(fio, tpx->lambda);
3031 gmx_fio_do_int(fio, tpx->bIr);
3032 gmx_fio_do_int(fio, tpx->bTop);
3033 gmx_fio_do_int(fio, tpx->bX);
3034 gmx_fio_do_int(fio, tpx->bV);
3035 gmx_fio_do_int(fio, tpx->bF);
3036 gmx_fio_do_int(fio, tpx->bBox);
3038 if ((fgen > tpx_generation))
3040 /* This can only happen if TopOnlyOK=TRUE */
3045 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3046 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3047 gmx_bool bXVallocated)
3053 int file_version, file_generation;
3057 gmx_bool bPeriodicMols;
3061 tpx.natoms = state->natoms;
3062 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3063 tpx.fep_state = state->fep_state;
3064 tpx.lambda = state->lambda[efptFEP];
3065 tpx.bIr = (ir != NULL);
3066 tpx.bTop = (mtop != NULL);
3067 tpx.bX = (state->x != NULL);
3068 tpx.bV = (state->v != NULL);
3069 tpx.bF = (f != NULL);
3073 TopOnlyOK = (ir == NULL);
3075 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3080 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3081 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3086 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3087 state->natoms = tpx.natoms;
3088 state->nalloc = tpx.natoms;
3094 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3098 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3100 do_test(fio, tpx.bBox, state->box);
3101 do_section(fio, eitemBOX, bRead);
3104 gmx_fio_ndo_rvec(fio, state->box, DIM);
3105 if (file_version >= 51)
3107 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3111 /* We initialize box_rel after reading the inputrec */
3112 clear_mat(state->box_rel);
3114 if (file_version >= 28)
3116 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3117 if (file_version < 56)
3120 gmx_fio_ndo_rvec(fio, mdum, DIM);
3125 if (state->ngtc > 0 && file_version >= 28)
3128 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3129 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3130 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3131 snew(dumv, state->ngtc);
3132 if (file_version < 69)
3134 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3136 /* These used to be the Berendsen tcoupl_lambda's */
3137 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3141 /* Prior to tpx version 26, the inputrec was here.
3142 * I moved it to enable partial forward-compatibility
3143 * for analysis/viewer programs.
3145 if (file_version < 26)
3147 do_test(fio, tpx.bIr, ir);
3148 do_section(fio, eitemIR, bRead);
3153 do_inputrec(fio, ir, bRead, file_version,
3154 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3157 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3162 do_inputrec(fio, &dum_ir, bRead, file_version,
3163 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3166 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3168 done_inputrec(&dum_ir);
3174 do_test(fio, tpx.bTop, mtop);
3175 do_section(fio, eitemTOP, bRead);
3180 do_mtop(fio, mtop, bRead, file_version);
3184 do_mtop(fio, &dum_top, bRead, file_version);
3185 done_mtop(&dum_top, TRUE);
3188 do_test(fio, tpx.bX, state->x);
3189 do_section(fio, eitemX, bRead);
3194 state->flags |= (1<<estX);
3196 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3199 do_test(fio, tpx.bV, state->v);
3200 do_section(fio, eitemV, bRead);
3205 state->flags |= (1<<estV);
3207 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3210 do_test(fio, tpx.bF, f);
3211 do_section(fio, eitemF, bRead);
3214 gmx_fio_ndo_rvec(fio, f, state->natoms);
3217 /* Starting with tpx version 26, we have the inputrec
3218 * at the end of the file, so we can ignore it
3219 * if the file is never than the software (but still the
3220 * same generation - see comments at the top of this file.
3225 bPeriodicMols = FALSE;
3226 if (file_version >= 26)
3228 do_test(fio, tpx.bIr, ir);
3229 do_section(fio, eitemIR, bRead);
3232 if (file_version >= 53)
3234 /* Removed the pbc info from do_inputrec, since we always want it */
3238 bPeriodicMols = ir->bPeriodicMols;
3240 gmx_fio_do_int(fio, ePBC);
3241 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3243 if (file_generation <= tpx_generation && ir)
3245 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3248 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3250 if (file_version < 51)
3252 set_box_rel(ir, state);
3254 if (file_version < 53)
3257 bPeriodicMols = ir->bPeriodicMols;
3260 if (bRead && ir && file_version >= 53)
3262 /* We need to do this after do_inputrec, since that initializes ir */
3264 ir->bPeriodicMols = bPeriodicMols;
3273 if (state->ngtc == 0)
3275 /* Reading old version without tcoupl state data: set it */
3276 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3278 if (tpx.bTop && mtop)
3280 if (file_version < 57)
3282 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3284 ir->eDisre = edrSimple;
3288 ir->eDisre = edrNone;
3291 set_disres_npair(mtop);
3295 if (tpx.bTop && mtop)
3297 gmx_mtop_finalize(mtop);
3300 if (file_version >= 57)
3304 env = getenv("GMX_NOCHARGEGROUPS");
3307 sscanf(env, "%d", &ienv);
3308 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3313 "Will make single atomic charge groups in non-solvent%s\n",
3314 ienv > 1 ? " and solvent" : "");
3315 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3317 fprintf(stderr, "\n");
3325 /************************************************************
3327 * The following routines are the exported ones
3329 ************************************************************/
3331 t_fileio *open_tpx(const char *fn, const char *mode)
3333 return gmx_fio_open(fn, mode);
3336 void close_tpx(t_fileio *fio)
3341 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3342 int *file_version, int *file_generation)
3346 fio = open_tpx(fn, "r");
3347 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3351 void write_tpx_state(const char *fn,
3352 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3356 fio = open_tpx(fn, "w");
3357 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3361 void read_tpx_state(const char *fn,
3362 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3366 fio = open_tpx(fn, "r");
3367 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3371 int read_tpx(const char *fn,
3372 t_inputrec *ir, matrix box, int *natoms,
3373 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3381 fio = open_tpx(fn, "r");
3382 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3384 *natoms = state.natoms;
3387 copy_mat(state.box, box);
3396 int read_tpx_top(const char *fn,
3397 t_inputrec *ir, matrix box, int *natoms,
3398 rvec *x, rvec *v, rvec *f, t_topology *top)
3404 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3406 *top = gmx_mtop_t_to_t_topology(&mtop);
3411 gmx_bool fn2bTPX(const char *file)
3413 switch (fn2ftp(file))
3424 static void done_gmx_groups_t(gmx_groups_t *g)
3428 for (i = 0; (i < egcNR); i++)
3430 if (NULL != g->grps[i].nm_ind)
3432 sfree(g->grps[i].nm_ind);
3433 g->grps[i].nm_ind = NULL;
3435 if (NULL != g->grpnr[i])
3441 /* The contents of this array is in symtab, don't free it here */
3445 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3446 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3449 int natoms, i, version, generation;
3450 gmx_bool bTop, bXNULL = FALSE;
3452 t_topology *topconv;
3455 bTop = fn2bTPX(infile);
3459 read_tpxheader(infile, &header, TRUE, &version, &generation);
3462 snew(*x, header.natoms);
3466 snew(*v, header.natoms);
3469 *ePBC = read_tpx(infile, NULL, box, &natoms,
3470 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3471 *top = gmx_mtop_t_to_t_topology(mtop);
3472 /* In this case we need to throw away the group data too */
3473 done_gmx_groups_t(&mtop->groups);
3475 strcpy(title, *top->name);
3476 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3480 get_stx_coordnum(infile, &natoms);
3481 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3492 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3500 aps = gmx_atomprop_init();
3501 for (i = 0; (i < natoms); i++)
3503 if (!gmx_atomprop_query(aps, epropMass,
3504 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3505 *top->atoms.atomname[i],
3506 &(top->atoms.atom[i].m)))
3510 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3511 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3512 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3513 *top->atoms.atomname[i]);
3517 gmx_atomprop_destroy(aps);
3519 top->idef.ntypes = -1;