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! */
48 #include "gromacs/utility/futil.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/topology/atomprop.h"
58 #include "gromacs/topology/block.h"
59 #include "gromacs/topology/mtop_util.h"
60 #include "gromacs/topology/symtab.h"
61 #include "gromacs/topology/topology.h"
62 #include "gromacs/utility/cstringutil.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/smalloc.h"
66 #define TPX_TAG_RELEASE "release"
68 /*! \brief Tag string for the file format written to run input files
69 * written by this version of the code.
71 * Change this if you want to change the run input format in a feature
72 * branch. This ensures that there will not be different run input
73 * formats around which cannot be distinguished, while not causing
74 * problems rebasing the feature branch onto upstream changes. When
75 * merging with mainstream GROMACS, set this tag string back to
76 * TPX_TAG_RELEASE, and instead add an element to tpxv and set
77 * tpx_version to that.
79 static const char *tpx_tag = TPX_TAG_RELEASE;
81 /*! \brief Enum of values that describe the contents of a tpr file
82 * whose format matches a version number
84 * The enum helps the code be more self-documenting and ensure merges
85 * do not silently resolve when two patches make the same bump. When
86 * adding new functionality, add a new element to the end of this
87 * enumeration, change the definition of tpx_version, and write code
88 * below that does the right thing according to the value of
91 tpxv_ComputationalElectrophysiology = 96, /**< support for ion/water position swaps (computational electrophysiology) */
92 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to gmx_int64_t */
93 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
94 tpxv_InteractiveMolecularDynamics, /**< interactive molecular dynamics (IMD) */
95 tpxv_RemoveObsoleteParameters1 /**< remove optimize_fft, dihre_fc, nstcheckpoint */
98 /*! \brief Version number of the file format written to run input
99 * files by this version of the code.
101 * The tpx_version number should be increased whenever the file format
102 * in the main development branch changes, generally to the highest
103 * value present in tpxv. Backward compatibility for reading old run
104 * input files is maintained by checking this version number against
105 * that of the file and then using the correct code path.
107 * When developing a feature branch that needs to change the run input
108 * file format, change tpx_tag instead. */
109 static const int tpx_version = tpxv_RemoveObsoleteParameters1;
112 /* This number should only be increased when you edit the TOPOLOGY section
113 * or the HEADER of the tpx format.
114 * This way we can maintain forward compatibility too for all analysis tools
115 * and/or external programs that only need to know the atom/residue names,
116 * charges, and bond connectivity.
118 * It first appeared in tpx version 26, when I also moved the inputrecord
119 * to the end of the tpx file, so we can just skip it if we only
122 * In particular, it must be increased when adding new elements to
123 * ftupd, so that old code can read new .tpr files.
125 static const int tpx_generation = 26;
127 /* This number should be the most recent backwards incompatible version
128 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
130 static const int tpx_incompatible_version = 9;
134 /* Struct used to maintain tpx compatibility when function types are added */
136 int fvnr; /* file version number in which the function type first appeared */
137 int ftype; /* function type */
141 * The entries should be ordered in:
142 * 1. ascending file version number
143 * 2. ascending function type number
145 /*static const t_ftupd ftupd[] = {
146 { 20, F_CUBICBONDS },
150 { 22, F_DISRESVIOL },
156 { 26, F_DIHRESVIOL },
157 { 30, F_CROSS_BOND_BONDS },
158 { 30, F_CROSS_BOND_ANGLES },
159 { 30, F_UREY_BRADLEY },
160 { 30, F_POLARIZATION },
164 * The entries should be ordered in:
165 * 1. ascending function type number
166 * 2. ascending file version number
168 /* question; what is the purpose of the commented code above? */
169 static const t_ftupd ftupd[] = {
170 { 20, F_CUBICBONDS },
175 { 43, F_TABBONDSNC },
176 { 70, F_RESTRBONDS },
177 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
178 { 76, F_LINEAR_ANGLES },
179 { 30, F_CROSS_BOND_BONDS },
180 { 30, F_CROSS_BOND_ANGLES },
181 { 30, F_UREY_BRADLEY },
182 { 34, F_QUARTIC_ANGLES },
184 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
185 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
194 { 72, F_NPSOLVATION },
196 { 41, F_LJC_PAIRS_NB },
199 { 32, F_COUL_RECIP },
202 { 30, F_POLARIZATION },
205 { 22, F_DISRESVIOL },
209 { 26, F_DIHRESVIOL },
214 { 46, F_ECONSERVED },
215 { 69, F_VTEMP_NOLONGERUSED},
217 { 54, F_DVDL_CONSTR },
218 { 76, F_ANHARM_POL },
221 { 79, F_DVDL_BONDED, },
222 { 79, F_DVDL_RESTRAINT },
223 { 79, F_DVDL_TEMPERATURE },
225 #define NFTUPD asize(ftupd)
227 /* Needed for backward compatibility */
230 static void _do_section(t_fileio *fio, int key, gmx_bool bRead, const char *src,
236 if (gmx_fio_getftp(fio) == efTPA)
240 gmx_fio_write_string(fio, itemstr[key]);
241 bDbg = gmx_fio_getdebug(fio);
242 gmx_fio_setdebug(fio, FALSE);
243 gmx_fio_write_string(fio, comment_str[key]);
244 gmx_fio_setdebug(fio, bDbg);
248 if (gmx_fio_getdebug(fio))
250 fprintf(stderr, "Looking for section %s (%s, %d)",
251 itemstr[key], src, line);
256 gmx_fio_do_string(fio, buf);
258 while ((gmx_strcasecmp(buf, itemstr[key]) != 0));
260 if (gmx_strcasecmp(buf, itemstr[key]) != 0)
262 gmx_fatal(FARGS, "\nCould not find section heading %s", itemstr[key]);
264 else if (gmx_fio_getdebug(fio))
266 fprintf(stderr, " and found it\n");
272 #define do_section(fio, key, bRead) _do_section(fio, key, bRead, __FILE__, __LINE__)
274 /**************************************************************
276 * Now the higer level routines that do io of the structures and arrays
278 **************************************************************/
279 static void do_pullgrp_tpx_pre95(t_fileio *fio,
288 gmx_fio_do_int(fio, pgrp->nat);
291 snew(pgrp->ind, pgrp->nat);
293 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
294 gmx_fio_do_int(fio, pgrp->nweight);
297 snew(pgrp->weight, pgrp->nweight);
299 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
300 gmx_fio_do_int(fio, pgrp->pbcatom);
301 gmx_fio_do_rvec(fio, pcrd->vec);
302 clear_rvec(pcrd->origin);
303 gmx_fio_do_rvec(fio, tmp);
305 gmx_fio_do_real(fio, pcrd->rate);
306 gmx_fio_do_real(fio, pcrd->k);
307 if (file_version >= 56)
309 gmx_fio_do_real(fio, pcrd->kB);
317 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
321 gmx_fio_do_int(fio, pgrp->nat);
324 snew(pgrp->ind, pgrp->nat);
326 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
327 gmx_fio_do_int(fio, pgrp->nweight);
330 snew(pgrp->weight, pgrp->nweight);
332 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
333 gmx_fio_do_int(fio, pgrp->pbcatom);
336 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd)
340 gmx_fio_do_int(fio, pcrd->group[0]);
341 gmx_fio_do_int(fio, pcrd->group[1]);
342 gmx_fio_do_rvec(fio, pcrd->origin);
343 gmx_fio_do_rvec(fio, pcrd->vec);
344 gmx_fio_do_real(fio, pcrd->init);
345 gmx_fio_do_real(fio, pcrd->rate);
346 gmx_fio_do_real(fio, pcrd->k);
347 gmx_fio_do_real(fio, pcrd->kB);
350 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
352 /* i is used in the ndo_double macro*/
356 int n_lambda = fepvals->n_lambda;
358 /* reset the lambda calculation window */
359 fepvals->lambda_start_n = 0;
360 fepvals->lambda_stop_n = n_lambda;
361 if (file_version >= 79)
367 snew(expand->init_lambda_weights, n_lambda);
369 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
370 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
373 gmx_fio_do_int(fio, expand->nstexpanded);
374 gmx_fio_do_int(fio, expand->elmcmove);
375 gmx_fio_do_int(fio, expand->elamstats);
376 gmx_fio_do_int(fio, expand->lmc_repeats);
377 gmx_fio_do_int(fio, expand->gibbsdeltalam);
378 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
379 gmx_fio_do_int(fio, expand->lmc_seed);
380 gmx_fio_do_real(fio, expand->mc_temp);
381 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
382 gmx_fio_do_int(fio, expand->nstTij);
383 gmx_fio_do_int(fio, expand->minvarmin);
384 gmx_fio_do_int(fio, expand->c_range);
385 gmx_fio_do_real(fio, expand->wl_scale);
386 gmx_fio_do_real(fio, expand->wl_ratio);
387 gmx_fio_do_real(fio, expand->init_wl_delta);
388 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
389 gmx_fio_do_int(fio, expand->elmceq);
390 gmx_fio_do_int(fio, expand->equil_steps);
391 gmx_fio_do_int(fio, expand->equil_samples);
392 gmx_fio_do_int(fio, expand->equil_n_at_lam);
393 gmx_fio_do_real(fio, expand->equil_wl_delta);
394 gmx_fio_do_real(fio, expand->equil_ratio);
398 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
401 if (file_version >= 79)
403 gmx_fio_do_int(fio, simtemp->eSimTempScale);
404 gmx_fio_do_real(fio, simtemp->simtemp_high);
405 gmx_fio_do_real(fio, simtemp->simtemp_low);
410 snew(simtemp->temperatures, n_lambda);
412 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
417 static void do_imd(t_fileio *fio, t_IMD *imd, gmx_bool bRead)
419 gmx_fio_do_int(fio, imd->nat);
422 snew(imd->ind, imd->nat);
424 gmx_fio_ndo_int(fio, imd->ind, imd->nat);
427 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
429 /* i is defined in the ndo_double macro; use g to iterate. */
434 /* free energy values */
436 if (file_version >= 79)
438 gmx_fio_do_int(fio, fepvals->init_fep_state);
439 gmx_fio_do_double(fio, fepvals->init_lambda);
440 gmx_fio_do_double(fio, fepvals->delta_lambda);
442 else if (file_version >= 59)
444 gmx_fio_do_double(fio, fepvals->init_lambda);
445 gmx_fio_do_double(fio, fepvals->delta_lambda);
449 gmx_fio_do_real(fio, rdum);
450 fepvals->init_lambda = rdum;
451 gmx_fio_do_real(fio, rdum);
452 fepvals->delta_lambda = rdum;
454 if (file_version >= 79)
456 gmx_fio_do_int(fio, fepvals->n_lambda);
459 snew(fepvals->all_lambda, efptNR);
461 for (g = 0; g < efptNR; g++)
463 if (fepvals->n_lambda > 0)
467 snew(fepvals->all_lambda[g], fepvals->n_lambda);
469 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
470 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
472 else if (fepvals->init_lambda >= 0)
474 fepvals->separate_dvdl[efptFEP] = TRUE;
478 else if (file_version >= 64)
480 gmx_fio_do_int(fio, fepvals->n_lambda);
485 snew(fepvals->all_lambda, efptNR);
486 /* still allocate the all_lambda array's contents. */
487 for (g = 0; g < efptNR; g++)
489 if (fepvals->n_lambda > 0)
491 snew(fepvals->all_lambda[g], fepvals->n_lambda);
495 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
497 if (fepvals->init_lambda >= 0)
501 fepvals->separate_dvdl[efptFEP] = TRUE;
505 /* copy the contents of the efptFEP lambda component to all
506 the other components */
507 for (g = 0; g < efptNR; g++)
509 for (h = 0; h < fepvals->n_lambda; h++)
513 fepvals->all_lambda[g][h] =
514 fepvals->all_lambda[efptFEP][h];
523 fepvals->n_lambda = 0;
524 fepvals->all_lambda = NULL;
525 if (fepvals->init_lambda >= 0)
527 fepvals->separate_dvdl[efptFEP] = TRUE;
530 if (file_version >= 13)
532 gmx_fio_do_real(fio, fepvals->sc_alpha);
536 fepvals->sc_alpha = 0;
538 if (file_version >= 38)
540 gmx_fio_do_int(fio, fepvals->sc_power);
544 fepvals->sc_power = 2;
546 if (file_version >= 79)
548 gmx_fio_do_real(fio, fepvals->sc_r_power);
552 fepvals->sc_r_power = 6.0;
554 if (file_version >= 15)
556 gmx_fio_do_real(fio, fepvals->sc_sigma);
560 fepvals->sc_sigma = 0.3;
564 if (file_version >= 71)
566 fepvals->sc_sigma_min = fepvals->sc_sigma;
570 fepvals->sc_sigma_min = 0;
573 if (file_version >= 79)
575 gmx_fio_do_int(fio, fepvals->bScCoul);
579 fepvals->bScCoul = TRUE;
581 if (file_version >= 64)
583 gmx_fio_do_int(fio, fepvals->nstdhdl);
587 fepvals->nstdhdl = 1;
590 if (file_version >= 73)
592 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
593 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
597 fepvals->separate_dhdl_file = esepdhdlfileYES;
598 fepvals->dhdl_derivatives = edhdlderivativesYES;
600 if (file_version >= 71)
602 gmx_fio_do_int(fio, fepvals->dh_hist_size);
603 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
607 fepvals->dh_hist_size = 0;
608 fepvals->dh_hist_spacing = 0.1;
610 if (file_version >= 79)
612 gmx_fio_do_int(fio, fepvals->bPrintEnergy);
616 fepvals->bPrintEnergy = FALSE;
619 /* handle lambda_neighbors */
620 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
622 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
623 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
624 (fepvals->init_lambda < 0) )
626 fepvals->lambda_start_n = (fepvals->init_fep_state -
627 fepvals->lambda_neighbors);
628 fepvals->lambda_stop_n = (fepvals->init_fep_state +
629 fepvals->lambda_neighbors + 1);
630 if (fepvals->lambda_start_n < 0)
632 fepvals->lambda_start_n = 0;;
634 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
636 fepvals->lambda_stop_n = fepvals->n_lambda;
641 fepvals->lambda_start_n = 0;
642 fepvals->lambda_stop_n = fepvals->n_lambda;
647 fepvals->lambda_start_n = 0;
648 fepvals->lambda_stop_n = fepvals->n_lambda;
652 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
656 if (file_version >= 95)
658 gmx_fio_do_int(fio, pull->ngroup);
660 gmx_fio_do_int(fio, pull->ncoord);
661 if (file_version < 95)
663 pull->ngroup = pull->ncoord + 1;
665 gmx_fio_do_int(fio, pull->eGeom);
666 gmx_fio_do_ivec(fio, pull->dim);
667 gmx_fio_do_real(fio, pull->cyl_r1);
668 gmx_fio_do_real(fio, pull->cyl_r0);
669 gmx_fio_do_real(fio, pull->constr_tol);
670 if (file_version >= 95)
672 gmx_fio_do_int(fio, pull->bPrintRef);
674 gmx_fio_do_int(fio, pull->nstxout);
675 gmx_fio_do_int(fio, pull->nstfout);
678 snew(pull->group, pull->ngroup);
679 snew(pull->coord, pull->ncoord);
681 if (file_version < 95)
683 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
684 if (pull->eGeom == epullgDIRPBC)
686 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
688 if (pull->eGeom > epullgDIRPBC)
693 for (g = 0; g < pull->ngroup; g++)
695 /* We read and ignore a pull coordinate for group 0 */
696 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[max(g-1, 0)],
697 bRead, file_version);
700 pull->coord[g-1].group[0] = 0;
701 pull->coord[g-1].group[1] = g;
705 pull->bPrintRef = (pull->group[0].nat > 0);
709 for (g = 0; g < pull->ngroup; g++)
711 do_pull_group(fio, &pull->group[g], bRead);
713 for (g = 0; g < pull->ncoord; g++)
715 do_pull_coord(fio, &pull->coord[g]);
721 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
725 gmx_fio_do_int(fio, rotg->eType);
726 gmx_fio_do_int(fio, rotg->bMassW);
727 gmx_fio_do_int(fio, rotg->nat);
730 snew(rotg->ind, rotg->nat);
732 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
735 snew(rotg->x_ref, rotg->nat);
737 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
738 gmx_fio_do_rvec(fio, rotg->vec);
739 gmx_fio_do_rvec(fio, rotg->pivot);
740 gmx_fio_do_real(fio, rotg->rate);
741 gmx_fio_do_real(fio, rotg->k);
742 gmx_fio_do_real(fio, rotg->slab_dist);
743 gmx_fio_do_real(fio, rotg->min_gaussian);
744 gmx_fio_do_real(fio, rotg->eps);
745 gmx_fio_do_int(fio, rotg->eFittype);
746 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
747 gmx_fio_do_real(fio, rotg->PotAngle_step);
750 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
754 gmx_fio_do_int(fio, rot->ngrp);
755 gmx_fio_do_int(fio, rot->nstrout);
756 gmx_fio_do_int(fio, rot->nstsout);
759 snew(rot->grp, rot->ngrp);
761 for (g = 0; g < rot->ngrp; g++)
763 do_rotgrp(fio, &rot->grp[g], bRead);
768 static void do_swapcoords(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead)
773 gmx_fio_do_int(fio, swap->nat);
774 gmx_fio_do_int(fio, swap->nat_sol);
775 for (j = 0; j < 2; j++)
777 gmx_fio_do_int(fio, swap->nat_split[j]);
778 gmx_fio_do_int(fio, swap->massw_split[j]);
780 gmx_fio_do_int(fio, swap->nstswap);
781 gmx_fio_do_int(fio, swap->nAverage);
782 gmx_fio_do_real(fio, swap->threshold);
783 gmx_fio_do_real(fio, swap->cyl0r);
784 gmx_fio_do_real(fio, swap->cyl0u);
785 gmx_fio_do_real(fio, swap->cyl0l);
786 gmx_fio_do_real(fio, swap->cyl1r);
787 gmx_fio_do_real(fio, swap->cyl1u);
788 gmx_fio_do_real(fio, swap->cyl1l);
792 snew(swap->ind, swap->nat);
793 snew(swap->ind_sol, swap->nat_sol);
794 for (j = 0; j < 2; j++)
796 snew(swap->ind_split[j], swap->nat_split[j]);
800 gmx_fio_ndo_int(fio, swap->ind, swap->nat);
801 gmx_fio_ndo_int(fio, swap->ind_sol, swap->nat_sol);
802 for (j = 0; j < 2; j++)
804 gmx_fio_ndo_int(fio, swap->ind_split[j], swap->nat_split[j]);
807 for (j = 0; j < eCompNR; j++)
809 gmx_fio_do_int(fio, swap->nanions[j]);
810 gmx_fio_do_int(fio, swap->ncations[j]);
816 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
817 int file_version, real *fudgeQQ)
819 int i, j, k, *tmp, idum = 0;
822 gmx_bool bSimAnn, bdum = 0;
823 real zerotemptime, finish_t, init_temp, finish_temp;
825 if (file_version != tpx_version)
827 /* Give a warning about features that are not accessible */
828 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
829 file_version, tpx_version);
837 if (file_version == 0)
842 /* Basic inputrec stuff */
843 gmx_fio_do_int(fio, ir->eI);
844 if (file_version >= 62)
846 gmx_fio_do_int64(fio, ir->nsteps);
850 gmx_fio_do_int(fio, idum);
854 if (file_version > 25)
856 if (file_version >= 62)
858 gmx_fio_do_int64(fio, ir->init_step);
862 gmx_fio_do_int(fio, idum);
863 ir->init_step = idum;
871 if (file_version >= 58)
873 gmx_fio_do_int(fio, ir->simulation_part);
877 ir->simulation_part = 1;
880 if (file_version >= 67)
882 gmx_fio_do_int(fio, ir->nstcalcenergy);
886 ir->nstcalcenergy = 1;
888 if (file_version < 53)
890 /* The pbc info has been moved out of do_inputrec,
891 * since we always want it, also without reading the inputrec.
893 gmx_fio_do_int(fio, ir->ePBC);
894 if ((file_version <= 15) && (ir->ePBC == 2))
898 if (file_version >= 45)
900 gmx_fio_do_int(fio, ir->bPeriodicMols);
907 ir->bPeriodicMols = TRUE;
911 ir->bPeriodicMols = FALSE;
915 if (file_version >= 81)
917 gmx_fio_do_int(fio, ir->cutoff_scheme);
918 if (file_version < 94)
920 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
925 ir->cutoff_scheme = ecutsGROUP;
927 gmx_fio_do_int(fio, ir->ns_type);
928 gmx_fio_do_int(fio, ir->nstlist);
929 gmx_fio_do_int(fio, idum); /* used to be ndelta; not used anymore */
930 if (file_version < 41)
932 gmx_fio_do_int(fio, idum);
933 gmx_fio_do_int(fio, idum);
935 if (file_version >= 45)
937 gmx_fio_do_real(fio, ir->rtpi);
943 gmx_fio_do_int(fio, ir->nstcomm);
944 if (file_version > 34)
946 gmx_fio_do_int(fio, ir->comm_mode);
948 else if (ir->nstcomm < 0)
950 ir->comm_mode = ecmANGULAR;
954 ir->comm_mode = ecmLINEAR;
956 ir->nstcomm = abs(ir->nstcomm);
958 /* ignore nstcheckpoint */
959 if (file_version > 25 && file_version < tpxv_RemoveObsoleteParameters1)
961 gmx_fio_do_int(fio, idum);
964 gmx_fio_do_int(fio, ir->nstcgsteep);
966 if (file_version >= 30)
968 gmx_fio_do_int(fio, ir->nbfgscorr);
975 gmx_fio_do_int(fio, ir->nstlog);
976 gmx_fio_do_int(fio, ir->nstxout);
977 gmx_fio_do_int(fio, ir->nstvout);
978 gmx_fio_do_int(fio, ir->nstfout);
979 gmx_fio_do_int(fio, ir->nstenergy);
980 gmx_fio_do_int(fio, ir->nstxout_compressed);
981 if (file_version >= 59)
983 gmx_fio_do_double(fio, ir->init_t);
984 gmx_fio_do_double(fio, ir->delta_t);
988 gmx_fio_do_real(fio, rdum);
990 gmx_fio_do_real(fio, rdum);
993 gmx_fio_do_real(fio, ir->x_compression_precision);
994 if (file_version < 19)
996 gmx_fio_do_int(fio, idum);
997 gmx_fio_do_int(fio, idum);
999 if (file_version < 18)
1001 gmx_fio_do_int(fio, idum);
1003 if (file_version >= 81)
1005 gmx_fio_do_real(fio, ir->verletbuf_tol);
1009 ir->verletbuf_tol = 0;
1011 gmx_fio_do_real(fio, ir->rlist);
1012 if (file_version >= 67)
1014 gmx_fio_do_real(fio, ir->rlistlong);
1016 if (file_version >= 82 && file_version != 90)
1018 gmx_fio_do_int(fio, ir->nstcalclr);
1022 /* Calculate at NS steps */
1023 ir->nstcalclr = ir->nstlist;
1025 gmx_fio_do_int(fio, ir->coulombtype);
1026 if (file_version < 32 && ir->coulombtype == eelRF)
1028 ir->coulombtype = eelRF_NEC;
1030 if (file_version >= 81)
1032 gmx_fio_do_int(fio, ir->coulomb_modifier);
1036 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1038 gmx_fio_do_real(fio, ir->rcoulomb_switch);
1039 gmx_fio_do_real(fio, ir->rcoulomb);
1040 gmx_fio_do_int(fio, ir->vdwtype);
1041 if (file_version >= 81)
1043 gmx_fio_do_int(fio, ir->vdw_modifier);
1047 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1049 gmx_fio_do_real(fio, ir->rvdw_switch);
1050 gmx_fio_do_real(fio, ir->rvdw);
1051 if (file_version < 67)
1053 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1055 gmx_fio_do_int(fio, ir->eDispCorr);
1056 gmx_fio_do_real(fio, ir->epsilon_r);
1057 if (file_version >= 37)
1059 gmx_fio_do_real(fio, ir->epsilon_rf);
1063 if (EEL_RF(ir->coulombtype))
1065 ir->epsilon_rf = ir->epsilon_r;
1066 ir->epsilon_r = 1.0;
1070 ir->epsilon_rf = 1.0;
1073 if (file_version >= 29)
1075 gmx_fio_do_real(fio, ir->tabext);
1082 if (file_version > 25)
1084 gmx_fio_do_int(fio, ir->gb_algorithm);
1085 gmx_fio_do_int(fio, ir->nstgbradii);
1086 gmx_fio_do_real(fio, ir->rgbradii);
1087 gmx_fio_do_real(fio, ir->gb_saltconc);
1088 gmx_fio_do_int(fio, ir->implicit_solvent);
1092 ir->gb_algorithm = egbSTILL;
1095 ir->gb_saltconc = 0;
1096 ir->implicit_solvent = eisNO;
1098 if (file_version >= 55)
1100 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1101 gmx_fio_do_real(fio, ir->gb_obc_alpha);
1102 gmx_fio_do_real(fio, ir->gb_obc_beta);
1103 gmx_fio_do_real(fio, ir->gb_obc_gamma);
1104 if (file_version >= 60)
1106 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1107 gmx_fio_do_int(fio, ir->sa_algorithm);
1111 ir->gb_dielectric_offset = 0.009;
1112 ir->sa_algorithm = esaAPPROX;
1114 gmx_fio_do_real(fio, ir->sa_surface_tension);
1116 /* Override sa_surface_tension if it is not changed in the mpd-file */
1117 if (ir->sa_surface_tension < 0)
1119 if (ir->gb_algorithm == egbSTILL)
1121 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1123 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1125 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1132 /* Better use sensible values than insane (0.0) ones... */
1133 ir->gb_epsilon_solvent = 80;
1134 ir->gb_obc_alpha = 1.0;
1135 ir->gb_obc_beta = 0.8;
1136 ir->gb_obc_gamma = 4.85;
1137 ir->sa_surface_tension = 2.092;
1141 if (file_version >= 81)
1143 gmx_fio_do_real(fio, ir->fourier_spacing);
1147 ir->fourier_spacing = 0.0;
1149 gmx_fio_do_int(fio, ir->nkx);
1150 gmx_fio_do_int(fio, ir->nky);
1151 gmx_fio_do_int(fio, ir->nkz);
1152 gmx_fio_do_int(fio, ir->pme_order);
1153 gmx_fio_do_real(fio, ir->ewald_rtol);
1155 if (file_version >= 93)
1157 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1161 ir->ewald_rtol_lj = ir->ewald_rtol;
1164 if (file_version >= 24)
1166 gmx_fio_do_int(fio, ir->ewald_geometry);
1170 ir->ewald_geometry = eewg3D;
1173 if (file_version <= 17)
1175 ir->epsilon_surface = 0;
1176 if (file_version == 17)
1178 gmx_fio_do_int(fio, idum);
1183 gmx_fio_do_real(fio, ir->epsilon_surface);
1186 /* ignore bOptFFT */
1187 if (file_version < tpxv_RemoveObsoleteParameters1)
1189 gmx_fio_do_gmx_bool(fio, bdum);
1192 if (file_version >= 93)
1194 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1196 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1197 gmx_fio_do_int(fio, ir->etc);
1198 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1199 * but the values 0 and 1 still mean no and
1200 * berendsen temperature coupling, respectively.
1202 if (file_version >= 79)
1204 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1206 if (file_version >= 71)
1208 gmx_fio_do_int(fio, ir->nsttcouple);
1212 ir->nsttcouple = ir->nstcalcenergy;
1214 if (file_version <= 15)
1216 gmx_fio_do_int(fio, idum);
1218 if (file_version <= 17)
1220 gmx_fio_do_int(fio, ir->epct);
1221 if (file_version <= 15)
1225 ir->epct = epctSURFACETENSION;
1227 gmx_fio_do_int(fio, idum);
1230 /* we have removed the NO alternative at the beginning */
1234 ir->epct = epctISOTROPIC;
1238 ir->epc = epcBERENDSEN;
1243 gmx_fio_do_int(fio, ir->epc);
1244 gmx_fio_do_int(fio, ir->epct);
1246 if (file_version >= 71)
1248 gmx_fio_do_int(fio, ir->nstpcouple);
1252 ir->nstpcouple = ir->nstcalcenergy;
1254 gmx_fio_do_real(fio, ir->tau_p);
1255 if (file_version <= 15)
1257 gmx_fio_do_rvec(fio, vdum);
1258 clear_mat(ir->ref_p);
1259 for (i = 0; i < DIM; i++)
1261 ir->ref_p[i][i] = vdum[i];
1266 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1267 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1268 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1270 if (file_version <= 15)
1272 gmx_fio_do_rvec(fio, vdum);
1273 clear_mat(ir->compress);
1274 for (i = 0; i < DIM; i++)
1276 ir->compress[i][i] = vdum[i];
1281 gmx_fio_do_rvec(fio, ir->compress[XX]);
1282 gmx_fio_do_rvec(fio, ir->compress[YY]);
1283 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1285 if (file_version >= 47)
1287 gmx_fio_do_int(fio, ir->refcoord_scaling);
1288 gmx_fio_do_rvec(fio, ir->posres_com);
1289 gmx_fio_do_rvec(fio, ir->posres_comB);
1293 ir->refcoord_scaling = erscNO;
1294 clear_rvec(ir->posres_com);
1295 clear_rvec(ir->posres_comB);
1297 if ((file_version > 25) && (file_version < 79))
1299 gmx_fio_do_int(fio, ir->andersen_seed);
1303 ir->andersen_seed = 0;
1305 if (file_version < 26)
1307 gmx_fio_do_gmx_bool(fio, bSimAnn);
1308 gmx_fio_do_real(fio, zerotemptime);
1311 if (file_version < 37)
1313 gmx_fio_do_real(fio, rdum);
1316 gmx_fio_do_real(fio, ir->shake_tol);
1317 if (file_version < 54)
1319 gmx_fio_do_real(fio, *fudgeQQ);
1322 gmx_fio_do_int(fio, ir->efep);
1323 if (file_version <= 14 && ir->efep != efepNO)
1327 do_fepvals(fio, ir->fepvals, bRead, file_version);
1329 if (file_version >= 79)
1331 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1334 ir->bSimTemp = TRUE;
1339 ir->bSimTemp = FALSE;
1343 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1346 if (file_version >= 79)
1348 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1351 ir->bExpanded = TRUE;
1355 ir->bExpanded = FALSE;
1360 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1362 if (file_version >= 57)
1364 gmx_fio_do_int(fio, ir->eDisre);
1366 gmx_fio_do_int(fio, ir->eDisreWeighting);
1367 if (file_version < 22)
1369 if (ir->eDisreWeighting == 0)
1371 ir->eDisreWeighting = edrwEqual;
1375 ir->eDisreWeighting = edrwConservative;
1378 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1379 gmx_fio_do_real(fio, ir->dr_fc);
1380 gmx_fio_do_real(fio, ir->dr_tau);
1381 gmx_fio_do_int(fio, ir->nstdisreout);
1382 if (file_version >= 22)
1384 gmx_fio_do_real(fio, ir->orires_fc);
1385 gmx_fio_do_real(fio, ir->orires_tau);
1386 gmx_fio_do_int(fio, ir->nstorireout);
1392 ir->nstorireout = 0;
1395 /* ignore dihre_fc */
1396 if (file_version >= 26 && file_version < 79)
1398 gmx_fio_do_real(fio, rdum);
1399 if (file_version < 56)
1401 gmx_fio_do_real(fio, rdum);
1402 gmx_fio_do_int(fio, idum);
1406 gmx_fio_do_real(fio, ir->em_stepsize);
1407 gmx_fio_do_real(fio, ir->em_tol);
1408 if (file_version >= 22)
1410 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1414 ir->bShakeSOR = TRUE;
1416 if (file_version >= 11)
1418 gmx_fio_do_int(fio, ir->niter);
1423 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1426 if (file_version >= 21)
1428 gmx_fio_do_real(fio, ir->fc_stepsize);
1432 ir->fc_stepsize = 0;
1434 gmx_fio_do_int(fio, ir->eConstrAlg);
1435 gmx_fio_do_int(fio, ir->nProjOrder);
1436 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1437 if (file_version <= 14)
1439 gmx_fio_do_int(fio, idum);
1441 if (file_version >= 26)
1443 gmx_fio_do_int(fio, ir->nLincsIter);
1448 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1451 if (file_version < 33)
1453 gmx_fio_do_real(fio, bd_temp);
1455 gmx_fio_do_real(fio, ir->bd_fric);
1456 if (file_version >= tpxv_Use64BitRandomSeed)
1458 gmx_fio_do_int64(fio, ir->ld_seed);
1462 gmx_fio_do_int(fio, idum);
1465 if (file_version >= 33)
1467 for (i = 0; i < DIM; i++)
1469 gmx_fio_do_rvec(fio, ir->deform[i]);
1474 for (i = 0; i < DIM; i++)
1476 clear_rvec(ir->deform[i]);
1479 if (file_version >= 14)
1481 gmx_fio_do_real(fio, ir->cos_accel);
1487 gmx_fio_do_int(fio, ir->userint1);
1488 gmx_fio_do_int(fio, ir->userint2);
1489 gmx_fio_do_int(fio, ir->userint3);
1490 gmx_fio_do_int(fio, ir->userint4);
1491 gmx_fio_do_real(fio, ir->userreal1);
1492 gmx_fio_do_real(fio, ir->userreal2);
1493 gmx_fio_do_real(fio, ir->userreal3);
1494 gmx_fio_do_real(fio, ir->userreal4);
1497 if (file_version >= 77)
1499 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1504 snew(ir->adress, 1);
1506 gmx_fio_do_int(fio, ir->adress->type);
1507 gmx_fio_do_real(fio, ir->adress->const_wf);
1508 gmx_fio_do_real(fio, ir->adress->ex_width);
1509 gmx_fio_do_real(fio, ir->adress->hy_width);
1510 gmx_fio_do_int(fio, ir->adress->icor);
1511 gmx_fio_do_int(fio, ir->adress->site);
1512 gmx_fio_do_rvec(fio, ir->adress->refs);
1513 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1514 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1515 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1516 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1520 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1522 if (ir->adress->n_tf_grps > 0)
1524 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1528 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1530 if (ir->adress->n_energy_grps > 0)
1532 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1538 ir->bAdress = FALSE;
1542 if (file_version >= 48)
1544 gmx_fio_do_int(fio, ir->ePull);
1545 if (ir->ePull != epullNO)
1551 do_pull(fio, ir->pull, bRead, file_version);
1556 ir->ePull = epullNO;
1559 /* Enforced rotation */
1560 if (file_version >= 74)
1562 gmx_fio_do_int(fio, ir->bRot);
1563 if (ir->bRot == TRUE)
1569 do_rot(fio, ir->rot, bRead);
1577 /* Interactive molecular dynamics */
1578 if (file_version >= tpxv_InteractiveMolecularDynamics)
1580 gmx_fio_do_int(fio, ir->bIMD);
1581 if (TRUE == ir->bIMD)
1587 do_imd(fio, ir->imd, bRead);
1592 /* We don't support IMD sessions for old .tpr files */
1597 gmx_fio_do_int(fio, ir->opts.ngtc);
1598 if (file_version >= 69)
1600 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1604 ir->opts.nhchainlength = 1;
1606 gmx_fio_do_int(fio, ir->opts.ngacc);
1607 gmx_fio_do_int(fio, ir->opts.ngfrz);
1608 gmx_fio_do_int(fio, ir->opts.ngener);
1612 snew(ir->opts.nrdf, ir->opts.ngtc);
1613 snew(ir->opts.ref_t, ir->opts.ngtc);
1614 snew(ir->opts.annealing, ir->opts.ngtc);
1615 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1616 snew(ir->opts.anneal_time, ir->opts.ngtc);
1617 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1618 snew(ir->opts.tau_t, ir->opts.ngtc);
1619 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1620 snew(ir->opts.acc, ir->opts.ngacc);
1621 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1623 if (ir->opts.ngtc > 0)
1625 if (bRead && file_version < 13)
1627 snew(tmp, ir->opts.ngtc);
1628 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1629 for (i = 0; i < ir->opts.ngtc; i++)
1631 ir->opts.nrdf[i] = tmp[i];
1637 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1639 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1640 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1641 if (file_version < 33 && ir->eI == eiBD)
1643 for (i = 0; i < ir->opts.ngtc; i++)
1645 ir->opts.tau_t[i] = bd_temp;
1649 if (ir->opts.ngfrz > 0)
1651 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1653 if (ir->opts.ngacc > 0)
1655 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1657 if (file_version >= 12)
1659 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1660 ir->opts.ngener*ir->opts.ngener);
1663 if (bRead && file_version < 26)
1665 for (i = 0; i < ir->opts.ngtc; i++)
1669 ir->opts.annealing[i] = eannSINGLE;
1670 ir->opts.anneal_npoints[i] = 2;
1671 snew(ir->opts.anneal_time[i], 2);
1672 snew(ir->opts.anneal_temp[i], 2);
1673 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1674 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1675 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1676 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1677 ir->opts.anneal_time[i][0] = ir->init_t;
1678 ir->opts.anneal_time[i][1] = finish_t;
1679 ir->opts.anneal_temp[i][0] = init_temp;
1680 ir->opts.anneal_temp[i][1] = finish_temp;
1684 ir->opts.annealing[i] = eannNO;
1685 ir->opts.anneal_npoints[i] = 0;
1691 /* file version 26 or later */
1692 /* First read the lists with annealing and npoints for each group */
1693 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1694 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1695 for (j = 0; j < (ir->opts.ngtc); j++)
1697 k = ir->opts.anneal_npoints[j];
1700 snew(ir->opts.anneal_time[j], k);
1701 snew(ir->opts.anneal_temp[j], k);
1703 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1704 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1708 if (file_version >= 45)
1710 gmx_fio_do_int(fio, ir->nwall);
1711 gmx_fio_do_int(fio, ir->wall_type);
1712 if (file_version >= 50)
1714 gmx_fio_do_real(fio, ir->wall_r_linpot);
1718 ir->wall_r_linpot = -1;
1720 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1721 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1722 gmx_fio_do_real(fio, ir->wall_density[0]);
1723 gmx_fio_do_real(fio, ir->wall_density[1]);
1724 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1730 ir->wall_atomtype[0] = -1;
1731 ir->wall_atomtype[1] = -1;
1732 ir->wall_density[0] = 0;
1733 ir->wall_density[1] = 0;
1734 ir->wall_ewald_zfac = 3;
1736 /* Cosine stuff for electric fields */
1737 for (j = 0; (j < DIM); j++)
1739 gmx_fio_do_int(fio, ir->ex[j].n);
1740 gmx_fio_do_int(fio, ir->et[j].n);
1743 snew(ir->ex[j].a, ir->ex[j].n);
1744 snew(ir->ex[j].phi, ir->ex[j].n);
1745 snew(ir->et[j].a, ir->et[j].n);
1746 snew(ir->et[j].phi, ir->et[j].n);
1748 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1749 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1750 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1751 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1755 if (file_version >= tpxv_ComputationalElectrophysiology)
1757 gmx_fio_do_int(fio, ir->eSwapCoords);
1758 if (ir->eSwapCoords != eswapNO)
1764 do_swapcoords(fio, ir->swap, bRead);
1769 if (file_version >= 39)
1771 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1772 gmx_fio_do_int(fio, ir->QMMMscheme);
1773 gmx_fio_do_real(fio, ir->scalefactor);
1774 gmx_fio_do_int(fio, ir->opts.ngQM);
1777 snew(ir->opts.QMmethod, ir->opts.ngQM);
1778 snew(ir->opts.QMbasis, ir->opts.ngQM);
1779 snew(ir->opts.QMcharge, ir->opts.ngQM);
1780 snew(ir->opts.QMmult, ir->opts.ngQM);
1781 snew(ir->opts.bSH, ir->opts.ngQM);
1782 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1783 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1784 snew(ir->opts.SAon, ir->opts.ngQM);
1785 snew(ir->opts.SAoff, ir->opts.ngQM);
1786 snew(ir->opts.SAsteps, ir->opts.ngQM);
1787 snew(ir->opts.bOPT, ir->opts.ngQM);
1788 snew(ir->opts.bTS, ir->opts.ngQM);
1790 if (ir->opts.ngQM > 0)
1792 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1793 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1794 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1795 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1796 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1797 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1798 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1799 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1800 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1801 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1802 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1803 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1805 /* end of QMMM stuff */
1810 static void do_harm(t_fileio *fio, t_iparams *iparams)
1812 gmx_fio_do_real(fio, iparams->harmonic.rA);
1813 gmx_fio_do_real(fio, iparams->harmonic.krA);
1814 gmx_fio_do_real(fio, iparams->harmonic.rB);
1815 gmx_fio_do_real(fio, iparams->harmonic.krB);
1818 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1819 gmx_bool bRead, int file_version)
1826 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1836 do_harm(fio, iparams);
1837 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1839 /* Correct incorrect storage of parameters */
1840 iparams->pdihs.phiB = iparams->pdihs.phiA;
1841 iparams->pdihs.cpB = iparams->pdihs.cpA;
1845 gmx_fio_do_real(fio, iparams->harmonic.rA);
1846 gmx_fio_do_real(fio, iparams->harmonic.krA);
1848 case F_LINEAR_ANGLES:
1849 gmx_fio_do_real(fio, iparams->linangle.klinA);
1850 gmx_fio_do_real(fio, iparams->linangle.aA);
1851 gmx_fio_do_real(fio, iparams->linangle.klinB);
1852 gmx_fio_do_real(fio, iparams->linangle.aB);
1855 gmx_fio_do_real(fio, iparams->fene.bm);
1856 gmx_fio_do_real(fio, iparams->fene.kb);
1860 gmx_fio_do_real(fio, iparams->restraint.lowA);
1861 gmx_fio_do_real(fio, iparams->restraint.up1A);
1862 gmx_fio_do_real(fio, iparams->restraint.up2A);
1863 gmx_fio_do_real(fio, iparams->restraint.kA);
1864 gmx_fio_do_real(fio, iparams->restraint.lowB);
1865 gmx_fio_do_real(fio, iparams->restraint.up1B);
1866 gmx_fio_do_real(fio, iparams->restraint.up2B);
1867 gmx_fio_do_real(fio, iparams->restraint.kB);
1873 gmx_fio_do_real(fio, iparams->tab.kA);
1874 gmx_fio_do_int(fio, iparams->tab.table);
1875 gmx_fio_do_real(fio, iparams->tab.kB);
1877 case F_CROSS_BOND_BONDS:
1878 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1879 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1880 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1882 case F_CROSS_BOND_ANGLES:
1883 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1884 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1885 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1886 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1888 case F_UREY_BRADLEY:
1889 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1890 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1891 gmx_fio_do_real(fio, iparams->u_b.r13A);
1892 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1893 if (file_version >= 79)
1895 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1896 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1897 gmx_fio_do_real(fio, iparams->u_b.r13B);
1898 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1902 iparams->u_b.thetaB = iparams->u_b.thetaA;
1903 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1904 iparams->u_b.r13B = iparams->u_b.r13A;
1905 iparams->u_b.kUBB = iparams->u_b.kUBA;
1908 case F_QUARTIC_ANGLES:
1909 gmx_fio_do_real(fio, iparams->qangle.theta);
1910 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1913 gmx_fio_do_real(fio, iparams->bham.a);
1914 gmx_fio_do_real(fio, iparams->bham.b);
1915 gmx_fio_do_real(fio, iparams->bham.c);
1918 gmx_fio_do_real(fio, iparams->morse.b0A);
1919 gmx_fio_do_real(fio, iparams->morse.cbA);
1920 gmx_fio_do_real(fio, iparams->morse.betaA);
1921 if (file_version >= 79)
1923 gmx_fio_do_real(fio, iparams->morse.b0B);
1924 gmx_fio_do_real(fio, iparams->morse.cbB);
1925 gmx_fio_do_real(fio, iparams->morse.betaB);
1929 iparams->morse.b0B = iparams->morse.b0A;
1930 iparams->morse.cbB = iparams->morse.cbA;
1931 iparams->morse.betaB = iparams->morse.betaA;
1935 gmx_fio_do_real(fio, iparams->cubic.b0);
1936 gmx_fio_do_real(fio, iparams->cubic.kb);
1937 gmx_fio_do_real(fio, iparams->cubic.kcub);
1941 case F_POLARIZATION:
1942 gmx_fio_do_real(fio, iparams->polarize.alpha);
1945 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1946 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1947 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1950 if (file_version < 31)
1952 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1954 gmx_fio_do_real(fio, iparams->wpol.al_x);
1955 gmx_fio_do_real(fio, iparams->wpol.al_y);
1956 gmx_fio_do_real(fio, iparams->wpol.al_z);
1957 gmx_fio_do_real(fio, iparams->wpol.rOH);
1958 gmx_fio_do_real(fio, iparams->wpol.rHH);
1959 gmx_fio_do_real(fio, iparams->wpol.rOD);
1962 gmx_fio_do_real(fio, iparams->thole.a);
1963 gmx_fio_do_real(fio, iparams->thole.alpha1);
1964 gmx_fio_do_real(fio, iparams->thole.alpha2);
1965 gmx_fio_do_real(fio, iparams->thole.rfac);
1968 gmx_fio_do_real(fio, iparams->lj.c6);
1969 gmx_fio_do_real(fio, iparams->lj.c12);
1972 gmx_fio_do_real(fio, iparams->lj14.c6A);
1973 gmx_fio_do_real(fio, iparams->lj14.c12A);
1974 gmx_fio_do_real(fio, iparams->lj14.c6B);
1975 gmx_fio_do_real(fio, iparams->lj14.c12B);
1978 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1979 gmx_fio_do_real(fio, iparams->ljc14.qi);
1980 gmx_fio_do_real(fio, iparams->ljc14.qj);
1981 gmx_fio_do_real(fio, iparams->ljc14.c6);
1982 gmx_fio_do_real(fio, iparams->ljc14.c12);
1984 case F_LJC_PAIRS_NB:
1985 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1986 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1987 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1988 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1994 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1995 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1996 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1998 /* Read the incorrectly stored multiplicity */
1999 gmx_fio_do_real(fio, iparams->harmonic.rB);
2000 gmx_fio_do_real(fio, iparams->harmonic.krB);
2001 iparams->pdihs.phiB = iparams->pdihs.phiA;
2002 iparams->pdihs.cpB = iparams->pdihs.cpA;
2006 gmx_fio_do_real(fio, iparams->pdihs.phiB);
2007 gmx_fio_do_real(fio, iparams->pdihs.cpB);
2008 gmx_fio_do_int(fio, iparams->pdihs.mult);
2012 gmx_fio_do_real(fio, iparams->pdihs.phiA);
2013 gmx_fio_do_real(fio, iparams->pdihs.cpA);
2016 gmx_fio_do_int(fio, iparams->disres.label);
2017 gmx_fio_do_int(fio, iparams->disres.type);
2018 gmx_fio_do_real(fio, iparams->disres.low);
2019 gmx_fio_do_real(fio, iparams->disres.up1);
2020 gmx_fio_do_real(fio, iparams->disres.up2);
2021 gmx_fio_do_real(fio, iparams->disres.kfac);
2024 gmx_fio_do_int(fio, iparams->orires.ex);
2025 gmx_fio_do_int(fio, iparams->orires.label);
2026 gmx_fio_do_int(fio, iparams->orires.power);
2027 gmx_fio_do_real(fio, iparams->orires.c);
2028 gmx_fio_do_real(fio, iparams->orires.obs);
2029 gmx_fio_do_real(fio, iparams->orires.kfac);
2032 if (file_version < 82)
2034 gmx_fio_do_int(fio, idum);
2035 gmx_fio_do_int(fio, idum);
2037 gmx_fio_do_real(fio, iparams->dihres.phiA);
2038 gmx_fio_do_real(fio, iparams->dihres.dphiA);
2039 gmx_fio_do_real(fio, iparams->dihres.kfacA);
2040 if (file_version >= 82)
2042 gmx_fio_do_real(fio, iparams->dihres.phiB);
2043 gmx_fio_do_real(fio, iparams->dihres.dphiB);
2044 gmx_fio_do_real(fio, iparams->dihres.kfacB);
2048 iparams->dihres.phiB = iparams->dihres.phiA;
2049 iparams->dihres.dphiB = iparams->dihres.dphiA;
2050 iparams->dihres.kfacB = iparams->dihres.kfacA;
2054 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
2055 gmx_fio_do_rvec(fio, iparams->posres.fcA);
2056 if (bRead && file_version < 27)
2058 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
2059 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
2063 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
2064 gmx_fio_do_rvec(fio, iparams->posres.fcB);
2068 gmx_fio_do_int(fio, iparams->fbposres.geom);
2069 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
2070 gmx_fio_do_real(fio, iparams->fbposres.r);
2071 gmx_fio_do_real(fio, iparams->fbposres.k);
2074 gmx_fio_ndo_real(fio, iparams->cbtdihs.cbtcA, NR_CBTDIHS);
2077 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2078 if (file_version >= 25)
2080 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2084 /* Fourier dihedrals are internally represented
2085 * as Ryckaert-Bellemans since those are faster to compute.
2087 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2088 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2092 gmx_fio_do_real(fio, iparams->constr.dA);
2093 gmx_fio_do_real(fio, iparams->constr.dB);
2096 gmx_fio_do_real(fio, iparams->settle.doh);
2097 gmx_fio_do_real(fio, iparams->settle.dhh);
2100 gmx_fio_do_real(fio, iparams->vsite.a);
2105 gmx_fio_do_real(fio, iparams->vsite.a);
2106 gmx_fio_do_real(fio, iparams->vsite.b);
2111 gmx_fio_do_real(fio, iparams->vsite.a);
2112 gmx_fio_do_real(fio, iparams->vsite.b);
2113 gmx_fio_do_real(fio, iparams->vsite.c);
2116 gmx_fio_do_int(fio, iparams->vsiten.n);
2117 gmx_fio_do_real(fio, iparams->vsiten.a);
2122 /* We got rid of some parameters in version 68 */
2123 if (bRead && file_version < 68)
2125 gmx_fio_do_real(fio, rdum);
2126 gmx_fio_do_real(fio, rdum);
2127 gmx_fio_do_real(fio, rdum);
2128 gmx_fio_do_real(fio, rdum);
2130 gmx_fio_do_real(fio, iparams->gb.sar);
2131 gmx_fio_do_real(fio, iparams->gb.st);
2132 gmx_fio_do_real(fio, iparams->gb.pi);
2133 gmx_fio_do_real(fio, iparams->gb.gbr);
2134 gmx_fio_do_real(fio, iparams->gb.bmlt);
2137 gmx_fio_do_int(fio, iparams->cmap.cmapA);
2138 gmx_fio_do_int(fio, iparams->cmap.cmapB);
2141 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
2142 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
2146 gmx_fio_unset_comment(fio);
2150 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
2157 gmx_fio_set_comment(fio, interaction_function[ftype].name);
2159 if (file_version < 44)
2161 for (i = 0; i < MAXNODES; i++)
2163 gmx_fio_do_int(fio, idum);
2166 gmx_fio_do_int(fio, ilist->nr);
2169 snew(ilist->iatoms, ilist->nr);
2171 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2174 gmx_fio_unset_comment(fio);
2178 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2179 gmx_bool bRead, int file_version)
2184 gmx_fio_do_int(fio, ffparams->atnr);
2185 if (file_version < 57)
2187 gmx_fio_do_int(fio, idum);
2189 gmx_fio_do_int(fio, ffparams->ntypes);
2192 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2193 ffparams->atnr, ffparams->ntypes);
2197 snew(ffparams->functype, ffparams->ntypes);
2198 snew(ffparams->iparams, ffparams->ntypes);
2200 /* Read/write all the function types */
2201 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2204 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2207 if (file_version >= 66)
2209 gmx_fio_do_double(fio, ffparams->reppow);
2213 ffparams->reppow = 12.0;
2216 if (file_version >= 57)
2218 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2221 /* Check whether all these function types are supported by the code.
2222 * In practice the code is backwards compatible, which means that the
2223 * numbering may have to be altered from old numbering to new numbering
2225 for (i = 0; (i < ffparams->ntypes); i++)
2229 /* Loop over file versions */
2230 for (k = 0; (k < NFTUPD); k++)
2232 /* Compare the read file_version to the update table */
2233 if ((file_version < ftupd[k].fvnr) &&
2234 (ffparams->functype[i] >= ftupd[k].ftype))
2236 ffparams->functype[i] += 1;
2239 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2240 i, ffparams->functype[i],
2241 interaction_function[ftupd[k].ftype].longname);
2248 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2252 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2257 static void add_settle_atoms(t_ilist *ilist)
2261 /* Settle used to only store the first atom: add the other two */
2262 srenew(ilist->iatoms, 2*ilist->nr);
2263 for (i = ilist->nr/2-1; i >= 0; i--)
2265 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2266 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2267 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2268 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2270 ilist->nr = 2*ilist->nr;
2273 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2276 int i, j, renum[F_NRE];
2280 for (j = 0; (j < F_NRE); j++)
2285 for (k = 0; k < NFTUPD; k++)
2287 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2296 ilist[j].iatoms = NULL;
2300 do_ilist(fio, &ilist[j], bRead, file_version, j);
2301 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2303 add_settle_atoms(&ilist[j]);
2307 if (bRead && gmx_debug_at)
2308 pr_ilist(debug,0,interaction_function[j].longname,
2309 functype,&ilist[j],TRUE);
2314 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2315 gmx_bool bRead, int file_version)
2317 do_ffparams(fio, ffparams, bRead, file_version);
2319 if (file_version >= 54)
2321 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2324 do_ilists(fio, molt->ilist, bRead, file_version);
2327 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2329 int i, idum, dum_nra, *dum_a;
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 if (file_version < 51)
2341 gmx_fio_do_int(fio, dum_nra);
2345 if ((block->nalloc_index > 0) && (NULL != block->index))
2347 sfree(block->index);
2349 block->nalloc_index = block->nr+1;
2350 snew(block->index, block->nalloc_index);
2352 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2354 if (file_version < 51 && dum_nra > 0)
2356 snew(dum_a, dum_nra);
2357 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2362 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2367 if (file_version < 44)
2369 for (i = 0; i < MAXNODES; i++)
2371 gmx_fio_do_int(fio, idum);
2374 gmx_fio_do_int(fio, block->nr);
2375 gmx_fio_do_int(fio, block->nra);
2378 block->nalloc_index = block->nr+1;
2379 snew(block->index, block->nalloc_index);
2380 block->nalloc_a = block->nra;
2381 snew(block->a, block->nalloc_a);
2383 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2384 gmx_fio_ndo_int(fio, block->a, block->nra);
2387 /* This is a primitive routine to make it possible to translate atomic numbers
2388 * to element names when reading TPR files, without making the Gromacs library
2389 * directory a dependency on mdrun (which is the case if we need elements.dat).
2392 atomicnumber_to_element(int atomicnumber)
2396 /* This does not have to be complete, so we only include elements likely
2397 * to occur in PDB files.
2399 switch (atomicnumber)
2401 case 1: p = "H"; break;
2402 case 5: p = "B"; break;
2403 case 6: p = "C"; break;
2404 case 7: p = "N"; break;
2405 case 8: p = "O"; break;
2406 case 9: p = "F"; break;
2407 case 11: p = "Na"; break;
2408 case 12: p = "Mg"; break;
2409 case 15: p = "P"; break;
2410 case 16: p = "S"; break;
2411 case 17: p = "Cl"; break;
2412 case 18: p = "Ar"; break;
2413 case 19: p = "K"; break;
2414 case 20: p = "Ca"; break;
2415 case 25: p = "Mn"; break;
2416 case 26: p = "Fe"; break;
2417 case 28: p = "Ni"; break;
2418 case 29: p = "Cu"; break;
2419 case 30: p = "Zn"; break;
2420 case 35: p = "Br"; break;
2421 case 47: p = "Ag"; break;
2422 default: p = ""; break;
2428 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2429 int file_version, gmx_groups_t *groups, int atnr)
2434 gmx_fio_do_real(fio, atom->m);
2435 gmx_fio_do_real(fio, atom->q);
2436 gmx_fio_do_real(fio, atom->mB);
2437 gmx_fio_do_real(fio, atom->qB);
2438 gmx_fio_do_ushort(fio, atom->type);
2439 gmx_fio_do_ushort(fio, atom->typeB);
2440 gmx_fio_do_int(fio, atom->ptype);
2441 gmx_fio_do_int(fio, atom->resind);
2442 if (file_version >= 52)
2444 gmx_fio_do_int(fio, atom->atomnumber);
2447 /* Set element string from atomic number if present.
2448 * This routine returns an empty string if the name is not found.
2450 strncpy(atom->elem, atomicnumber_to_element(atom->atomnumber), 4);
2451 /* avoid warnings about potentially unterminated string */
2452 atom->elem[3] = '\0';
2457 atom->atomnumber = NOTSET;
2459 if (file_version < 23)
2463 else if (file_version < 39)
2472 if (file_version < 57)
2474 unsigned char uchar[egcNR];
2475 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2476 for (i = myngrp; (i < ngrp); i++)
2480 /* Copy the old data format to the groups struct */
2481 for (i = 0; i < ngrp; i++)
2483 groups->grpnr[i][atnr] = uchar[i];
2488 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2493 if (file_version < 23)
2497 else if (file_version < 39)
2506 for (j = 0; (j < ngrp); j++)
2510 gmx_fio_do_int(fio, grps[j].nr);
2513 snew(grps[j].nm_ind, grps[j].nr);
2515 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2520 snew(grps[j].nm_ind, grps[j].nr);
2525 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2531 gmx_fio_do_int(fio, ls);
2532 *nm = get_symtab_handle(symtab, ls);
2536 ls = lookup_symtab(symtab, *nm);
2537 gmx_fio_do_int(fio, ls);
2541 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2546 for (j = 0; (j < nstr); j++)
2548 do_symstr(fio, &(nm[j]), bRead, symtab);
2552 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2553 t_symtab *symtab, int file_version)
2557 for (j = 0; (j < n); j++)
2559 do_symstr(fio, &(ri[j].name), bRead, symtab);
2560 if (file_version >= 63)
2562 gmx_fio_do_int(fio, ri[j].nr);
2563 gmx_fio_do_uchar(fio, ri[j].ic);
2573 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2575 gmx_groups_t *groups)
2579 gmx_fio_do_int(fio, atoms->nr);
2580 gmx_fio_do_int(fio, atoms->nres);
2581 if (file_version < 57)
2583 gmx_fio_do_int(fio, groups->ngrpname);
2584 for (i = 0; i < egcNR; i++)
2586 groups->ngrpnr[i] = atoms->nr;
2587 snew(groups->grpnr[i], groups->ngrpnr[i]);
2592 snew(atoms->atom, atoms->nr);
2593 snew(atoms->atomname, atoms->nr);
2594 snew(atoms->atomtype, atoms->nr);
2595 snew(atoms->atomtypeB, atoms->nr);
2596 snew(atoms->resinfo, atoms->nres);
2597 if (file_version < 57)
2599 snew(groups->grpname, groups->ngrpname);
2601 atoms->pdbinfo = NULL;
2603 for (i = 0; (i < atoms->nr); i++)
2605 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2607 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2608 if (bRead && (file_version <= 20))
2610 for (i = 0; i < atoms->nr; i++)
2612 atoms->atomtype[i] = put_symtab(symtab, "?");
2613 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2618 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2619 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2621 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2623 if (file_version < 57)
2625 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2627 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2631 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2632 gmx_bool bRead, t_symtab *symtab,
2637 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2638 gmx_fio_do_int(fio, groups->ngrpname);
2641 snew(groups->grpname, groups->ngrpname);
2643 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2644 for (g = 0; g < egcNR; g++)
2646 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2647 if (groups->ngrpnr[g] == 0)
2651 groups->grpnr[g] = NULL;
2658 snew(groups->grpnr[g], groups->ngrpnr[g]);
2660 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2665 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2670 if (file_version > 25)
2672 gmx_fio_do_int(fio, atomtypes->nr);
2676 snew(atomtypes->radius, j);
2677 snew(atomtypes->vol, j);
2678 snew(atomtypes->surftens, j);
2679 snew(atomtypes->atomnumber, j);
2680 snew(atomtypes->gb_radius, j);
2681 snew(atomtypes->S_hct, j);
2683 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2684 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2685 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2686 if (file_version >= 40)
2688 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2690 if (file_version >= 60)
2692 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2693 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2698 /* File versions prior to 26 cannot do GBSA,
2699 * so they dont use this structure
2702 atomtypes->radius = NULL;
2703 atomtypes->vol = NULL;
2704 atomtypes->surftens = NULL;
2705 atomtypes->atomnumber = NULL;
2706 atomtypes->gb_radius = NULL;
2707 atomtypes->S_hct = NULL;
2711 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2717 gmx_fio_do_int(fio, symtab->nr);
2721 snew(symtab->symbuf, 1);
2722 symbuf = symtab->symbuf;
2723 symbuf->bufsize = nr;
2724 snew(symbuf->buf, nr);
2725 for (i = 0; (i < nr); i++)
2727 gmx_fio_do_string(fio, buf);
2728 symbuf->buf[i] = strdup(buf);
2733 symbuf = symtab->symbuf;
2734 while (symbuf != NULL)
2736 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2738 gmx_fio_do_string(fio, symbuf->buf[i]);
2741 symbuf = symbuf->next;
2745 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2750 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2752 int i, j, ngrid, gs, nelem;
2754 gmx_fio_do_int(fio, cmap_grid->ngrid);
2755 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2757 ngrid = cmap_grid->ngrid;
2758 gs = cmap_grid->grid_spacing;
2763 snew(cmap_grid->cmapdata, ngrid);
2765 for (i = 0; i < cmap_grid->ngrid; i++)
2767 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2771 for (i = 0; i < cmap_grid->ngrid; i++)
2773 for (j = 0; j < nelem; j++)
2775 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2776 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2777 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2778 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2784 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2786 int m, a, a0, a1, r;
2790 /* We always assign a new chain number, but save the chain id characters
2791 * for larger molecules.
2793 #define CHAIN_MIN_ATOMS 15
2797 for (m = 0; m < mols->nr; m++)
2799 a0 = mols->index[m];
2800 a1 = mols->index[m+1];
2801 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2810 for (a = a0; a < a1; a++)
2812 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2813 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2818 /* Blank out the chain id if there was only one chain */
2821 for (r = 0; r < atoms->nres; r++)
2823 atoms->resinfo[r].chainid = ' ';
2828 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2829 t_symtab *symtab, int file_version,
2830 gmx_groups_t *groups)
2834 if (file_version >= 57)
2836 do_symstr(fio, &(molt->name), bRead, symtab);
2839 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2841 if (bRead && gmx_debug_at)
2843 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2846 if (file_version >= 57)
2848 do_ilists(fio, molt->ilist, bRead, file_version);
2850 do_block(fio, &molt->cgs, bRead, file_version);
2851 if (bRead && gmx_debug_at)
2853 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2857 /* This used to be in the atoms struct */
2858 do_blocka(fio, &molt->excls, bRead, file_version);
2861 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2865 gmx_fio_do_int(fio, molb->type);
2866 gmx_fio_do_int(fio, molb->nmol);
2867 gmx_fio_do_int(fio, molb->natoms_mol);
2868 /* Position restraint coordinates */
2869 gmx_fio_do_int(fio, molb->nposres_xA);
2870 if (molb->nposres_xA > 0)
2874 snew(molb->posres_xA, molb->nposres_xA);
2876 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2878 gmx_fio_do_int(fio, molb->nposres_xB);
2879 if (molb->nposres_xB > 0)
2883 snew(molb->posres_xB, molb->nposres_xB);
2885 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2890 static t_block mtop_mols(gmx_mtop_t *mtop)
2896 for (mb = 0; mb < mtop->nmolblock; mb++)
2898 mols.nr += mtop->molblock[mb].nmol;
2900 mols.nalloc_index = mols.nr + 1;
2901 snew(mols.index, mols.nalloc_index);
2906 for (mb = 0; mb < mtop->nmolblock; mb++)
2908 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2910 a += mtop->molblock[mb].natoms_mol;
2919 static void add_posres_molblock(gmx_mtop_t *mtop)
2924 gmx_molblock_t *molb;
2927 /* posres reference positions are stored in ip->posres (if present) and
2928 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2929 posres.pos0A are identical to fbposres.pos0. */
2930 il = &mtop->moltype[0].ilist[F_POSRES];
2931 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2932 if (il->nr == 0 && ilfb->nr == 0)
2938 for (i = 0; i < il->nr; i += 2)
2940 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2941 am = max(am, il->iatoms[i+1]);
2942 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2943 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2944 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2949 /* This loop is required if we have only flat-bottomed posres:
2951 - bFE == FALSE (no B-state for flat-bottomed posres) */
2954 for (i = 0; i < ilfb->nr; i += 2)
2956 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2957 am = max(am, ilfb->iatoms[i+1]);
2960 /* Make the posres coordinate block end at a molecule end */
2962 while (am >= mtop->mols.index[mol+1])
2966 molb = &mtop->molblock[0];
2967 molb->nposres_xA = mtop->mols.index[mol+1];
2968 snew(molb->posres_xA, molb->nposres_xA);
2971 molb->nposres_xB = molb->nposres_xA;
2972 snew(molb->posres_xB, molb->nposres_xB);
2976 molb->nposres_xB = 0;
2978 for (i = 0; i < il->nr; i += 2)
2980 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2981 a = il->iatoms[i+1];
2982 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2983 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2984 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2987 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2988 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2989 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2994 /* If only flat-bottomed posres are present, take reference pos from them.
2995 Here: bFE == FALSE */
2996 for (i = 0; i < ilfb->nr; i += 2)
2998 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2999 a = ilfb->iatoms[i+1];
3000 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
3001 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
3002 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
3007 static void set_disres_npair(gmx_mtop_t *mtop)
3014 ip = mtop->ffparams.iparams;
3016 for (mt = 0; mt < mtop->nmoltype; mt++)
3018 il = &mtop->moltype[mt].ilist[F_DISRES];
3023 for (i = 0; i < il->nr; i += 3)
3026 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
3028 ip[a[i]].disres.npair = npair;
3036 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
3046 do_symtab(fio, &(mtop->symtab), bRead);
3049 pr_symtab(debug, 0, "symtab", &mtop->symtab);
3052 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
3054 if (file_version >= 57)
3056 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
3058 gmx_fio_do_int(fio, mtop->nmoltype);
3066 snew(mtop->moltype, mtop->nmoltype);
3067 if (file_version < 57)
3069 mtop->moltype[0].name = mtop->name;
3072 for (mt = 0; mt < mtop->nmoltype; mt++)
3074 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
3078 if (file_version >= 57)
3080 gmx_fio_do_int(fio, mtop->nmolblock);
3084 mtop->nmolblock = 1;
3088 snew(mtop->molblock, mtop->nmolblock);
3090 if (file_version >= 57)
3092 for (mb = 0; mb < mtop->nmolblock; mb++)
3094 do_molblock(fio, &mtop->molblock[mb], bRead);
3096 gmx_fio_do_int(fio, mtop->natoms);
3100 mtop->molblock[0].type = 0;
3101 mtop->molblock[0].nmol = 1;
3102 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
3103 mtop->molblock[0].nposres_xA = 0;
3104 mtop->molblock[0].nposres_xB = 0;
3107 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
3110 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
3113 if (file_version < 57)
3115 /* Debug statements are inside do_idef */
3116 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
3117 mtop->natoms = mtop->moltype[0].atoms.nr;
3120 if (file_version >= 65)
3122 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
3126 mtop->ffparams.cmap_grid.ngrid = 0;
3127 mtop->ffparams.cmap_grid.grid_spacing = 0;
3128 mtop->ffparams.cmap_grid.cmapdata = NULL;
3131 if (file_version >= 57)
3133 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
3136 if (file_version < 57)
3138 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
3139 if (bRead && gmx_debug_at)
3141 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
3143 do_block(fio, &mtop->mols, bRead, file_version);
3144 /* Add the posres coordinates to the molblock */
3145 add_posres_molblock(mtop);
3149 if (file_version >= 57)
3151 done_block(&mtop->mols);
3152 mtop->mols = mtop_mols(mtop);
3156 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
3160 if (file_version < 51)
3162 /* Here used to be the shake blocks */
3163 do_blocka(fio, &dumb, bRead, file_version);
3176 close_symtab(&(mtop->symtab));
3180 /* If TopOnlyOK is TRUE then we can read even future versions
3181 * of tpx files, provided the file_generation hasn't changed.
3182 * If it is FALSE, we need the inputrecord too, and bail out
3183 * if the file is newer than the program.
3185 * The version and generation if the topology (see top of this file)
3186 * are returned in the two last arguments.
3188 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
3190 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
3191 gmx_bool TopOnlyOK, int *file_version,
3192 int *file_generation)
3195 char file_tag[STRLEN];
3202 gmx_fio_checktype(fio);
3203 gmx_fio_setdebug(fio, bDebugMode());
3205 /* NEW! XDR tpb file */
3206 precision = sizeof(real);
3209 gmx_fio_do_string(fio, buf);
3210 if (strncmp(buf, "VERSION", 7))
3212 gmx_fatal(FARGS, "Can not read file %s,\n"
3213 " this file is from a Gromacs version which is older than 2.0\n"
3214 " Make a new one with grompp or use a gro or pdb file, if possible",
3215 gmx_fio_getname(fio));
3217 gmx_fio_do_int(fio, precision);
3218 bDouble = (precision == sizeof(double));
3219 if ((precision != sizeof(float)) && !bDouble)
3221 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3222 "instead of %d or %d",
3223 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3225 gmx_fio_setprecision(fio, bDouble);
3226 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3227 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3231 gmx_fio_write_string(fio, GromacsVersion());
3232 bDouble = (precision == sizeof(double));
3233 gmx_fio_setprecision(fio, bDouble);
3234 gmx_fio_do_int(fio, precision);
3236 sprintf(file_tag, "%s", tpx_tag);
3237 fgen = tpx_generation;
3240 /* Check versions! */
3241 gmx_fio_do_int(fio, fver);
3243 /* This is for backward compatibility with development versions 77-79
3244 * where the tag was, mistakenly, placed before the generation,
3245 * which would cause a segv instead of a proper error message
3246 * when reading the topology only from tpx with <77 code.
3248 if (fver >= 77 && fver <= 79)
3250 gmx_fio_do_string(fio, file_tag);
3255 gmx_fio_do_int(fio, fgen);
3264 gmx_fio_do_string(fio, file_tag);
3270 /* Versions before 77 don't have the tag, set it to release */
3271 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3274 if (strcmp(file_tag, tpx_tag) != 0)
3276 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3279 /* We only support reading tpx files with the same tag as the code
3280 * or tpx files with the release tag and with lower version number.
3282 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
3284 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3285 gmx_fio_getname(fio), fver, file_tag,
3286 tpx_version, tpx_tag);
3291 if (file_version != NULL)
3293 *file_version = fver;
3295 if (file_generation != NULL)
3297 *file_generation = fgen;
3301 if ((fver <= tpx_incompatible_version) ||
3302 ((fver > tpx_version) && !TopOnlyOK) ||
3303 (fgen > tpx_generation) ||
3304 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3306 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3307 gmx_fio_getname(fio), fver, tpx_version);
3310 do_section(fio, eitemHEADER, bRead);
3311 gmx_fio_do_int(fio, tpx->natoms);
3314 gmx_fio_do_int(fio, tpx->ngtc);
3322 gmx_fio_do_int(fio, idum);
3323 gmx_fio_do_real(fio, rdum);
3325 /*a better decision will eventually (5.0 or later) need to be made
3326 on how to treat the alchemical state of the system, which can now
3327 vary through a simulation, and cannot be completely described
3328 though a single lambda variable, or even a single state
3329 index. Eventually, should probably be a vector. MRS*/
3332 gmx_fio_do_int(fio, tpx->fep_state);
3334 gmx_fio_do_real(fio, tpx->lambda);
3335 gmx_fio_do_int(fio, tpx->bIr);
3336 gmx_fio_do_int(fio, tpx->bTop);
3337 gmx_fio_do_int(fio, tpx->bX);
3338 gmx_fio_do_int(fio, tpx->bV);
3339 gmx_fio_do_int(fio, tpx->bF);
3340 gmx_fio_do_int(fio, tpx->bBox);
3342 if ((fgen > tpx_generation))
3344 /* This can only happen if TopOnlyOK=TRUE */
3349 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3350 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3351 gmx_bool bXVallocated)
3357 int file_version, file_generation;
3361 gmx_bool bPeriodicMols;
3365 tpx.natoms = state->natoms;
3366 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3367 tpx.fep_state = state->fep_state;
3368 tpx.lambda = state->lambda[efptFEP];
3369 tpx.bIr = (ir != NULL);
3370 tpx.bTop = (mtop != NULL);
3371 tpx.bX = (state->x != NULL);
3372 tpx.bV = (state->v != NULL);
3373 tpx.bF = (f != NULL);
3377 TopOnlyOK = (ir == NULL);
3379 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3384 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3385 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3390 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3391 state->natoms = tpx.natoms;
3392 state->nalloc = tpx.natoms;
3398 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3402 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3404 do_test(fio, tpx.bBox, state->box);
3405 do_section(fio, eitemBOX, bRead);
3408 gmx_fio_ndo_rvec(fio, state->box, DIM);
3409 if (file_version >= 51)
3411 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3415 /* We initialize box_rel after reading the inputrec */
3416 clear_mat(state->box_rel);
3418 if (file_version >= 28)
3420 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3421 if (file_version < 56)
3424 gmx_fio_ndo_rvec(fio, mdum, DIM);
3429 if (state->ngtc > 0 && file_version >= 28)
3432 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3433 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3434 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3435 snew(dumv, state->ngtc);
3436 if (file_version < 69)
3438 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3440 /* These used to be the Berendsen tcoupl_lambda's */
3441 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3445 /* Prior to tpx version 26, the inputrec was here.
3446 * I moved it to enable partial forward-compatibility
3447 * for analysis/viewer programs.
3449 if (file_version < 26)
3451 do_test(fio, tpx.bIr, ir);
3452 do_section(fio, eitemIR, bRead);
3457 do_inputrec(fio, ir, bRead, file_version,
3458 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3461 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3466 do_inputrec(fio, &dum_ir, bRead, file_version,
3467 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3470 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3472 done_inputrec(&dum_ir);
3478 do_test(fio, tpx.bTop, mtop);
3479 do_section(fio, eitemTOP, bRead);
3484 do_mtop(fio, mtop, bRead, file_version);
3488 do_mtop(fio, &dum_top, bRead, file_version);
3489 done_mtop(&dum_top, TRUE);
3492 do_test(fio, tpx.bX, state->x);
3493 do_section(fio, eitemX, bRead);
3498 state->flags |= (1<<estX);
3500 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3503 do_test(fio, tpx.bV, state->v);
3504 do_section(fio, eitemV, bRead);
3509 state->flags |= (1<<estV);
3511 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3514 do_test(fio, tpx.bF, f);
3515 do_section(fio, eitemF, bRead);
3518 gmx_fio_ndo_rvec(fio, f, state->natoms);
3521 /* Starting with tpx version 26, we have the inputrec
3522 * at the end of the file, so we can ignore it
3523 * if the file is never than the software (but still the
3524 * same generation - see comments at the top of this file.
3529 bPeriodicMols = FALSE;
3530 if (file_version >= 26)
3532 do_test(fio, tpx.bIr, ir);
3533 do_section(fio, eitemIR, bRead);
3536 if (file_version >= 53)
3538 /* Removed the pbc info from do_inputrec, since we always want it */
3542 bPeriodicMols = ir->bPeriodicMols;
3544 gmx_fio_do_int(fio, ePBC);
3545 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3547 if (file_generation <= tpx_generation && ir)
3549 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3552 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3554 if (file_version < 51)
3556 set_box_rel(ir, state);
3558 if (file_version < 53)
3561 bPeriodicMols = ir->bPeriodicMols;
3564 if (bRead && ir && file_version >= 53)
3566 /* We need to do this after do_inputrec, since that initializes ir */
3568 ir->bPeriodicMols = bPeriodicMols;
3577 if (state->ngtc == 0)
3579 /* Reading old version without tcoupl state data: set it */
3580 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3582 if (tpx.bTop && mtop)
3584 if (file_version < 57)
3586 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3588 ir->eDisre = edrSimple;
3592 ir->eDisre = edrNone;
3595 set_disres_npair(mtop);
3599 if (tpx.bTop && mtop)
3601 gmx_mtop_finalize(mtop);
3604 if (file_version >= 57)
3608 env = getenv("GMX_NOCHARGEGROUPS");
3611 sscanf(env, "%d", &ienv);
3612 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3617 "Will make single atomic charge groups in non-solvent%s\n",
3618 ienv > 1 ? " and solvent" : "");
3619 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3621 fprintf(stderr, "\n");
3629 /************************************************************
3631 * The following routines are the exported ones
3633 ************************************************************/
3635 t_fileio *open_tpx(const char *fn, const char *mode)
3637 return gmx_fio_open(fn, mode);
3640 void close_tpx(t_fileio *fio)
3645 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3646 int *file_version, int *file_generation)
3650 fio = open_tpx(fn, "r");
3651 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3655 void write_tpx_state(const char *fn,
3656 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3660 fio = open_tpx(fn, "w");
3661 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3665 void read_tpx_state(const char *fn,
3666 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3670 fio = open_tpx(fn, "r");
3671 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3675 int read_tpx(const char *fn,
3676 t_inputrec *ir, matrix box, int *natoms,
3677 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3685 fio = open_tpx(fn, "r");
3686 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3688 *natoms = state.natoms;
3691 copy_mat(state.box, box);
3700 int read_tpx_top(const char *fn,
3701 t_inputrec *ir, matrix box, int *natoms,
3702 rvec *x, rvec *v, rvec *f, t_topology *top)
3708 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3710 *top = gmx_mtop_t_to_t_topology(&mtop);
3715 gmx_bool fn2bTPX(const char *file)
3717 switch (fn2ftp(file))
3728 static void done_gmx_groups_t(gmx_groups_t *g)
3732 for (i = 0; (i < egcNR); i++)
3734 if (NULL != g->grps[i].nm_ind)
3736 sfree(g->grps[i].nm_ind);
3737 g->grps[i].nm_ind = NULL;
3739 if (NULL != g->grpnr[i])
3745 /* The contents of this array is in symtab, don't free it here */
3749 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3750 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3753 int natoms, i, version, generation;
3754 gmx_bool bTop, bXNULL = FALSE;
3756 t_topology *topconv;
3759 bTop = fn2bTPX(infile);
3763 read_tpxheader(infile, &header, TRUE, &version, &generation);
3766 snew(*x, header.natoms);
3770 snew(*v, header.natoms);
3773 *ePBC = read_tpx(infile, NULL, box, &natoms,
3774 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3775 *top = gmx_mtop_t_to_t_topology(mtop);
3776 /* In this case we need to throw away the group data too */
3777 done_gmx_groups_t(&mtop->groups);
3779 strcpy(title, *top->name);
3780 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3784 get_stx_coordnum(infile, &natoms);
3785 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3796 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3804 aps = gmx_atomprop_init();
3805 for (i = 0; (i < natoms); i++)
3807 if (!gmx_atomprop_query(aps, epropMass,
3808 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3809 *top->atoms.atomname[i],
3810 &(top->atoms.atom[i].m)))
3814 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3815 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3816 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3817 *top->atoms.atomname[i]);
3821 gmx_atomprop_destroy(aps);
3823 top->idef.ntypes = -1;