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,2014, 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"
60 #include "mtop_util.h"
62 #define TPX_TAG_RELEASE "release"
64 /*! \brief Tag string for the file format written to run input files
65 * written by this version of the code.
67 * Change this if you want to change the run input format in a feature
68 * branch. This ensures that there will not be different run input
69 * formats around which cannot be distinguished, while not causing
70 * problems rebasing the feature branch onto upstream changes. When
71 * merging with mainstream GROMACS, set this tag string back to
72 * TPX_TAG_RELEASE, and instead add an element to tpxv and set
73 * tpx_version to that.
75 static const char *tpx_tag = TPX_TAG_RELEASE;
77 /*! \brief Enum of values that describe the contents of a tpr file
78 * whose format matches a version number
80 * The enum helps the code be more self-documenting and ensure merges
81 * do not silently resolve when two patches make the same bump. When
82 * adding new functionality, add a new element to the end of this
83 * enumeration, change the definition of tpx_version, and write code
84 * below that does the right thing according to the value of
87 tpxv_ComputationalElectrophysiology = 96, /**< support for ion/water position swaps (computational electrophysiology) */
88 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to gmx_int64_t */
89 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials /**< potentials for supporting coarse-grained force fields */
92 /*! \brief Version number of the file format written to run input
93 * files by this version of the code.
95 * The tpx_version number should be increased whenever the file format
96 * in the main development branch changes, generally to the highest
97 * value present in tpxv. Backward compatibility for reading old run
98 * input files is maintained by checking this version number against
99 * that of the file and then using the correct code path.
101 * When developing a feature branch that needs to change the run input
102 * file format, change tpx_tag instead. */
103 static const int tpx_version = tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials;
105 /* This number should only be increased when you edit the TOPOLOGY section
106 * or the HEADER of the tpx format.
107 * This way we can maintain forward compatibility too for all analysis tools
108 * and/or external programs that only need to know the atom/residue names,
109 * charges, and bond connectivity.
111 * It first appeared in tpx version 26, when I also moved the inputrecord
112 * to the end of the tpx file, so we can just skip it if we only
115 * In particular, it must be increased when adding new elements to
116 * ftupd, so that old code can read new .tpr files.
118 static const int tpx_generation = 26;
120 /* This number should be the most recent backwards incompatible version
121 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
123 static const int tpx_incompatible_version = 9;
127 /* Struct used to maintain tpx compatibility when function types are added */
129 int fvnr; /* file version number in which the function type first appeared */
130 int ftype; /* function type */
134 * The entries should be ordered in:
135 * 1. ascending file version number
136 * 2. ascending function type number
138 /*static const t_ftupd ftupd[] = {
139 { 20, F_CUBICBONDS },
143 { 22, F_DISRESVIOL },
149 { 26, F_DIHRESVIOL },
150 { 30, F_CROSS_BOND_BONDS },
151 { 30, F_CROSS_BOND_ANGLES },
152 { 30, F_UREY_BRADLEY },
153 { 30, F_POLARIZATION },
157 * The entries should be ordered in:
158 * 1. ascending function type number
159 * 2. ascending file version number
161 /* question; what is the purpose of the commented code above? */
162 static const t_ftupd ftupd[] = {
163 { 20, F_CUBICBONDS },
168 { 43, F_TABBONDSNC },
169 { 70, F_RESTRBONDS },
170 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
171 { 76, F_LINEAR_ANGLES },
172 { 30, F_CROSS_BOND_BONDS },
173 { 30, F_CROSS_BOND_ANGLES },
174 { 30, F_UREY_BRADLEY },
175 { 34, F_QUARTIC_ANGLES },
177 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
178 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
187 { 72, F_NPSOLVATION },
189 { 41, F_LJC_PAIRS_NB },
192 { 32, F_COUL_RECIP },
195 { 30, F_POLARIZATION },
198 { 22, F_DISRESVIOL },
202 { 26, F_DIHRESVIOL },
207 { 46, F_ECONSERVED },
208 { 69, F_VTEMP_NOLONGERUSED},
210 { 54, F_DVDL_CONSTR },
211 { 76, F_ANHARM_POL },
214 { 79, F_DVDL_BONDED, },
215 { 79, F_DVDL_RESTRAINT },
216 { 79, F_DVDL_TEMPERATURE },
218 #define NFTUPD asize(ftupd)
220 /* Needed for backward compatibility */
223 static void _do_section(t_fileio *fio, int key, gmx_bool bRead, const char *src,
229 if (gmx_fio_getftp(fio) == efTPA)
233 gmx_fio_write_string(fio, itemstr[key]);
234 bDbg = gmx_fio_getdebug(fio);
235 gmx_fio_setdebug(fio, FALSE);
236 gmx_fio_write_string(fio, comment_str[key]);
237 gmx_fio_setdebug(fio, bDbg);
241 if (gmx_fio_getdebug(fio))
243 fprintf(stderr, "Looking for section %s (%s, %d)",
244 itemstr[key], src, line);
249 gmx_fio_do_string(fio, buf);
251 while ((gmx_strcasecmp(buf, itemstr[key]) != 0));
253 if (gmx_strcasecmp(buf, itemstr[key]) != 0)
255 gmx_fatal(FARGS, "\nCould not find section heading %s", itemstr[key]);
257 else if (gmx_fio_getdebug(fio))
259 fprintf(stderr, " and found it\n");
265 #define do_section(fio, key, bRead) _do_section(fio, key, bRead, __FILE__, __LINE__)
267 /**************************************************************
269 * Now the higer level routines that do io of the structures and arrays
271 **************************************************************/
272 static void do_pullgrp_tpx_pre95(t_fileio *fio,
281 gmx_fio_do_int(fio, pgrp->nat);
284 snew(pgrp->ind, pgrp->nat);
286 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
287 gmx_fio_do_int(fio, pgrp->nweight);
290 snew(pgrp->weight, pgrp->nweight);
292 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
293 gmx_fio_do_int(fio, pgrp->pbcatom);
294 gmx_fio_do_rvec(fio, pcrd->vec);
295 clear_rvec(pcrd->origin);
296 gmx_fio_do_rvec(fio, tmp);
298 gmx_fio_do_real(fio, pcrd->rate);
299 gmx_fio_do_real(fio, pcrd->k);
300 if (file_version >= 56)
302 gmx_fio_do_real(fio, pcrd->kB);
310 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
314 gmx_fio_do_int(fio, pgrp->nat);
317 snew(pgrp->ind, pgrp->nat);
319 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
320 gmx_fio_do_int(fio, pgrp->nweight);
323 snew(pgrp->weight, pgrp->nweight);
325 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
326 gmx_fio_do_int(fio, pgrp->pbcatom);
329 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd)
333 gmx_fio_do_int(fio, pcrd->group[0]);
334 gmx_fio_do_int(fio, pcrd->group[1]);
335 gmx_fio_do_rvec(fio, pcrd->origin);
336 gmx_fio_do_rvec(fio, pcrd->vec);
337 gmx_fio_do_real(fio, pcrd->init);
338 gmx_fio_do_real(fio, pcrd->rate);
339 gmx_fio_do_real(fio, pcrd->k);
340 gmx_fio_do_real(fio, pcrd->kB);
343 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
345 /* i is used in the ndo_double macro*/
349 int n_lambda = fepvals->n_lambda;
351 /* reset the lambda calculation window */
352 fepvals->lambda_start_n = 0;
353 fepvals->lambda_stop_n = n_lambda;
354 if (file_version >= 79)
360 snew(expand->init_lambda_weights, n_lambda);
362 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
363 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
366 gmx_fio_do_int(fio, expand->nstexpanded);
367 gmx_fio_do_int(fio, expand->elmcmove);
368 gmx_fio_do_int(fio, expand->elamstats);
369 gmx_fio_do_int(fio, expand->lmc_repeats);
370 gmx_fio_do_int(fio, expand->gibbsdeltalam);
371 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
372 gmx_fio_do_int(fio, expand->lmc_seed);
373 gmx_fio_do_real(fio, expand->mc_temp);
374 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
375 gmx_fio_do_int(fio, expand->nstTij);
376 gmx_fio_do_int(fio, expand->minvarmin);
377 gmx_fio_do_int(fio, expand->c_range);
378 gmx_fio_do_real(fio, expand->wl_scale);
379 gmx_fio_do_real(fio, expand->wl_ratio);
380 gmx_fio_do_real(fio, expand->init_wl_delta);
381 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
382 gmx_fio_do_int(fio, expand->elmceq);
383 gmx_fio_do_int(fio, expand->equil_steps);
384 gmx_fio_do_int(fio, expand->equil_samples);
385 gmx_fio_do_int(fio, expand->equil_n_at_lam);
386 gmx_fio_do_real(fio, expand->equil_wl_delta);
387 gmx_fio_do_real(fio, expand->equil_ratio);
391 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
394 if (file_version >= 79)
396 gmx_fio_do_int(fio, simtemp->eSimTempScale);
397 gmx_fio_do_real(fio, simtemp->simtemp_high);
398 gmx_fio_do_real(fio, simtemp->simtemp_low);
403 snew(simtemp->temperatures, n_lambda);
405 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
410 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
412 /* i is defined in the ndo_double macro; use g to iterate. */
417 /* free energy values */
419 if (file_version >= 79)
421 gmx_fio_do_int(fio, fepvals->init_fep_state);
422 gmx_fio_do_double(fio, fepvals->init_lambda);
423 gmx_fio_do_double(fio, fepvals->delta_lambda);
425 else if (file_version >= 59)
427 gmx_fio_do_double(fio, fepvals->init_lambda);
428 gmx_fio_do_double(fio, fepvals->delta_lambda);
432 gmx_fio_do_real(fio, rdum);
433 fepvals->init_lambda = rdum;
434 gmx_fio_do_real(fio, rdum);
435 fepvals->delta_lambda = rdum;
437 if (file_version >= 79)
439 gmx_fio_do_int(fio, fepvals->n_lambda);
442 snew(fepvals->all_lambda, efptNR);
444 for (g = 0; g < efptNR; g++)
446 if (fepvals->n_lambda > 0)
450 snew(fepvals->all_lambda[g], fepvals->n_lambda);
452 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
453 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
455 else if (fepvals->init_lambda >= 0)
457 fepvals->separate_dvdl[efptFEP] = TRUE;
461 else if (file_version >= 64)
463 gmx_fio_do_int(fio, fepvals->n_lambda);
468 snew(fepvals->all_lambda, efptNR);
469 /* still allocate the all_lambda array's contents. */
470 for (g = 0; g < efptNR; g++)
472 if (fepvals->n_lambda > 0)
474 snew(fepvals->all_lambda[g], fepvals->n_lambda);
478 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
480 if (fepvals->init_lambda >= 0)
484 fepvals->separate_dvdl[efptFEP] = TRUE;
488 /* copy the contents of the efptFEP lambda component to all
489 the other components */
490 for (g = 0; g < efptNR; g++)
492 for (h = 0; h < fepvals->n_lambda; h++)
496 fepvals->all_lambda[g][h] =
497 fepvals->all_lambda[efptFEP][h];
506 fepvals->n_lambda = 0;
507 fepvals->all_lambda = NULL;
508 if (fepvals->init_lambda >= 0)
510 fepvals->separate_dvdl[efptFEP] = TRUE;
513 if (file_version >= 13)
515 gmx_fio_do_real(fio, fepvals->sc_alpha);
519 fepvals->sc_alpha = 0;
521 if (file_version >= 38)
523 gmx_fio_do_int(fio, fepvals->sc_power);
527 fepvals->sc_power = 2;
529 if (file_version >= 79)
531 gmx_fio_do_real(fio, fepvals->sc_r_power);
535 fepvals->sc_r_power = 6.0;
537 if (file_version >= 15)
539 gmx_fio_do_real(fio, fepvals->sc_sigma);
543 fepvals->sc_sigma = 0.3;
547 if (file_version >= 71)
549 fepvals->sc_sigma_min = fepvals->sc_sigma;
553 fepvals->sc_sigma_min = 0;
556 if (file_version >= 79)
558 gmx_fio_do_int(fio, fepvals->bScCoul);
562 fepvals->bScCoul = TRUE;
564 if (file_version >= 64)
566 gmx_fio_do_int(fio, fepvals->nstdhdl);
570 fepvals->nstdhdl = 1;
573 if (file_version >= 73)
575 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
576 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
580 fepvals->separate_dhdl_file = esepdhdlfileYES;
581 fepvals->dhdl_derivatives = edhdlderivativesYES;
583 if (file_version >= 71)
585 gmx_fio_do_int(fio, fepvals->dh_hist_size);
586 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
590 fepvals->dh_hist_size = 0;
591 fepvals->dh_hist_spacing = 0.1;
593 if (file_version >= 79)
595 gmx_fio_do_int(fio, fepvals->bPrintEnergy);
599 fepvals->bPrintEnergy = FALSE;
602 /* handle lambda_neighbors */
603 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
605 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
606 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
607 (fepvals->init_lambda < 0) )
609 fepvals->lambda_start_n = (fepvals->init_fep_state -
610 fepvals->lambda_neighbors);
611 fepvals->lambda_stop_n = (fepvals->init_fep_state +
612 fepvals->lambda_neighbors + 1);
613 if (fepvals->lambda_start_n < 0)
615 fepvals->lambda_start_n = 0;;
617 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
619 fepvals->lambda_stop_n = fepvals->n_lambda;
624 fepvals->lambda_start_n = 0;
625 fepvals->lambda_stop_n = fepvals->n_lambda;
630 fepvals->lambda_start_n = 0;
631 fepvals->lambda_stop_n = fepvals->n_lambda;
635 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
639 if (file_version >= 95)
641 gmx_fio_do_int(fio, pull->ngroup);
643 gmx_fio_do_int(fio, pull->ncoord);
644 if (file_version < 95)
646 pull->ngroup = pull->ncoord + 1;
648 gmx_fio_do_int(fio, pull->eGeom);
649 gmx_fio_do_ivec(fio, pull->dim);
650 gmx_fio_do_real(fio, pull->cyl_r1);
651 gmx_fio_do_real(fio, pull->cyl_r0);
652 gmx_fio_do_real(fio, pull->constr_tol);
653 if (file_version >= 95)
655 gmx_fio_do_int(fio, pull->bPrintRef);
657 gmx_fio_do_int(fio, pull->nstxout);
658 gmx_fio_do_int(fio, pull->nstfout);
661 snew(pull->group, pull->ngroup);
662 snew(pull->coord, pull->ncoord);
664 if (file_version < 95)
666 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
667 if (pull->eGeom == epullgDIRPBC)
669 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
671 if (pull->eGeom > epullgDIRPBC)
676 for (g = 0; g < pull->ngroup; g++)
678 /* We read and ignore a pull coordinate for group 0 */
679 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[max(g-1, 0)],
680 bRead, file_version);
683 pull->coord[g-1].group[0] = 0;
684 pull->coord[g-1].group[1] = g;
688 pull->bPrintRef = (pull->group[0].nat > 0);
692 for (g = 0; g < pull->ngroup; g++)
694 do_pull_group(fio, &pull->group[g], bRead);
696 for (g = 0; g < pull->ncoord; g++)
698 do_pull_coord(fio, &pull->coord[g]);
704 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
708 gmx_fio_do_int(fio, rotg->eType);
709 gmx_fio_do_int(fio, rotg->bMassW);
710 gmx_fio_do_int(fio, rotg->nat);
713 snew(rotg->ind, rotg->nat);
715 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
718 snew(rotg->x_ref, rotg->nat);
720 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
721 gmx_fio_do_rvec(fio, rotg->vec);
722 gmx_fio_do_rvec(fio, rotg->pivot);
723 gmx_fio_do_real(fio, rotg->rate);
724 gmx_fio_do_real(fio, rotg->k);
725 gmx_fio_do_real(fio, rotg->slab_dist);
726 gmx_fio_do_real(fio, rotg->min_gaussian);
727 gmx_fio_do_real(fio, rotg->eps);
728 gmx_fio_do_int(fio, rotg->eFittype);
729 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
730 gmx_fio_do_real(fio, rotg->PotAngle_step);
733 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
737 gmx_fio_do_int(fio, rot->ngrp);
738 gmx_fio_do_int(fio, rot->nstrout);
739 gmx_fio_do_int(fio, rot->nstsout);
742 snew(rot->grp, rot->ngrp);
744 for (g = 0; g < rot->ngrp; g++)
746 do_rotgrp(fio, &rot->grp[g], bRead);
751 static void do_swapcoords(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead)
756 gmx_fio_do_int(fio, swap->nat);
757 gmx_fio_do_int(fio, swap->nat_sol);
758 for (j = 0; j < 2; j++)
760 gmx_fio_do_int(fio, swap->nat_split[j]);
761 gmx_fio_do_int(fio, swap->massw_split[j]);
763 gmx_fio_do_int(fio, swap->nstswap);
764 gmx_fio_do_int(fio, swap->nAverage);
765 gmx_fio_do_real(fio, swap->threshold);
766 gmx_fio_do_real(fio, swap->cyl0r);
767 gmx_fio_do_real(fio, swap->cyl0u);
768 gmx_fio_do_real(fio, swap->cyl0l);
769 gmx_fio_do_real(fio, swap->cyl1r);
770 gmx_fio_do_real(fio, swap->cyl1u);
771 gmx_fio_do_real(fio, swap->cyl1l);
775 snew(swap->ind, swap->nat);
776 snew(swap->ind_sol, swap->nat_sol);
777 for (j = 0; j < 2; j++)
779 snew(swap->ind_split[j], swap->nat_split[j]);
783 gmx_fio_ndo_int(fio, swap->ind, swap->nat);
784 gmx_fio_ndo_int(fio, swap->ind_sol, swap->nat_sol);
785 for (j = 0; j < 2; j++)
787 gmx_fio_ndo_int(fio, swap->ind_split[j], swap->nat_split[j]);
790 for (j = 0; j < eCompNR; j++)
792 gmx_fio_do_int(fio, swap->nanions[j]);
793 gmx_fio_do_int(fio, swap->ncations[j]);
799 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
800 int file_version, real *fudgeQQ)
802 int i, j, k, *tmp, idum = 0;
806 real zerotemptime, finish_t, init_temp, finish_temp;
808 if (file_version != tpx_version)
810 /* Give a warning about features that are not accessible */
811 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
812 file_version, tpx_version);
820 if (file_version == 0)
825 /* Basic inputrec stuff */
826 gmx_fio_do_int(fio, ir->eI);
827 if (file_version >= 62)
829 gmx_fio_do_int64(fio, ir->nsteps);
833 gmx_fio_do_int(fio, idum);
836 if (file_version > 25)
838 if (file_version >= 62)
840 gmx_fio_do_int64(fio, ir->init_step);
844 gmx_fio_do_int(fio, idum);
845 ir->init_step = idum;
853 if (file_version >= 58)
855 gmx_fio_do_int(fio, ir->simulation_part);
859 ir->simulation_part = 1;
862 if (file_version >= 67)
864 gmx_fio_do_int(fio, ir->nstcalcenergy);
868 ir->nstcalcenergy = 1;
870 if (file_version < 53)
872 /* The pbc info has been moved out of do_inputrec,
873 * since we always want it, also without reading the inputrec.
875 gmx_fio_do_int(fio, ir->ePBC);
876 if ((file_version <= 15) && (ir->ePBC == 2))
880 if (file_version >= 45)
882 gmx_fio_do_int(fio, ir->bPeriodicMols);
889 ir->bPeriodicMols = TRUE;
893 ir->bPeriodicMols = FALSE;
897 if (file_version >= 81)
899 gmx_fio_do_int(fio, ir->cutoff_scheme);
900 if (file_version < 94)
902 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
907 ir->cutoff_scheme = ecutsGROUP;
909 gmx_fio_do_int(fio, ir->ns_type);
910 gmx_fio_do_int(fio, ir->nstlist);
911 gmx_fio_do_int(fio, ir->ndelta);
912 if (file_version < 41)
914 gmx_fio_do_int(fio, idum);
915 gmx_fio_do_int(fio, idum);
917 if (file_version >= 45)
919 gmx_fio_do_real(fio, ir->rtpi);
925 gmx_fio_do_int(fio, ir->nstcomm);
926 if (file_version > 34)
928 gmx_fio_do_int(fio, ir->comm_mode);
930 else if (ir->nstcomm < 0)
932 ir->comm_mode = ecmANGULAR;
936 ir->comm_mode = ecmLINEAR;
938 ir->nstcomm = abs(ir->nstcomm);
940 if (file_version > 25)
942 gmx_fio_do_int(fio, ir->nstcheckpoint);
946 ir->nstcheckpoint = 0;
949 gmx_fio_do_int(fio, ir->nstcgsteep);
951 if (file_version >= 30)
953 gmx_fio_do_int(fio, ir->nbfgscorr);
960 gmx_fio_do_int(fio, ir->nstlog);
961 gmx_fio_do_int(fio, ir->nstxout);
962 gmx_fio_do_int(fio, ir->nstvout);
963 gmx_fio_do_int(fio, ir->nstfout);
964 gmx_fio_do_int(fio, ir->nstenergy);
965 gmx_fio_do_int(fio, ir->nstxout_compressed);
966 if (file_version >= 59)
968 gmx_fio_do_double(fio, ir->init_t);
969 gmx_fio_do_double(fio, ir->delta_t);
973 gmx_fio_do_real(fio, rdum);
975 gmx_fio_do_real(fio, rdum);
978 gmx_fio_do_real(fio, ir->x_compression_precision);
979 if (file_version < 19)
981 gmx_fio_do_int(fio, idum);
982 gmx_fio_do_int(fio, idum);
984 if (file_version < 18)
986 gmx_fio_do_int(fio, idum);
988 if (file_version >= 81)
990 gmx_fio_do_real(fio, ir->verletbuf_tol);
994 ir->verletbuf_tol = 0;
996 gmx_fio_do_real(fio, ir->rlist);
997 if (file_version >= 67)
999 gmx_fio_do_real(fio, ir->rlistlong);
1001 if (file_version >= 82 && file_version != 90)
1003 gmx_fio_do_int(fio, ir->nstcalclr);
1007 /* Calculate at NS steps */
1008 ir->nstcalclr = ir->nstlist;
1010 gmx_fio_do_int(fio, ir->coulombtype);
1011 if (file_version < 32 && ir->coulombtype == eelRF)
1013 ir->coulombtype = eelRF_NEC;
1015 if (file_version >= 81)
1017 gmx_fio_do_int(fio, ir->coulomb_modifier);
1021 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1023 gmx_fio_do_real(fio, ir->rcoulomb_switch);
1024 gmx_fio_do_real(fio, ir->rcoulomb);
1025 gmx_fio_do_int(fio, ir->vdwtype);
1026 if (file_version >= 81)
1028 gmx_fio_do_int(fio, ir->vdw_modifier);
1032 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1034 gmx_fio_do_real(fio, ir->rvdw_switch);
1035 gmx_fio_do_real(fio, ir->rvdw);
1036 if (file_version < 67)
1038 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1040 gmx_fio_do_int(fio, ir->eDispCorr);
1041 gmx_fio_do_real(fio, ir->epsilon_r);
1042 if (file_version >= 37)
1044 gmx_fio_do_real(fio, ir->epsilon_rf);
1048 if (EEL_RF(ir->coulombtype))
1050 ir->epsilon_rf = ir->epsilon_r;
1051 ir->epsilon_r = 1.0;
1055 ir->epsilon_rf = 1.0;
1058 if (file_version >= 29)
1060 gmx_fio_do_real(fio, ir->tabext);
1067 if (file_version > 25)
1069 gmx_fio_do_int(fio, ir->gb_algorithm);
1070 gmx_fio_do_int(fio, ir->nstgbradii);
1071 gmx_fio_do_real(fio, ir->rgbradii);
1072 gmx_fio_do_real(fio, ir->gb_saltconc);
1073 gmx_fio_do_int(fio, ir->implicit_solvent);
1077 ir->gb_algorithm = egbSTILL;
1080 ir->gb_saltconc = 0;
1081 ir->implicit_solvent = eisNO;
1083 if (file_version >= 55)
1085 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1086 gmx_fio_do_real(fio, ir->gb_obc_alpha);
1087 gmx_fio_do_real(fio, ir->gb_obc_beta);
1088 gmx_fio_do_real(fio, ir->gb_obc_gamma);
1089 if (file_version >= 60)
1091 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1092 gmx_fio_do_int(fio, ir->sa_algorithm);
1096 ir->gb_dielectric_offset = 0.009;
1097 ir->sa_algorithm = esaAPPROX;
1099 gmx_fio_do_real(fio, ir->sa_surface_tension);
1101 /* Override sa_surface_tension if it is not changed in the mpd-file */
1102 if (ir->sa_surface_tension < 0)
1104 if (ir->gb_algorithm == egbSTILL)
1106 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1108 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1110 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1117 /* Better use sensible values than insane (0.0) ones... */
1118 ir->gb_epsilon_solvent = 80;
1119 ir->gb_obc_alpha = 1.0;
1120 ir->gb_obc_beta = 0.8;
1121 ir->gb_obc_gamma = 4.85;
1122 ir->sa_surface_tension = 2.092;
1126 if (file_version >= 81)
1128 gmx_fio_do_real(fio, ir->fourier_spacing);
1132 ir->fourier_spacing = 0.0;
1134 gmx_fio_do_int(fio, ir->nkx);
1135 gmx_fio_do_int(fio, ir->nky);
1136 gmx_fio_do_int(fio, ir->nkz);
1137 gmx_fio_do_int(fio, ir->pme_order);
1138 gmx_fio_do_real(fio, ir->ewald_rtol);
1140 if (file_version >= 93)
1142 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1146 ir->ewald_rtol_lj = ir->ewald_rtol;
1149 if (file_version >= 24)
1151 gmx_fio_do_int(fio, ir->ewald_geometry);
1155 ir->ewald_geometry = eewg3D;
1158 if (file_version <= 17)
1160 ir->epsilon_surface = 0;
1161 if (file_version == 17)
1163 gmx_fio_do_int(fio, idum);
1168 gmx_fio_do_real(fio, ir->epsilon_surface);
1171 gmx_fio_do_gmx_bool(fio, ir->bOptFFT);
1173 if (file_version >= 93)
1175 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1177 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1178 gmx_fio_do_int(fio, ir->etc);
1179 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1180 * but the values 0 and 1 still mean no and
1181 * berendsen temperature coupling, respectively.
1183 if (file_version >= 79)
1185 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1187 if (file_version >= 71)
1189 gmx_fio_do_int(fio, ir->nsttcouple);
1193 ir->nsttcouple = ir->nstcalcenergy;
1195 if (file_version <= 15)
1197 gmx_fio_do_int(fio, idum);
1199 if (file_version <= 17)
1201 gmx_fio_do_int(fio, ir->epct);
1202 if (file_version <= 15)
1206 ir->epct = epctSURFACETENSION;
1208 gmx_fio_do_int(fio, idum);
1211 /* we have removed the NO alternative at the beginning */
1215 ir->epct = epctISOTROPIC;
1219 ir->epc = epcBERENDSEN;
1224 gmx_fio_do_int(fio, ir->epc);
1225 gmx_fio_do_int(fio, ir->epct);
1227 if (file_version >= 71)
1229 gmx_fio_do_int(fio, ir->nstpcouple);
1233 ir->nstpcouple = ir->nstcalcenergy;
1235 gmx_fio_do_real(fio, ir->tau_p);
1236 if (file_version <= 15)
1238 gmx_fio_do_rvec(fio, vdum);
1239 clear_mat(ir->ref_p);
1240 for (i = 0; i < DIM; i++)
1242 ir->ref_p[i][i] = vdum[i];
1247 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1248 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1249 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1251 if (file_version <= 15)
1253 gmx_fio_do_rvec(fio, vdum);
1254 clear_mat(ir->compress);
1255 for (i = 0; i < DIM; i++)
1257 ir->compress[i][i] = vdum[i];
1262 gmx_fio_do_rvec(fio, ir->compress[XX]);
1263 gmx_fio_do_rvec(fio, ir->compress[YY]);
1264 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1266 if (file_version >= 47)
1268 gmx_fio_do_int(fio, ir->refcoord_scaling);
1269 gmx_fio_do_rvec(fio, ir->posres_com);
1270 gmx_fio_do_rvec(fio, ir->posres_comB);
1274 ir->refcoord_scaling = erscNO;
1275 clear_rvec(ir->posres_com);
1276 clear_rvec(ir->posres_comB);
1278 if ((file_version > 25) && (file_version < 79))
1280 gmx_fio_do_int(fio, ir->andersen_seed);
1284 ir->andersen_seed = 0;
1286 if (file_version < 26)
1288 gmx_fio_do_gmx_bool(fio, bSimAnn);
1289 gmx_fio_do_real(fio, zerotemptime);
1292 if (file_version < 37)
1294 gmx_fio_do_real(fio, rdum);
1297 gmx_fio_do_real(fio, ir->shake_tol);
1298 if (file_version < 54)
1300 gmx_fio_do_real(fio, *fudgeQQ);
1303 gmx_fio_do_int(fio, ir->efep);
1304 if (file_version <= 14 && ir->efep != efepNO)
1308 do_fepvals(fio, ir->fepvals, bRead, file_version);
1310 if (file_version >= 79)
1312 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1315 ir->bSimTemp = TRUE;
1320 ir->bSimTemp = FALSE;
1324 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1327 if (file_version >= 79)
1329 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1332 ir->bExpanded = TRUE;
1336 ir->bExpanded = FALSE;
1341 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1343 if (file_version >= 57)
1345 gmx_fio_do_int(fio, ir->eDisre);
1347 gmx_fio_do_int(fio, ir->eDisreWeighting);
1348 if (file_version < 22)
1350 if (ir->eDisreWeighting == 0)
1352 ir->eDisreWeighting = edrwEqual;
1356 ir->eDisreWeighting = edrwConservative;
1359 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1360 gmx_fio_do_real(fio, ir->dr_fc);
1361 gmx_fio_do_real(fio, ir->dr_tau);
1362 gmx_fio_do_int(fio, ir->nstdisreout);
1363 if (file_version >= 22)
1365 gmx_fio_do_real(fio, ir->orires_fc);
1366 gmx_fio_do_real(fio, ir->orires_tau);
1367 gmx_fio_do_int(fio, ir->nstorireout);
1373 ir->nstorireout = 0;
1375 if (file_version >= 26 && file_version < 79)
1377 gmx_fio_do_real(fio, ir->dihre_fc);
1378 if (file_version < 56)
1380 gmx_fio_do_real(fio, rdum);
1381 gmx_fio_do_int(fio, idum);
1389 gmx_fio_do_real(fio, ir->em_stepsize);
1390 gmx_fio_do_real(fio, ir->em_tol);
1391 if (file_version >= 22)
1393 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1397 ir->bShakeSOR = TRUE;
1399 if (file_version >= 11)
1401 gmx_fio_do_int(fio, ir->niter);
1406 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1409 if (file_version >= 21)
1411 gmx_fio_do_real(fio, ir->fc_stepsize);
1415 ir->fc_stepsize = 0;
1417 gmx_fio_do_int(fio, ir->eConstrAlg);
1418 gmx_fio_do_int(fio, ir->nProjOrder);
1419 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1420 if (file_version <= 14)
1422 gmx_fio_do_int(fio, idum);
1424 if (file_version >= 26)
1426 gmx_fio_do_int(fio, ir->nLincsIter);
1431 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1434 if (file_version < 33)
1436 gmx_fio_do_real(fio, bd_temp);
1438 gmx_fio_do_real(fio, ir->bd_fric);
1439 if (file_version >= tpxv_Use64BitRandomSeed)
1441 gmx_fio_do_int64(fio, ir->ld_seed);
1445 gmx_fio_do_int(fio, idum);
1448 if (file_version >= 33)
1450 for (i = 0; i < DIM; i++)
1452 gmx_fio_do_rvec(fio, ir->deform[i]);
1457 for (i = 0; i < DIM; i++)
1459 clear_rvec(ir->deform[i]);
1462 if (file_version >= 14)
1464 gmx_fio_do_real(fio, ir->cos_accel);
1470 gmx_fio_do_int(fio, ir->userint1);
1471 gmx_fio_do_int(fio, ir->userint2);
1472 gmx_fio_do_int(fio, ir->userint3);
1473 gmx_fio_do_int(fio, ir->userint4);
1474 gmx_fio_do_real(fio, ir->userreal1);
1475 gmx_fio_do_real(fio, ir->userreal2);
1476 gmx_fio_do_real(fio, ir->userreal3);
1477 gmx_fio_do_real(fio, ir->userreal4);
1480 if (file_version >= 77)
1482 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1487 snew(ir->adress, 1);
1489 gmx_fio_do_int(fio, ir->adress->type);
1490 gmx_fio_do_real(fio, ir->adress->const_wf);
1491 gmx_fio_do_real(fio, ir->adress->ex_width);
1492 gmx_fio_do_real(fio, ir->adress->hy_width);
1493 gmx_fio_do_int(fio, ir->adress->icor);
1494 gmx_fio_do_int(fio, ir->adress->site);
1495 gmx_fio_do_rvec(fio, ir->adress->refs);
1496 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1497 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1498 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1499 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1503 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1505 if (ir->adress->n_tf_grps > 0)
1507 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1511 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1513 if (ir->adress->n_energy_grps > 0)
1515 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1521 ir->bAdress = FALSE;
1525 if (file_version >= 48)
1527 gmx_fio_do_int(fio, ir->ePull);
1528 if (ir->ePull != epullNO)
1534 do_pull(fio, ir->pull, bRead, file_version);
1539 ir->ePull = epullNO;
1542 /* Enforced rotation */
1543 if (file_version >= 74)
1545 gmx_fio_do_int(fio, ir->bRot);
1546 if (ir->bRot == TRUE)
1552 do_rot(fio, ir->rot, bRead);
1561 gmx_fio_do_int(fio, ir->opts.ngtc);
1562 if (file_version >= 69)
1564 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1568 ir->opts.nhchainlength = 1;
1570 gmx_fio_do_int(fio, ir->opts.ngacc);
1571 gmx_fio_do_int(fio, ir->opts.ngfrz);
1572 gmx_fio_do_int(fio, ir->opts.ngener);
1576 snew(ir->opts.nrdf, ir->opts.ngtc);
1577 snew(ir->opts.ref_t, ir->opts.ngtc);
1578 snew(ir->opts.annealing, ir->opts.ngtc);
1579 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1580 snew(ir->opts.anneal_time, ir->opts.ngtc);
1581 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1582 snew(ir->opts.tau_t, ir->opts.ngtc);
1583 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1584 snew(ir->opts.acc, ir->opts.ngacc);
1585 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1587 if (ir->opts.ngtc > 0)
1589 if (bRead && file_version < 13)
1591 snew(tmp, ir->opts.ngtc);
1592 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1593 for (i = 0; i < ir->opts.ngtc; i++)
1595 ir->opts.nrdf[i] = tmp[i];
1601 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1603 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1604 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1605 if (file_version < 33 && ir->eI == eiBD)
1607 for (i = 0; i < ir->opts.ngtc; i++)
1609 ir->opts.tau_t[i] = bd_temp;
1613 if (ir->opts.ngfrz > 0)
1615 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1617 if (ir->opts.ngacc > 0)
1619 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1621 if (file_version >= 12)
1623 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1624 ir->opts.ngener*ir->opts.ngener);
1627 if (bRead && file_version < 26)
1629 for (i = 0; i < ir->opts.ngtc; i++)
1633 ir->opts.annealing[i] = eannSINGLE;
1634 ir->opts.anneal_npoints[i] = 2;
1635 snew(ir->opts.anneal_time[i], 2);
1636 snew(ir->opts.anneal_temp[i], 2);
1637 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1638 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1639 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1640 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1641 ir->opts.anneal_time[i][0] = ir->init_t;
1642 ir->opts.anneal_time[i][1] = finish_t;
1643 ir->opts.anneal_temp[i][0] = init_temp;
1644 ir->opts.anneal_temp[i][1] = finish_temp;
1648 ir->opts.annealing[i] = eannNO;
1649 ir->opts.anneal_npoints[i] = 0;
1655 /* file version 26 or later */
1656 /* First read the lists with annealing and npoints for each group */
1657 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1658 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1659 for (j = 0; j < (ir->opts.ngtc); j++)
1661 k = ir->opts.anneal_npoints[j];
1664 snew(ir->opts.anneal_time[j], k);
1665 snew(ir->opts.anneal_temp[j], k);
1667 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1668 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1672 if (file_version >= 45)
1674 gmx_fio_do_int(fio, ir->nwall);
1675 gmx_fio_do_int(fio, ir->wall_type);
1676 if (file_version >= 50)
1678 gmx_fio_do_real(fio, ir->wall_r_linpot);
1682 ir->wall_r_linpot = -1;
1684 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1685 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1686 gmx_fio_do_real(fio, ir->wall_density[0]);
1687 gmx_fio_do_real(fio, ir->wall_density[1]);
1688 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1694 ir->wall_atomtype[0] = -1;
1695 ir->wall_atomtype[1] = -1;
1696 ir->wall_density[0] = 0;
1697 ir->wall_density[1] = 0;
1698 ir->wall_ewald_zfac = 3;
1700 /* Cosine stuff for electric fields */
1701 for (j = 0; (j < DIM); j++)
1703 gmx_fio_do_int(fio, ir->ex[j].n);
1704 gmx_fio_do_int(fio, ir->et[j].n);
1707 snew(ir->ex[j].a, ir->ex[j].n);
1708 snew(ir->ex[j].phi, ir->ex[j].n);
1709 snew(ir->et[j].a, ir->et[j].n);
1710 snew(ir->et[j].phi, ir->et[j].n);
1712 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1713 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1714 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1715 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1719 if (file_version >= tpxv_ComputationalElectrophysiology)
1721 gmx_fio_do_int(fio, ir->eSwapCoords);
1722 if (ir->eSwapCoords != eswapNO)
1728 do_swapcoords(fio, ir->swap, bRead);
1733 if (file_version >= 39)
1735 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1736 gmx_fio_do_int(fio, ir->QMMMscheme);
1737 gmx_fio_do_real(fio, ir->scalefactor);
1738 gmx_fio_do_int(fio, ir->opts.ngQM);
1741 snew(ir->opts.QMmethod, ir->opts.ngQM);
1742 snew(ir->opts.QMbasis, ir->opts.ngQM);
1743 snew(ir->opts.QMcharge, ir->opts.ngQM);
1744 snew(ir->opts.QMmult, ir->opts.ngQM);
1745 snew(ir->opts.bSH, ir->opts.ngQM);
1746 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1747 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1748 snew(ir->opts.SAon, ir->opts.ngQM);
1749 snew(ir->opts.SAoff, ir->opts.ngQM);
1750 snew(ir->opts.SAsteps, ir->opts.ngQM);
1751 snew(ir->opts.bOPT, ir->opts.ngQM);
1752 snew(ir->opts.bTS, ir->opts.ngQM);
1754 if (ir->opts.ngQM > 0)
1756 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1757 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1758 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1759 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1760 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1761 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1762 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1763 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1764 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1765 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1766 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1767 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1769 /* end of QMMM stuff */
1774 static void do_harm(t_fileio *fio, t_iparams *iparams)
1776 gmx_fio_do_real(fio, iparams->harmonic.rA);
1777 gmx_fio_do_real(fio, iparams->harmonic.krA);
1778 gmx_fio_do_real(fio, iparams->harmonic.rB);
1779 gmx_fio_do_real(fio, iparams->harmonic.krB);
1782 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1783 gmx_bool bRead, int file_version)
1790 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1800 do_harm(fio, iparams);
1801 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1803 /* Correct incorrect storage of parameters */
1804 iparams->pdihs.phiB = iparams->pdihs.phiA;
1805 iparams->pdihs.cpB = iparams->pdihs.cpA;
1809 gmx_fio_do_real(fio, iparams->harmonic.rA);
1810 gmx_fio_do_real(fio, iparams->harmonic.krA);
1812 case F_LINEAR_ANGLES:
1813 gmx_fio_do_real(fio, iparams->linangle.klinA);
1814 gmx_fio_do_real(fio, iparams->linangle.aA);
1815 gmx_fio_do_real(fio, iparams->linangle.klinB);
1816 gmx_fio_do_real(fio, iparams->linangle.aB);
1819 gmx_fio_do_real(fio, iparams->fene.bm);
1820 gmx_fio_do_real(fio, iparams->fene.kb);
1824 gmx_fio_do_real(fio, iparams->restraint.lowA);
1825 gmx_fio_do_real(fio, iparams->restraint.up1A);
1826 gmx_fio_do_real(fio, iparams->restraint.up2A);
1827 gmx_fio_do_real(fio, iparams->restraint.kA);
1828 gmx_fio_do_real(fio, iparams->restraint.lowB);
1829 gmx_fio_do_real(fio, iparams->restraint.up1B);
1830 gmx_fio_do_real(fio, iparams->restraint.up2B);
1831 gmx_fio_do_real(fio, iparams->restraint.kB);
1837 gmx_fio_do_real(fio, iparams->tab.kA);
1838 gmx_fio_do_int(fio, iparams->tab.table);
1839 gmx_fio_do_real(fio, iparams->tab.kB);
1841 case F_CROSS_BOND_BONDS:
1842 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1843 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1844 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1846 case F_CROSS_BOND_ANGLES:
1847 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1848 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1849 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1850 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1852 case F_UREY_BRADLEY:
1853 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1854 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1855 gmx_fio_do_real(fio, iparams->u_b.r13A);
1856 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1857 if (file_version >= 79)
1859 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1860 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1861 gmx_fio_do_real(fio, iparams->u_b.r13B);
1862 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1866 iparams->u_b.thetaB = iparams->u_b.thetaA;
1867 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1868 iparams->u_b.r13B = iparams->u_b.r13A;
1869 iparams->u_b.kUBB = iparams->u_b.kUBA;
1872 case F_QUARTIC_ANGLES:
1873 gmx_fio_do_real(fio, iparams->qangle.theta);
1874 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1877 gmx_fio_do_real(fio, iparams->bham.a);
1878 gmx_fio_do_real(fio, iparams->bham.b);
1879 gmx_fio_do_real(fio, iparams->bham.c);
1882 gmx_fio_do_real(fio, iparams->morse.b0A);
1883 gmx_fio_do_real(fio, iparams->morse.cbA);
1884 gmx_fio_do_real(fio, iparams->morse.betaA);
1885 if (file_version >= 79)
1887 gmx_fio_do_real(fio, iparams->morse.b0B);
1888 gmx_fio_do_real(fio, iparams->morse.cbB);
1889 gmx_fio_do_real(fio, iparams->morse.betaB);
1893 iparams->morse.b0B = iparams->morse.b0A;
1894 iparams->morse.cbB = iparams->morse.cbA;
1895 iparams->morse.betaB = iparams->morse.betaA;
1899 gmx_fio_do_real(fio, iparams->cubic.b0);
1900 gmx_fio_do_real(fio, iparams->cubic.kb);
1901 gmx_fio_do_real(fio, iparams->cubic.kcub);
1905 case F_POLARIZATION:
1906 gmx_fio_do_real(fio, iparams->polarize.alpha);
1909 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1910 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1911 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1914 if (file_version < 31)
1916 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1918 gmx_fio_do_real(fio, iparams->wpol.al_x);
1919 gmx_fio_do_real(fio, iparams->wpol.al_y);
1920 gmx_fio_do_real(fio, iparams->wpol.al_z);
1921 gmx_fio_do_real(fio, iparams->wpol.rOH);
1922 gmx_fio_do_real(fio, iparams->wpol.rHH);
1923 gmx_fio_do_real(fio, iparams->wpol.rOD);
1926 gmx_fio_do_real(fio, iparams->thole.a);
1927 gmx_fio_do_real(fio, iparams->thole.alpha1);
1928 gmx_fio_do_real(fio, iparams->thole.alpha2);
1929 gmx_fio_do_real(fio, iparams->thole.rfac);
1932 gmx_fio_do_real(fio, iparams->lj.c6);
1933 gmx_fio_do_real(fio, iparams->lj.c12);
1936 gmx_fio_do_real(fio, iparams->lj14.c6A);
1937 gmx_fio_do_real(fio, iparams->lj14.c12A);
1938 gmx_fio_do_real(fio, iparams->lj14.c6B);
1939 gmx_fio_do_real(fio, iparams->lj14.c12B);
1942 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1943 gmx_fio_do_real(fio, iparams->ljc14.qi);
1944 gmx_fio_do_real(fio, iparams->ljc14.qj);
1945 gmx_fio_do_real(fio, iparams->ljc14.c6);
1946 gmx_fio_do_real(fio, iparams->ljc14.c12);
1948 case F_LJC_PAIRS_NB:
1949 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1950 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1951 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1952 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1958 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1959 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1960 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1962 /* Read the incorrectly stored multiplicity */
1963 gmx_fio_do_real(fio, iparams->harmonic.rB);
1964 gmx_fio_do_real(fio, iparams->harmonic.krB);
1965 iparams->pdihs.phiB = iparams->pdihs.phiA;
1966 iparams->pdihs.cpB = iparams->pdihs.cpA;
1970 gmx_fio_do_real(fio, iparams->pdihs.phiB);
1971 gmx_fio_do_real(fio, iparams->pdihs.cpB);
1972 gmx_fio_do_int(fio, iparams->pdihs.mult);
1976 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1977 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1980 gmx_fio_do_int(fio, iparams->disres.label);
1981 gmx_fio_do_int(fio, iparams->disres.type);
1982 gmx_fio_do_real(fio, iparams->disres.low);
1983 gmx_fio_do_real(fio, iparams->disres.up1);
1984 gmx_fio_do_real(fio, iparams->disres.up2);
1985 gmx_fio_do_real(fio, iparams->disres.kfac);
1988 gmx_fio_do_int(fio, iparams->orires.ex);
1989 gmx_fio_do_int(fio, iparams->orires.label);
1990 gmx_fio_do_int(fio, iparams->orires.power);
1991 gmx_fio_do_real(fio, iparams->orires.c);
1992 gmx_fio_do_real(fio, iparams->orires.obs);
1993 gmx_fio_do_real(fio, iparams->orires.kfac);
1996 if (file_version < 82)
1998 gmx_fio_do_int(fio, idum);
1999 gmx_fio_do_int(fio, idum);
2001 gmx_fio_do_real(fio, iparams->dihres.phiA);
2002 gmx_fio_do_real(fio, iparams->dihres.dphiA);
2003 gmx_fio_do_real(fio, iparams->dihres.kfacA);
2004 if (file_version >= 82)
2006 gmx_fio_do_real(fio, iparams->dihres.phiB);
2007 gmx_fio_do_real(fio, iparams->dihres.dphiB);
2008 gmx_fio_do_real(fio, iparams->dihres.kfacB);
2012 iparams->dihres.phiB = iparams->dihres.phiA;
2013 iparams->dihres.dphiB = iparams->dihres.dphiA;
2014 iparams->dihres.kfacB = iparams->dihres.kfacA;
2018 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
2019 gmx_fio_do_rvec(fio, iparams->posres.fcA);
2020 if (bRead && file_version < 27)
2022 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
2023 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
2027 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
2028 gmx_fio_do_rvec(fio, iparams->posres.fcB);
2032 gmx_fio_do_int(fio, iparams->fbposres.geom);
2033 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
2034 gmx_fio_do_real(fio, iparams->fbposres.r);
2035 gmx_fio_do_real(fio, iparams->fbposres.k);
2038 gmx_fio_ndo_real(fio, iparams->cbtdihs.cbtcA, NR_CBTDIHS);
2041 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2042 if (file_version >= 25)
2044 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2048 /* Fourier dihedrals are internally represented
2049 * as Ryckaert-Bellemans since those are faster to compute.
2051 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2052 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2056 gmx_fio_do_real(fio, iparams->constr.dA);
2057 gmx_fio_do_real(fio, iparams->constr.dB);
2060 gmx_fio_do_real(fio, iparams->settle.doh);
2061 gmx_fio_do_real(fio, iparams->settle.dhh);
2064 gmx_fio_do_real(fio, iparams->vsite.a);
2069 gmx_fio_do_real(fio, iparams->vsite.a);
2070 gmx_fio_do_real(fio, iparams->vsite.b);
2075 gmx_fio_do_real(fio, iparams->vsite.a);
2076 gmx_fio_do_real(fio, iparams->vsite.b);
2077 gmx_fio_do_real(fio, iparams->vsite.c);
2080 gmx_fio_do_int(fio, iparams->vsiten.n);
2081 gmx_fio_do_real(fio, iparams->vsiten.a);
2086 /* We got rid of some parameters in version 68 */
2087 if (bRead && file_version < 68)
2089 gmx_fio_do_real(fio, rdum);
2090 gmx_fio_do_real(fio, rdum);
2091 gmx_fio_do_real(fio, rdum);
2092 gmx_fio_do_real(fio, rdum);
2094 gmx_fio_do_real(fio, iparams->gb.sar);
2095 gmx_fio_do_real(fio, iparams->gb.st);
2096 gmx_fio_do_real(fio, iparams->gb.pi);
2097 gmx_fio_do_real(fio, iparams->gb.gbr);
2098 gmx_fio_do_real(fio, iparams->gb.bmlt);
2101 gmx_fio_do_int(fio, iparams->cmap.cmapA);
2102 gmx_fio_do_int(fio, iparams->cmap.cmapB);
2105 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
2106 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
2110 gmx_fio_unset_comment(fio);
2114 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
2121 gmx_fio_set_comment(fio, interaction_function[ftype].name);
2123 if (file_version < 44)
2125 for (i = 0; i < MAXNODES; i++)
2127 gmx_fio_do_int(fio, idum);
2130 gmx_fio_do_int(fio, ilist->nr);
2133 snew(ilist->iatoms, ilist->nr);
2135 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2138 gmx_fio_unset_comment(fio);
2142 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2143 gmx_bool bRead, int file_version)
2148 gmx_fio_do_int(fio, ffparams->atnr);
2149 if (file_version < 57)
2151 gmx_fio_do_int(fio, idum);
2153 gmx_fio_do_int(fio, ffparams->ntypes);
2156 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2157 ffparams->atnr, ffparams->ntypes);
2161 snew(ffparams->functype, ffparams->ntypes);
2162 snew(ffparams->iparams, ffparams->ntypes);
2164 /* Read/write all the function types */
2165 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2168 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2171 if (file_version >= 66)
2173 gmx_fio_do_double(fio, ffparams->reppow);
2177 ffparams->reppow = 12.0;
2180 if (file_version >= 57)
2182 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2185 /* Check whether all these function types are supported by the code.
2186 * In practice the code is backwards compatible, which means that the
2187 * numbering may have to be altered from old numbering to new numbering
2189 for (i = 0; (i < ffparams->ntypes); i++)
2193 /* Loop over file versions */
2194 for (k = 0; (k < NFTUPD); k++)
2196 /* Compare the read file_version to the update table */
2197 if ((file_version < ftupd[k].fvnr) &&
2198 (ffparams->functype[i] >= ftupd[k].ftype))
2200 ffparams->functype[i] += 1;
2203 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2204 i, ffparams->functype[i],
2205 interaction_function[ftupd[k].ftype].longname);
2212 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2216 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2221 static void add_settle_atoms(t_ilist *ilist)
2225 /* Settle used to only store the first atom: add the other two */
2226 srenew(ilist->iatoms, 2*ilist->nr);
2227 for (i = ilist->nr/2-1; i >= 0; i--)
2229 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2230 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2231 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2232 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2234 ilist->nr = 2*ilist->nr;
2237 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2240 int i, j, renum[F_NRE];
2244 for (j = 0; (j < F_NRE); j++)
2249 for (k = 0; k < NFTUPD; k++)
2251 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2260 ilist[j].iatoms = NULL;
2264 do_ilist(fio, &ilist[j], bRead, file_version, j);
2265 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2267 add_settle_atoms(&ilist[j]);
2271 if (bRead && gmx_debug_at)
2272 pr_ilist(debug,0,interaction_function[j].longname,
2273 functype,&ilist[j],TRUE);
2278 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2279 gmx_bool bRead, int file_version)
2281 do_ffparams(fio, ffparams, bRead, file_version);
2283 if (file_version >= 54)
2285 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2288 do_ilists(fio, molt->ilist, bRead, file_version);
2291 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2293 int i, idum, dum_nra, *dum_a;
2295 if (file_version < 44)
2297 for (i = 0; i < MAXNODES; i++)
2299 gmx_fio_do_int(fio, idum);
2302 gmx_fio_do_int(fio, block->nr);
2303 if (file_version < 51)
2305 gmx_fio_do_int(fio, dum_nra);
2309 if ((block->nalloc_index > 0) && (NULL != block->index))
2311 sfree(block->index);
2313 block->nalloc_index = block->nr+1;
2314 snew(block->index, block->nalloc_index);
2316 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2318 if (file_version < 51 && dum_nra > 0)
2320 snew(dum_a, dum_nra);
2321 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2326 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2331 if (file_version < 44)
2333 for (i = 0; i < MAXNODES; i++)
2335 gmx_fio_do_int(fio, idum);
2338 gmx_fio_do_int(fio, block->nr);
2339 gmx_fio_do_int(fio, block->nra);
2342 block->nalloc_index = block->nr+1;
2343 snew(block->index, block->nalloc_index);
2344 block->nalloc_a = block->nra;
2345 snew(block->a, block->nalloc_a);
2347 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2348 gmx_fio_ndo_int(fio, block->a, block->nra);
2351 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2352 int file_version, gmx_groups_t *groups, int atnr)
2356 gmx_fio_do_real(fio, atom->m);
2357 gmx_fio_do_real(fio, atom->q);
2358 gmx_fio_do_real(fio, atom->mB);
2359 gmx_fio_do_real(fio, atom->qB);
2360 gmx_fio_do_ushort(fio, atom->type);
2361 gmx_fio_do_ushort(fio, atom->typeB);
2362 gmx_fio_do_int(fio, atom->ptype);
2363 gmx_fio_do_int(fio, atom->resind);
2364 if (file_version >= 52)
2366 gmx_fio_do_int(fio, atom->atomnumber);
2370 atom->atomnumber = NOTSET;
2372 if (file_version < 23)
2376 else if (file_version < 39)
2385 if (file_version < 57)
2387 unsigned char uchar[egcNR];
2388 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2389 for (i = myngrp; (i < ngrp); i++)
2393 /* Copy the old data format to the groups struct */
2394 for (i = 0; i < ngrp; i++)
2396 groups->grpnr[i][atnr] = uchar[i];
2401 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2406 if (file_version < 23)
2410 else if (file_version < 39)
2419 for (j = 0; (j < ngrp); j++)
2423 gmx_fio_do_int(fio, grps[j].nr);
2426 snew(grps[j].nm_ind, grps[j].nr);
2428 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2433 snew(grps[j].nm_ind, grps[j].nr);
2438 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2444 gmx_fio_do_int(fio, ls);
2445 *nm = get_symtab_handle(symtab, ls);
2449 ls = lookup_symtab(symtab, *nm);
2450 gmx_fio_do_int(fio, ls);
2454 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2459 for (j = 0; (j < nstr); j++)
2461 do_symstr(fio, &(nm[j]), bRead, symtab);
2465 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2466 t_symtab *symtab, int file_version)
2470 for (j = 0; (j < n); j++)
2472 do_symstr(fio, &(ri[j].name), bRead, symtab);
2473 if (file_version >= 63)
2475 gmx_fio_do_int(fio, ri[j].nr);
2476 gmx_fio_do_uchar(fio, ri[j].ic);
2486 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2488 gmx_groups_t *groups)
2492 gmx_fio_do_int(fio, atoms->nr);
2493 gmx_fio_do_int(fio, atoms->nres);
2494 if (file_version < 57)
2496 gmx_fio_do_int(fio, groups->ngrpname);
2497 for (i = 0; i < egcNR; i++)
2499 groups->ngrpnr[i] = atoms->nr;
2500 snew(groups->grpnr[i], groups->ngrpnr[i]);
2505 snew(atoms->atom, atoms->nr);
2506 snew(atoms->atomname, atoms->nr);
2507 snew(atoms->atomtype, atoms->nr);
2508 snew(atoms->atomtypeB, atoms->nr);
2509 snew(atoms->resinfo, atoms->nres);
2510 if (file_version < 57)
2512 snew(groups->grpname, groups->ngrpname);
2514 atoms->pdbinfo = NULL;
2516 for (i = 0; (i < atoms->nr); i++)
2518 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2520 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2521 if (bRead && (file_version <= 20))
2523 for (i = 0; i < atoms->nr; i++)
2525 atoms->atomtype[i] = put_symtab(symtab, "?");
2526 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2531 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2532 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2534 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2536 if (file_version < 57)
2538 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2540 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2544 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2545 gmx_bool bRead, t_symtab *symtab,
2550 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2551 gmx_fio_do_int(fio, groups->ngrpname);
2554 snew(groups->grpname, groups->ngrpname);
2556 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2557 for (g = 0; g < egcNR; g++)
2559 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2560 if (groups->ngrpnr[g] == 0)
2564 groups->grpnr[g] = NULL;
2571 snew(groups->grpnr[g], groups->ngrpnr[g]);
2573 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2578 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2583 if (file_version > 25)
2585 gmx_fio_do_int(fio, atomtypes->nr);
2589 snew(atomtypes->radius, j);
2590 snew(atomtypes->vol, j);
2591 snew(atomtypes->surftens, j);
2592 snew(atomtypes->atomnumber, j);
2593 snew(atomtypes->gb_radius, j);
2594 snew(atomtypes->S_hct, j);
2596 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2597 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2598 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2599 if (file_version >= 40)
2601 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2603 if (file_version >= 60)
2605 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2606 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2611 /* File versions prior to 26 cannot do GBSA,
2612 * so they dont use this structure
2615 atomtypes->radius = NULL;
2616 atomtypes->vol = NULL;
2617 atomtypes->surftens = NULL;
2618 atomtypes->atomnumber = NULL;
2619 atomtypes->gb_radius = NULL;
2620 atomtypes->S_hct = NULL;
2624 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2630 gmx_fio_do_int(fio, symtab->nr);
2634 snew(symtab->symbuf, 1);
2635 symbuf = symtab->symbuf;
2636 symbuf->bufsize = nr;
2637 snew(symbuf->buf, nr);
2638 for (i = 0; (i < nr); i++)
2640 gmx_fio_do_string(fio, buf);
2641 symbuf->buf[i] = strdup(buf);
2646 symbuf = symtab->symbuf;
2647 while (symbuf != NULL)
2649 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2651 gmx_fio_do_string(fio, symbuf->buf[i]);
2654 symbuf = symbuf->next;
2658 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2663 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2665 int i, j, ngrid, gs, nelem;
2667 gmx_fio_do_int(fio, cmap_grid->ngrid);
2668 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2670 ngrid = cmap_grid->ngrid;
2671 gs = cmap_grid->grid_spacing;
2676 snew(cmap_grid->cmapdata, ngrid);
2678 for (i = 0; i < cmap_grid->ngrid; i++)
2680 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2684 for (i = 0; i < cmap_grid->ngrid; i++)
2686 for (j = 0; j < nelem; j++)
2688 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2689 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2690 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2691 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2697 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2699 int m, a, a0, a1, r;
2703 /* We always assign a new chain number, but save the chain id characters
2704 * for larger molecules.
2706 #define CHAIN_MIN_ATOMS 15
2710 for (m = 0; m < mols->nr; m++)
2712 a0 = mols->index[m];
2713 a1 = mols->index[m+1];
2714 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2723 for (a = a0; a < a1; a++)
2725 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2726 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2731 /* Blank out the chain id if there was only one chain */
2734 for (r = 0; r < atoms->nres; r++)
2736 atoms->resinfo[r].chainid = ' ';
2741 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2742 t_symtab *symtab, int file_version,
2743 gmx_groups_t *groups)
2747 if (file_version >= 57)
2749 do_symstr(fio, &(molt->name), bRead, symtab);
2752 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2754 if (bRead && gmx_debug_at)
2756 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2759 if (file_version >= 57)
2761 do_ilists(fio, molt->ilist, bRead, file_version);
2763 do_block(fio, &molt->cgs, bRead, file_version);
2764 if (bRead && gmx_debug_at)
2766 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2770 /* This used to be in the atoms struct */
2771 do_blocka(fio, &molt->excls, bRead, file_version);
2774 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2778 gmx_fio_do_int(fio, molb->type);
2779 gmx_fio_do_int(fio, molb->nmol);
2780 gmx_fio_do_int(fio, molb->natoms_mol);
2781 /* Position restraint coordinates */
2782 gmx_fio_do_int(fio, molb->nposres_xA);
2783 if (molb->nposres_xA > 0)
2787 snew(molb->posres_xA, molb->nposres_xA);
2789 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2791 gmx_fio_do_int(fio, molb->nposres_xB);
2792 if (molb->nposres_xB > 0)
2796 snew(molb->posres_xB, molb->nposres_xB);
2798 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2803 static t_block mtop_mols(gmx_mtop_t *mtop)
2809 for (mb = 0; mb < mtop->nmolblock; mb++)
2811 mols.nr += mtop->molblock[mb].nmol;
2813 mols.nalloc_index = mols.nr + 1;
2814 snew(mols.index, mols.nalloc_index);
2819 for (mb = 0; mb < mtop->nmolblock; mb++)
2821 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2823 a += mtop->molblock[mb].natoms_mol;
2832 static void add_posres_molblock(gmx_mtop_t *mtop)
2837 gmx_molblock_t *molb;
2840 /* posres reference positions are stored in ip->posres (if present) and
2841 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2842 posres.pos0A are identical to fbposres.pos0. */
2843 il = &mtop->moltype[0].ilist[F_POSRES];
2844 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2845 if (il->nr == 0 && ilfb->nr == 0)
2851 for (i = 0; i < il->nr; i += 2)
2853 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2854 am = max(am, il->iatoms[i+1]);
2855 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2856 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2857 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2862 /* This loop is required if we have only flat-bottomed posres:
2864 - bFE == FALSE (no B-state for flat-bottomed posres) */
2867 for (i = 0; i < ilfb->nr; i += 2)
2869 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2870 am = max(am, ilfb->iatoms[i+1]);
2873 /* Make the posres coordinate block end at a molecule end */
2875 while (am >= mtop->mols.index[mol+1])
2879 molb = &mtop->molblock[0];
2880 molb->nposres_xA = mtop->mols.index[mol+1];
2881 snew(molb->posres_xA, molb->nposres_xA);
2884 molb->nposres_xB = molb->nposres_xA;
2885 snew(molb->posres_xB, molb->nposres_xB);
2889 molb->nposres_xB = 0;
2891 for (i = 0; i < il->nr; i += 2)
2893 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2894 a = il->iatoms[i+1];
2895 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2896 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2897 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2900 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2901 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2902 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2907 /* If only flat-bottomed posres are present, take reference pos from them.
2908 Here: bFE == FALSE */
2909 for (i = 0; i < ilfb->nr; i += 2)
2911 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2912 a = ilfb->iatoms[i+1];
2913 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2914 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2915 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2920 static void set_disres_npair(gmx_mtop_t *mtop)
2927 ip = mtop->ffparams.iparams;
2929 for (mt = 0; mt < mtop->nmoltype; mt++)
2931 il = &mtop->moltype[mt].ilist[F_DISRES];
2936 for (i = 0; i < il->nr; i += 3)
2939 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2941 ip[a[i]].disres.npair = npair;
2949 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2959 do_symtab(fio, &(mtop->symtab), bRead);
2962 pr_symtab(debug, 0, "symtab", &mtop->symtab);
2965 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2967 if (file_version >= 57)
2969 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
2971 gmx_fio_do_int(fio, mtop->nmoltype);
2979 snew(mtop->moltype, mtop->nmoltype);
2980 if (file_version < 57)
2982 mtop->moltype[0].name = mtop->name;
2985 for (mt = 0; mt < mtop->nmoltype; mt++)
2987 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
2991 if (file_version >= 57)
2993 gmx_fio_do_int(fio, mtop->nmolblock);
2997 mtop->nmolblock = 1;
3001 snew(mtop->molblock, mtop->nmolblock);
3003 if (file_version >= 57)
3005 for (mb = 0; mb < mtop->nmolblock; mb++)
3007 do_molblock(fio, &mtop->molblock[mb], bRead);
3009 gmx_fio_do_int(fio, mtop->natoms);
3013 mtop->molblock[0].type = 0;
3014 mtop->molblock[0].nmol = 1;
3015 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
3016 mtop->molblock[0].nposres_xA = 0;
3017 mtop->molblock[0].nposres_xB = 0;
3020 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
3023 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
3026 if (file_version < 57)
3028 /* Debug statements are inside do_idef */
3029 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
3030 mtop->natoms = mtop->moltype[0].atoms.nr;
3033 if (file_version >= 65)
3035 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
3039 mtop->ffparams.cmap_grid.ngrid = 0;
3040 mtop->ffparams.cmap_grid.grid_spacing = 0;
3041 mtop->ffparams.cmap_grid.cmapdata = NULL;
3044 if (file_version >= 57)
3046 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
3049 if (file_version < 57)
3051 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
3052 if (bRead && gmx_debug_at)
3054 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
3056 do_block(fio, &mtop->mols, bRead, file_version);
3057 /* Add the posres coordinates to the molblock */
3058 add_posres_molblock(mtop);
3062 if (file_version >= 57)
3064 done_block(&mtop->mols);
3065 mtop->mols = mtop_mols(mtop);
3069 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
3073 if (file_version < 51)
3075 /* Here used to be the shake blocks */
3076 do_blocka(fio, &dumb, bRead, file_version);
3089 close_symtab(&(mtop->symtab));
3093 /* If TopOnlyOK is TRUE then we can read even future versions
3094 * of tpx files, provided the file_generation hasn't changed.
3095 * If it is FALSE, we need the inputrecord too, and bail out
3096 * if the file is newer than the program.
3098 * The version and generation if the topology (see top of this file)
3099 * are returned in the two last arguments.
3101 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
3103 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
3104 gmx_bool TopOnlyOK, int *file_version,
3105 int *file_generation)
3108 char file_tag[STRLEN];
3115 gmx_fio_checktype(fio);
3116 gmx_fio_setdebug(fio, bDebugMode());
3118 /* NEW! XDR tpb file */
3119 precision = sizeof(real);
3122 gmx_fio_do_string(fio, buf);
3123 if (strncmp(buf, "VERSION", 7))
3125 gmx_fatal(FARGS, "Can not read file %s,\n"
3126 " this file is from a Gromacs version which is older than 2.0\n"
3127 " Make a new one with grompp or use a gro or pdb file, if possible",
3128 gmx_fio_getname(fio));
3130 gmx_fio_do_int(fio, precision);
3131 bDouble = (precision == sizeof(double));
3132 if ((precision != sizeof(float)) && !bDouble)
3134 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3135 "instead of %d or %d",
3136 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3138 gmx_fio_setprecision(fio, bDouble);
3139 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3140 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3144 gmx_fio_write_string(fio, GromacsVersion());
3145 bDouble = (precision == sizeof(double));
3146 gmx_fio_setprecision(fio, bDouble);
3147 gmx_fio_do_int(fio, precision);
3149 sprintf(file_tag, "%s", tpx_tag);
3150 fgen = tpx_generation;
3153 /* Check versions! */
3154 gmx_fio_do_int(fio, fver);
3156 /* This is for backward compatibility with development versions 77-79
3157 * where the tag was, mistakenly, placed before the generation,
3158 * which would cause a segv instead of a proper error message
3159 * when reading the topology only from tpx with <77 code.
3161 if (fver >= 77 && fver <= 79)
3163 gmx_fio_do_string(fio, file_tag);
3168 gmx_fio_do_int(fio, fgen);
3177 gmx_fio_do_string(fio, file_tag);
3183 /* Versions before 77 don't have the tag, set it to release */
3184 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3187 if (strcmp(file_tag, tpx_tag) != 0)
3189 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3192 /* We only support reading tpx files with the same tag as the code
3193 * or tpx files with the release tag and with lower version number.
3195 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
3197 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3198 gmx_fio_getname(fio), fver, file_tag,
3199 tpx_version, tpx_tag);
3204 if (file_version != NULL)
3206 *file_version = fver;
3208 if (file_generation != NULL)
3210 *file_generation = fgen;
3214 if ((fver <= tpx_incompatible_version) ||
3215 ((fver > tpx_version) && !TopOnlyOK) ||
3216 (fgen > tpx_generation) ||
3217 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3219 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3220 gmx_fio_getname(fio), fver, tpx_version);
3223 do_section(fio, eitemHEADER, bRead);
3224 gmx_fio_do_int(fio, tpx->natoms);
3227 gmx_fio_do_int(fio, tpx->ngtc);
3235 gmx_fio_do_int(fio, idum);
3236 gmx_fio_do_real(fio, rdum);
3238 /*a better decision will eventually (5.0 or later) need to be made
3239 on how to treat the alchemical state of the system, which can now
3240 vary through a simulation, and cannot be completely described
3241 though a single lambda variable, or even a single state
3242 index. Eventually, should probably be a vector. MRS*/
3245 gmx_fio_do_int(fio, tpx->fep_state);
3247 gmx_fio_do_real(fio, tpx->lambda);
3248 gmx_fio_do_int(fio, tpx->bIr);
3249 gmx_fio_do_int(fio, tpx->bTop);
3250 gmx_fio_do_int(fio, tpx->bX);
3251 gmx_fio_do_int(fio, tpx->bV);
3252 gmx_fio_do_int(fio, tpx->bF);
3253 gmx_fio_do_int(fio, tpx->bBox);
3255 if ((fgen > tpx_generation))
3257 /* This can only happen if TopOnlyOK=TRUE */
3262 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3263 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3264 gmx_bool bXVallocated)
3270 int file_version, file_generation;
3274 gmx_bool bPeriodicMols;
3278 tpx.natoms = state->natoms;
3279 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3280 tpx.fep_state = state->fep_state;
3281 tpx.lambda = state->lambda[efptFEP];
3282 tpx.bIr = (ir != NULL);
3283 tpx.bTop = (mtop != NULL);
3284 tpx.bX = (state->x != NULL);
3285 tpx.bV = (state->v != NULL);
3286 tpx.bF = (f != NULL);
3290 TopOnlyOK = (ir == NULL);
3292 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3297 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3298 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3303 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3304 state->natoms = tpx.natoms;
3305 state->nalloc = tpx.natoms;
3311 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3315 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3317 do_test(fio, tpx.bBox, state->box);
3318 do_section(fio, eitemBOX, bRead);
3321 gmx_fio_ndo_rvec(fio, state->box, DIM);
3322 if (file_version >= 51)
3324 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3328 /* We initialize box_rel after reading the inputrec */
3329 clear_mat(state->box_rel);
3331 if (file_version >= 28)
3333 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3334 if (file_version < 56)
3337 gmx_fio_ndo_rvec(fio, mdum, DIM);
3342 if (state->ngtc > 0 && file_version >= 28)
3345 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3346 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3347 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3348 snew(dumv, state->ngtc);
3349 if (file_version < 69)
3351 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3353 /* These used to be the Berendsen tcoupl_lambda's */
3354 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3358 /* Prior to tpx version 26, the inputrec was here.
3359 * I moved it to enable partial forward-compatibility
3360 * for analysis/viewer programs.
3362 if (file_version < 26)
3364 do_test(fio, tpx.bIr, ir);
3365 do_section(fio, eitemIR, bRead);
3370 do_inputrec(fio, ir, bRead, file_version,
3371 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3374 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3379 do_inputrec(fio, &dum_ir, bRead, file_version,
3380 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3383 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3385 done_inputrec(&dum_ir);
3391 do_test(fio, tpx.bTop, mtop);
3392 do_section(fio, eitemTOP, bRead);
3397 do_mtop(fio, mtop, bRead, file_version);
3401 do_mtop(fio, &dum_top, bRead, file_version);
3402 done_mtop(&dum_top, TRUE);
3405 do_test(fio, tpx.bX, state->x);
3406 do_section(fio, eitemX, bRead);
3411 state->flags |= (1<<estX);
3413 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3416 do_test(fio, tpx.bV, state->v);
3417 do_section(fio, eitemV, bRead);
3422 state->flags |= (1<<estV);
3424 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3427 do_test(fio, tpx.bF, f);
3428 do_section(fio, eitemF, bRead);
3431 gmx_fio_ndo_rvec(fio, f, state->natoms);
3434 /* Starting with tpx version 26, we have the inputrec
3435 * at the end of the file, so we can ignore it
3436 * if the file is never than the software (but still the
3437 * same generation - see comments at the top of this file.
3442 bPeriodicMols = FALSE;
3443 if (file_version >= 26)
3445 do_test(fio, tpx.bIr, ir);
3446 do_section(fio, eitemIR, bRead);
3449 if (file_version >= 53)
3451 /* Removed the pbc info from do_inputrec, since we always want it */
3455 bPeriodicMols = ir->bPeriodicMols;
3457 gmx_fio_do_int(fio, ePBC);
3458 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3460 if (file_generation <= tpx_generation && ir)
3462 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3465 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3467 if (file_version < 51)
3469 set_box_rel(ir, state);
3471 if (file_version < 53)
3474 bPeriodicMols = ir->bPeriodicMols;
3477 if (bRead && ir && file_version >= 53)
3479 /* We need to do this after do_inputrec, since that initializes ir */
3481 ir->bPeriodicMols = bPeriodicMols;
3490 if (state->ngtc == 0)
3492 /* Reading old version without tcoupl state data: set it */
3493 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3495 if (tpx.bTop && mtop)
3497 if (file_version < 57)
3499 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3501 ir->eDisre = edrSimple;
3505 ir->eDisre = edrNone;
3508 set_disres_npair(mtop);
3512 if (tpx.bTop && mtop)
3514 gmx_mtop_finalize(mtop);
3517 if (file_version >= 57)
3521 env = getenv("GMX_NOCHARGEGROUPS");
3524 sscanf(env, "%d", &ienv);
3525 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3530 "Will make single atomic charge groups in non-solvent%s\n",
3531 ienv > 1 ? " and solvent" : "");
3532 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3534 fprintf(stderr, "\n");
3542 /************************************************************
3544 * The following routines are the exported ones
3546 ************************************************************/
3548 t_fileio *open_tpx(const char *fn, const char *mode)
3550 return gmx_fio_open(fn, mode);
3553 void close_tpx(t_fileio *fio)
3558 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3559 int *file_version, int *file_generation)
3563 fio = open_tpx(fn, "r");
3564 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3568 void write_tpx_state(const char *fn,
3569 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3573 fio = open_tpx(fn, "w");
3574 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3578 void read_tpx_state(const char *fn,
3579 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3583 fio = open_tpx(fn, "r");
3584 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3588 int read_tpx(const char *fn,
3589 t_inputrec *ir, matrix box, int *natoms,
3590 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3598 fio = open_tpx(fn, "r");
3599 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3601 *natoms = state.natoms;
3604 copy_mat(state.box, box);
3613 int read_tpx_top(const char *fn,
3614 t_inputrec *ir, matrix box, int *natoms,
3615 rvec *x, rvec *v, rvec *f, t_topology *top)
3621 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3623 *top = gmx_mtop_t_to_t_topology(&mtop);
3628 gmx_bool fn2bTPX(const char *file)
3630 switch (fn2ftp(file))
3641 static void done_gmx_groups_t(gmx_groups_t *g)
3645 for (i = 0; (i < egcNR); i++)
3647 if (NULL != g->grps[i].nm_ind)
3649 sfree(g->grps[i].nm_ind);
3650 g->grps[i].nm_ind = NULL;
3652 if (NULL != g->grpnr[i])
3658 /* The contents of this array is in symtab, don't free it here */
3662 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3663 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3666 int natoms, i, version, generation;
3667 gmx_bool bTop, bXNULL = FALSE;
3669 t_topology *topconv;
3672 bTop = fn2bTPX(infile);
3676 read_tpxheader(infile, &header, TRUE, &version, &generation);
3679 snew(*x, header.natoms);
3683 snew(*v, header.natoms);
3686 *ePBC = read_tpx(infile, NULL, box, &natoms,
3687 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3688 *top = gmx_mtop_t_to_t_topology(mtop);
3689 /* In this case we need to throw away the group data too */
3690 done_gmx_groups_t(&mtop->groups);
3692 strcpy(title, *top->name);
3693 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3697 get_stx_coordnum(infile, &natoms);
3698 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3709 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3717 aps = gmx_atomprop_init();
3718 for (i = 0; (i < natoms); i++)
3720 if (!gmx_atomprop_query(aps, epropMass,
3721 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3722 *top->atoms.atomname[i],
3723 &(top->atoms.atom[i].m)))
3727 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3728 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3729 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3730 *top->atoms.atomname[i]);
3734 gmx_atomprop_destroy(aps);
3736 top->idef.ntypes = -1;