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,2015,2016,2017,2018, 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.
39 /* This file is completely threadsafe - keep it that way! */
50 #include "gromacs/compat/make_unique.h"
51 #include "gromacs/fileio/filetypes.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/gmxfio-xdr.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdtypes/awh-history.h"
57 #include "gromacs/mdtypes/awh-params.h"
58 #include "gromacs/mdtypes/inputrec.h"
59 #include "gromacs/mdtypes/md_enums.h"
60 #include "gromacs/mdtypes/pull-params.h"
61 #include "gromacs/mdtypes/state.h"
62 #include "gromacs/pbcutil/boxutilities.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/topology/block.h"
65 #include "gromacs/topology/ifunc.h"
66 #include "gromacs/topology/mtop_util.h"
67 #include "gromacs/topology/symtab.h"
68 #include "gromacs/topology/topology.h"
69 #include "gromacs/utility/arraysize.h"
70 #include "gromacs/utility/baseversion.h"
71 #include "gromacs/utility/cstringutil.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/futil.h"
74 #include "gromacs/utility/gmxassert.h"
75 #include "gromacs/utility/keyvaluetreebuilder.h"
76 #include "gromacs/utility/keyvaluetreeserializer.h"
77 #include "gromacs/utility/smalloc.h"
78 #include "gromacs/utility/snprintf.h"
79 #include "gromacs/utility/txtdump.h"
81 #define TPX_TAG_RELEASE "release"
83 /*! \brief Tag string for the file format written to run input files
84 * written by this version of the code.
86 * Change this if you want to change the run input format in a feature
87 * branch. This ensures that there will not be different run input
88 * formats around which cannot be distinguished, while not causing
89 * problems rebasing the feature branch onto upstream changes. When
90 * merging with mainstream GROMACS, set this tag string back to
91 * TPX_TAG_RELEASE, and instead add an element to tpxv.
93 static const char *tpx_tag = TPX_TAG_RELEASE;
95 /*! \brief Enum of values that describe the contents of a tpr file
96 * whose format matches a version number
98 * The enum helps the code be more self-documenting and ensure merges
99 * do not silently resolve when two patches make the same bump. When
100 * adding new functionality, add a new element just above tpxv_Count
101 * in this enumeration, and write code below that does the right thing
102 * according to the value of file_version.
105 tpxv_ComputationalElectrophysiology = 96, /**< support for ion/water position swaps (computational electrophysiology) */
106 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to int64_t */
107 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
108 tpxv_InteractiveMolecularDynamics, /**< interactive molecular dynamics (IMD) */
109 tpxv_RemoveObsoleteParameters1, /**< remove optimize_fft, dihre_fc, nstcheckpoint */
110 tpxv_PullCoordTypeGeom, /**< add pull type and geometry per group and flat-bottom */
111 tpxv_PullGeomDirRel, /**< add pull geometry direction-relative */
112 tpxv_IntermolecularBondeds, /**< permit inter-molecular bonded interactions in the topology */
113 tpxv_CompElWithSwapLayerOffset, /**< added parameters for improved CompEl setups */
114 tpxv_CompElPolyatomicIonsAndMultipleIonTypes, /**< CompEl now can handle polyatomic ions and more than two types of ions */
115 tpxv_RemoveAdress, /**< removed support for AdResS */
116 tpxv_PullCoordNGroup, /**< add ngroup to pull coord */
117 tpxv_RemoveTwinRange, /**< removed support for twin-range interactions */
118 tpxv_ReplacePullPrintCOM12, /**< Replaced print-com-1, 2 with pull-print-com */
119 tpxv_PullExternalPotential, /**< Added pull type external potential */
120 tpxv_GenericParamsForElectricField, /**< Introduced KeyValueTree and moved electric field parameters */
121 tpxv_AcceleratedWeightHistogram, /**< sampling with accelerated weight histogram method (AWH) */
122 tpxv_RemoveImplicitSolvation, /**< removed support for implicit solvation */
123 tpxv_PullPrevStepCOMAsReference, /**< Enabled using the COM of the pull group of the last frame as reference for PBC */
124 tpxv_MimicQMMM, /**< Inroduced support for MiMiC QM/MM interface */
125 tpxv_PullAverage, /**< Added possibility to output average pull force and position */
126 tpxv_Count /**< the total number of tpxv versions */
129 /*! \brief Version number of the file format written to run input
130 * files by this version of the code.
132 * The tpx_version increases whenever the file format in the main
133 * development branch changes, due to an extension of the tpxv enum above.
134 * Backward compatibility for reading old run input files is maintained
135 * by checking this version number against that of the file and then using
136 * the correct code path.
138 * When developing a feature branch that needs to change the run input
139 * file format, change tpx_tag instead. */
140 static const int tpx_version = tpxv_Count - 1;
143 /* This number should only be increased when you edit the TOPOLOGY section
144 * or the HEADER of the tpx format.
145 * This way we can maintain forward compatibility too for all analysis tools
146 * and/or external programs that only need to know the atom/residue names,
147 * charges, and bond connectivity.
149 * It first appeared in tpx version 26, when I also moved the inputrecord
150 * to the end of the tpx file, so we can just skip it if we only
153 * In particular, it must be increased when adding new elements to
154 * ftupd, so that old code can read new .tpr files.
156 static const int tpx_generation = 26;
158 /* This number should be the most recent backwards incompatible version
159 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
161 static const int tpx_incompatible_version = 57; // GMX4.0 has version 58
165 /* Struct used to maintain tpx compatibility when function types are added */
167 int fvnr; /* file version number in which the function type first appeared */
168 int ftype; /* function type */
172 * TODO The following three lines make little sense, please clarify if
173 * you've had to work out how ftupd works.
175 * The entries should be ordered in:
176 * 1. ascending function type number
177 * 2. ascending file version number
179 * Because we support reading of old .tpr file versions (even when
180 * mdrun can no longer run the simulation), we need to be able to read
181 * obsolete t_interaction_function types. Any data read from such
182 * fields is discarded. Their names have _NOLONGERUSED appended to
183 * them to make things clear.
185 static const t_ftupd ftupd[] = {
186 { 70, F_RESTRBONDS },
187 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
188 { 76, F_LINEAR_ANGLES },
189 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
190 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
192 { 60, F_GB12_NOLONGERUSED },
193 { 61, F_GB13_NOLONGERUSED },
194 { 61, F_GB14_NOLONGERUSED },
195 { 72, F_GBPOL_NOLONGERUSED },
196 { 72, F_NPSOLVATION_NOLONGERUSED },
199 { 69, F_VTEMP_NOLONGERUSED},
201 { 76, F_ANHARM_POL },
204 { 79, F_DVDL_BONDED, },
205 { 79, F_DVDL_RESTRAINT },
206 { 79, F_DVDL_TEMPERATURE },
208 #define NFTUPD asize(ftupd)
210 /* Needed for backward compatibility */
213 /**************************************************************
215 * Now the higer level routines that do io of the structures and arrays
217 **************************************************************/
218 static void do_pullgrp_tpx_pre95(t_fileio *fio,
225 gmx_fio_do_int(fio, pgrp->nat);
228 snew(pgrp->ind, pgrp->nat);
230 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
231 gmx_fio_do_int(fio, pgrp->nweight);
234 snew(pgrp->weight, pgrp->nweight);
236 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
237 gmx_fio_do_int(fio, pgrp->pbcatom);
238 gmx_fio_do_rvec(fio, pcrd->vec);
239 clear_rvec(pcrd->origin);
240 gmx_fio_do_rvec(fio, tmp);
242 gmx_fio_do_real(fio, pcrd->rate);
243 gmx_fio_do_real(fio, pcrd->k);
244 gmx_fio_do_real(fio, pcrd->kB);
247 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
249 gmx_fio_do_int(fio, pgrp->nat);
252 snew(pgrp->ind, pgrp->nat);
254 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
255 gmx_fio_do_int(fio, pgrp->nweight);
258 snew(pgrp->weight, pgrp->nweight);
260 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
261 gmx_fio_do_int(fio, pgrp->pbcatom);
264 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd,
265 gmx_bool bRead, int file_version,
266 int ePullOld, int eGeomOld, ivec dimOld)
268 if (file_version >= tpxv_PullCoordNGroup)
270 gmx_fio_do_int(fio, pcrd->eType);
271 if (file_version >= tpxv_PullExternalPotential)
273 if (pcrd->eType == epullEXTERNAL)
279 gmx_fio_do_string(fio, buf);
280 pcrd->externalPotentialProvider = gmx_strdup(buf);
284 gmx_fio_do_string(fio, pcrd->externalPotentialProvider);
289 pcrd->externalPotentialProvider = nullptr;
296 pcrd->externalPotentialProvider = nullptr;
299 /* Note that we try to support adding new geometries without
300 * changing the tpx version. This requires checks when printing the
301 * geometry string and a check and fatal_error in init_pull.
303 gmx_fio_do_int(fio, pcrd->eGeom);
304 gmx_fio_do_int(fio, pcrd->ngroup);
305 if (pcrd->ngroup <= c_pullCoordNgroupMax)
307 gmx_fio_ndo_int(fio, pcrd->group, pcrd->ngroup);
311 /* More groups in file than supported, this must be a new geometry
312 * that is not supported by our current code. Since we will not
313 * use the groups for this coord (checks in the pull and WHAM code
314 * ensure this), we can ignore the groups and set ngroup=0.
317 snew(dum, pcrd->ngroup);
318 gmx_fio_ndo_int(fio, dum, pcrd->ngroup);
323 gmx_fio_do_ivec(fio, pcrd->dim);
328 gmx_fio_do_int(fio, pcrd->group[0]);
329 gmx_fio_do_int(fio, pcrd->group[1]);
330 if (file_version >= tpxv_PullCoordTypeGeom)
332 pcrd->ngroup = (pcrd->eGeom == epullgDIRRELATIVE ? 4 : 2);
333 gmx_fio_do_int(fio, pcrd->eType);
334 gmx_fio_do_int(fio, pcrd->eGeom);
335 if (pcrd->ngroup == 4)
337 gmx_fio_do_int(fio, pcrd->group[2]);
338 gmx_fio_do_int(fio, pcrd->group[3]);
340 gmx_fio_do_ivec(fio, pcrd->dim);
344 pcrd->eType = ePullOld;
345 pcrd->eGeom = eGeomOld;
346 copy_ivec(dimOld, pcrd->dim);
349 gmx_fio_do_rvec(fio, pcrd->origin);
350 gmx_fio_do_rvec(fio, pcrd->vec);
351 if (file_version >= tpxv_PullCoordTypeGeom)
353 gmx_fio_do_gmx_bool(fio, pcrd->bStart);
357 /* This parameter is only printed, but not actually used by mdrun */
358 pcrd->bStart = FALSE;
360 gmx_fio_do_real(fio, pcrd->init);
361 gmx_fio_do_real(fio, pcrd->rate);
362 gmx_fio_do_real(fio, pcrd->k);
363 gmx_fio_do_real(fio, pcrd->kB);
366 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
368 int n_lambda = fepvals->n_lambda;
370 /* reset the lambda calculation window */
371 fepvals->lambda_start_n = 0;
372 fepvals->lambda_stop_n = n_lambda;
373 if (file_version >= 79)
379 snew(expand->init_lambda_weights, n_lambda);
381 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
382 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
385 gmx_fio_do_int(fio, expand->nstexpanded);
386 gmx_fio_do_int(fio, expand->elmcmove);
387 gmx_fio_do_int(fio, expand->elamstats);
388 gmx_fio_do_int(fio, expand->lmc_repeats);
389 gmx_fio_do_int(fio, expand->gibbsdeltalam);
390 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
391 gmx_fio_do_int(fio, expand->lmc_seed);
392 gmx_fio_do_real(fio, expand->mc_temp);
393 gmx_fio_do_gmx_bool(fio, expand->bSymmetrizedTMatrix);
394 gmx_fio_do_int(fio, expand->nstTij);
395 gmx_fio_do_int(fio, expand->minvarmin);
396 gmx_fio_do_int(fio, expand->c_range);
397 gmx_fio_do_real(fio, expand->wl_scale);
398 gmx_fio_do_real(fio, expand->wl_ratio);
399 gmx_fio_do_real(fio, expand->init_wl_delta);
400 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
401 gmx_fio_do_int(fio, expand->elmceq);
402 gmx_fio_do_int(fio, expand->equil_steps);
403 gmx_fio_do_int(fio, expand->equil_samples);
404 gmx_fio_do_int(fio, expand->equil_n_at_lam);
405 gmx_fio_do_real(fio, expand->equil_wl_delta);
406 gmx_fio_do_real(fio, expand->equil_ratio);
410 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
413 if (file_version >= 79)
415 gmx_fio_do_int(fio, simtemp->eSimTempScale);
416 gmx_fio_do_real(fio, simtemp->simtemp_high);
417 gmx_fio_do_real(fio, simtemp->simtemp_low);
422 snew(simtemp->temperatures, n_lambda);
424 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
429 static void do_imd(t_fileio *fio, t_IMD *imd, gmx_bool bRead)
431 gmx_fio_do_int(fio, imd->nat);
434 snew(imd->ind, imd->nat);
436 gmx_fio_ndo_int(fio, imd->ind, imd->nat);
439 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
441 /* i is defined in the ndo_double macro; use g to iterate. */
445 /* free energy values */
447 if (file_version >= 79)
449 gmx_fio_do_int(fio, fepvals->init_fep_state);
450 gmx_fio_do_double(fio, fepvals->init_lambda);
451 gmx_fio_do_double(fio, fepvals->delta_lambda);
453 else if (file_version >= 59)
455 gmx_fio_do_double(fio, fepvals->init_lambda);
456 gmx_fio_do_double(fio, fepvals->delta_lambda);
460 gmx_fio_do_real(fio, rdum);
461 fepvals->init_lambda = rdum;
462 gmx_fio_do_real(fio, rdum);
463 fepvals->delta_lambda = rdum;
465 if (file_version >= 79)
467 gmx_fio_do_int(fio, fepvals->n_lambda);
470 snew(fepvals->all_lambda, efptNR);
472 for (g = 0; g < efptNR; g++)
474 if (fepvals->n_lambda > 0)
478 snew(fepvals->all_lambda[g], fepvals->n_lambda);
480 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
481 gmx_fio_ndo_gmx_bool(fio, fepvals->separate_dvdl, efptNR);
483 else if (fepvals->init_lambda >= 0)
485 fepvals->separate_dvdl[efptFEP] = TRUE;
489 else if (file_version >= 64)
491 gmx_fio_do_int(fio, fepvals->n_lambda);
496 snew(fepvals->all_lambda, efptNR);
497 /* still allocate the all_lambda array's contents. */
498 for (g = 0; g < efptNR; g++)
500 if (fepvals->n_lambda > 0)
502 snew(fepvals->all_lambda[g], fepvals->n_lambda);
506 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
508 if (fepvals->init_lambda >= 0)
512 fepvals->separate_dvdl[efptFEP] = TRUE;
516 /* copy the contents of the efptFEP lambda component to all
517 the other components */
518 for (g = 0; g < efptNR; g++)
520 for (h = 0; h < fepvals->n_lambda; h++)
524 fepvals->all_lambda[g][h] =
525 fepvals->all_lambda[efptFEP][h];
534 fepvals->n_lambda = 0;
535 fepvals->all_lambda = nullptr;
536 if (fepvals->init_lambda >= 0)
538 fepvals->separate_dvdl[efptFEP] = TRUE;
541 gmx_fio_do_real(fio, fepvals->sc_alpha);
542 gmx_fio_do_int(fio, fepvals->sc_power);
543 if (file_version >= 79)
545 gmx_fio_do_real(fio, fepvals->sc_r_power);
549 fepvals->sc_r_power = 6.0;
551 gmx_fio_do_real(fio, fepvals->sc_sigma);
554 if (file_version >= 71)
556 fepvals->sc_sigma_min = fepvals->sc_sigma;
560 fepvals->sc_sigma_min = 0;
563 if (file_version >= 79)
565 gmx_fio_do_gmx_bool(fio, fepvals->bScCoul);
569 fepvals->bScCoul = TRUE;
571 if (file_version >= 64)
573 gmx_fio_do_int(fio, fepvals->nstdhdl);
577 fepvals->nstdhdl = 1;
580 if (file_version >= 73)
582 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
583 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
587 fepvals->separate_dhdl_file = esepdhdlfileYES;
588 fepvals->dhdl_derivatives = edhdlderivativesYES;
590 if (file_version >= 71)
592 gmx_fio_do_int(fio, fepvals->dh_hist_size);
593 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
597 fepvals->dh_hist_size = 0;
598 fepvals->dh_hist_spacing = 0.1;
600 if (file_version >= 79)
602 gmx_fio_do_int(fio, fepvals->edHdLPrintEnergy);
606 fepvals->edHdLPrintEnergy = edHdLPrintEnergyNO;
609 /* handle lambda_neighbors */
610 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
612 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
613 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
614 (fepvals->init_lambda < 0) )
616 fepvals->lambda_start_n = (fepvals->init_fep_state -
617 fepvals->lambda_neighbors);
618 fepvals->lambda_stop_n = (fepvals->init_fep_state +
619 fepvals->lambda_neighbors + 1);
620 if (fepvals->lambda_start_n < 0)
622 fepvals->lambda_start_n = 0;;
624 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
626 fepvals->lambda_stop_n = fepvals->n_lambda;
631 fepvals->lambda_start_n = 0;
632 fepvals->lambda_stop_n = fepvals->n_lambda;
637 fepvals->lambda_start_n = 0;
638 fepvals->lambda_stop_n = fepvals->n_lambda;
642 static void do_awhBias(t_fileio *fio, gmx::AwhBiasParams *awhBiasParams, gmx_bool bRead,
643 int gmx_unused file_version)
645 gmx_fio_do_int(fio, awhBiasParams->eTarget);
646 gmx_fio_do_double(fio, awhBiasParams->targetBetaScaling);
647 gmx_fio_do_double(fio, awhBiasParams->targetCutoff);
648 gmx_fio_do_int(fio, awhBiasParams->eGrowth);
649 gmx_fio_do_int(fio, awhBiasParams->bUserData);
650 gmx_fio_do_double(fio, awhBiasParams->errorInitial);
651 gmx_fio_do_int(fio, awhBiasParams->ndim);
652 gmx_fio_do_int(fio, awhBiasParams->shareGroup);
653 gmx_fio_do_gmx_bool(fio, awhBiasParams->equilibrateHistogram);
657 snew(awhBiasParams->dimParams, awhBiasParams->ndim);
660 for (int d = 0; d < awhBiasParams->ndim; d++)
662 gmx::AwhDimParams *dimParams = &awhBiasParams->dimParams[d];
664 gmx_fio_do_int(fio, dimParams->eCoordProvider);
665 gmx_fio_do_int(fio, dimParams->coordIndex);
666 gmx_fio_do_double(fio, dimParams->origin);
667 gmx_fio_do_double(fio, dimParams->end);
668 gmx_fio_do_double(fio, dimParams->period);
669 gmx_fio_do_double(fio, dimParams->forceConstant);
670 gmx_fio_do_double(fio, dimParams->diffusion);
671 gmx_fio_do_double(fio, dimParams->coordValueInit);
672 gmx_fio_do_double(fio, dimParams->coverDiameter);
676 static void do_awh(t_fileio *fio, gmx::AwhParams *awhParams, gmx_bool bRead,
677 int gmx_unused file_version)
679 gmx_fio_do_int(fio, awhParams->numBias);
680 gmx_fio_do_int(fio, awhParams->nstOut);
681 gmx_fio_do_int64(fio, awhParams->seed);
682 gmx_fio_do_int(fio, awhParams->nstSampleCoord);
683 gmx_fio_do_int(fio, awhParams->numSamplesUpdateFreeEnergy);
684 gmx_fio_do_int(fio, awhParams->ePotential);
685 gmx_fio_do_gmx_bool(fio, awhParams->shareBiasMultisim);
687 if (awhParams->numBias > 0)
691 snew(awhParams->awhBiasParams, awhParams->numBias);
694 for (int k = 0; k < awhParams->numBias; k++)
696 do_awhBias(fio, &awhParams->awhBiasParams[k], bRead, file_version);
701 static void do_pull(t_fileio *fio, pull_params_t *pull, gmx_bool bRead,
702 int file_version, int ePullOld)
708 if (file_version >= 95)
710 gmx_fio_do_int(fio, pull->ngroup);
712 gmx_fio_do_int(fio, pull->ncoord);
713 if (file_version < 95)
715 pull->ngroup = pull->ncoord + 1;
717 if (file_version < tpxv_PullCoordTypeGeom)
721 gmx_fio_do_int(fio, eGeomOld);
722 gmx_fio_do_ivec(fio, dimOld);
723 /* The inner cylinder radius, now removed */
724 gmx_fio_do_real(fio, dum);
726 gmx_fio_do_real(fio, pull->cylinder_r);
727 gmx_fio_do_real(fio, pull->constr_tol);
728 if (file_version >= 95)
730 gmx_fio_do_gmx_bool(fio, pull->bPrintCOM);
731 /* With file_version < 95 this value is set below */
733 if (file_version >= tpxv_ReplacePullPrintCOM12)
735 gmx_fio_do_gmx_bool(fio, pull->bPrintRefValue);
736 gmx_fio_do_gmx_bool(fio, pull->bPrintComp);
738 else if (file_version >= tpxv_PullCoordTypeGeom)
741 gmx_fio_do_int(fio, idum); /* used to be bPrintCOM2 */
742 gmx_fio_do_gmx_bool(fio, pull->bPrintRefValue);
743 gmx_fio_do_gmx_bool(fio, pull->bPrintComp);
747 pull->bPrintRefValue = FALSE;
748 pull->bPrintComp = TRUE;
750 gmx_fio_do_int(fio, pull->nstxout);
751 gmx_fio_do_int(fio, pull->nstfout);
752 if (file_version >= tpxv_PullPrevStepCOMAsReference)
754 gmx_fio_do_gmx_bool(fio, pull->bSetPbcRefToPrevStepCOM);
758 pull->bSetPbcRefToPrevStepCOM = FALSE;
762 snew(pull->group, pull->ngroup);
763 snew(pull->coord, pull->ncoord);
765 if (file_version < 95)
767 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
768 if (eGeomOld == epullgDIRPBC)
770 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
772 if (eGeomOld > epullgDIRPBC)
777 for (g = 0; g < pull->ngroup; g++)
779 /* We read and ignore a pull coordinate for group 0 */
780 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[std::max(g-1, 0)],
784 pull->coord[g-1].group[0] = 0;
785 pull->coord[g-1].group[1] = g;
789 pull->bPrintCOM = (pull->group[0].nat > 0);
793 for (g = 0; g < pull->ngroup; g++)
795 do_pull_group(fio, &pull->group[g], bRead);
797 for (g = 0; g < pull->ncoord; g++)
799 do_pull_coord(fio, &pull->coord[g],
800 bRead, file_version, ePullOld, eGeomOld, dimOld);
803 if (file_version >= tpxv_PullAverage)
807 v = pull->bXOutAverage;
808 gmx_fio_do_gmx_bool(fio, v);
809 pull->bXOutAverage = v;
810 v = pull->bFOutAverage;
811 gmx_fio_do_gmx_bool(fio, v);
812 pull->bFOutAverage = v;
817 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
819 gmx_fio_do_int(fio, rotg->eType);
820 gmx_fio_do_int(fio, rotg->bMassW);
821 gmx_fio_do_int(fio, rotg->nat);
824 snew(rotg->ind, rotg->nat);
826 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
829 snew(rotg->x_ref, rotg->nat);
831 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
832 gmx_fio_do_rvec(fio, rotg->inputVec);
833 gmx_fio_do_rvec(fio, rotg->pivot);
834 gmx_fio_do_real(fio, rotg->rate);
835 gmx_fio_do_real(fio, rotg->k);
836 gmx_fio_do_real(fio, rotg->slab_dist);
837 gmx_fio_do_real(fio, rotg->min_gaussian);
838 gmx_fio_do_real(fio, rotg->eps);
839 gmx_fio_do_int(fio, rotg->eFittype);
840 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
841 gmx_fio_do_real(fio, rotg->PotAngle_step);
844 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
848 gmx_fio_do_int(fio, rot->ngrp);
849 gmx_fio_do_int(fio, rot->nstrout);
850 gmx_fio_do_int(fio, rot->nstsout);
853 snew(rot->grp, rot->ngrp);
855 for (g = 0; g < rot->ngrp; g++)
857 do_rotgrp(fio, &rot->grp[g], bRead);
862 static void do_swapgroup(t_fileio *fio, t_swapGroup *g, gmx_bool bRead)
865 /* Name of the group or molecule */
870 gmx_fio_do_string(fio, buf);
871 g->molname = gmx_strdup(buf);
875 gmx_fio_do_string(fio, g->molname);
878 /* Number of atoms in the group */
879 gmx_fio_do_int(fio, g->nat);
881 /* The group's atom indices */
884 snew(g->ind, g->nat);
886 gmx_fio_ndo_int(fio, g->ind, g->nat);
888 /* Requested counts for compartments A and B */
889 gmx_fio_ndo_int(fio, g->nmolReq, eCompNR);
892 static void do_swapcoords_tpx(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead, int file_version)
894 /* Enums for better readability of the code */
899 eChannel0 = 0, eChannel1
903 if (file_version >= tpxv_CompElPolyatomicIonsAndMultipleIonTypes)
905 /* The total number of swap groups is the sum of the fixed groups
906 * (split0, split1, solvent), and the user-defined groups (2+ types of ions)
908 gmx_fio_do_int(fio, swap->ngrp);
911 snew(swap->grp, swap->ngrp);
913 for (int ig = 0; ig < swap->ngrp; ig++)
915 do_swapgroup(fio, &swap->grp[ig], bRead);
917 gmx_fio_do_gmx_bool(fio, swap->massw_split[eChannel0]);
918 gmx_fio_do_gmx_bool(fio, swap->massw_split[eChannel1]);
919 gmx_fio_do_int(fio, swap->nstswap);
920 gmx_fio_do_int(fio, swap->nAverage);
921 gmx_fio_do_real(fio, swap->threshold);
922 gmx_fio_do_real(fio, swap->cyl0r);
923 gmx_fio_do_real(fio, swap->cyl0u);
924 gmx_fio_do_real(fio, swap->cyl0l);
925 gmx_fio_do_real(fio, swap->cyl1r);
926 gmx_fio_do_real(fio, swap->cyl1u);
927 gmx_fio_do_real(fio, swap->cyl1l);
931 /*** Support reading older CompEl .tpr files ***/
933 /* In the original CompEl .tpr files, we always have 5 groups: */
935 snew(swap->grp, swap->ngrp);
937 swap->grp[eGrpSplit0 ].molname = gmx_strdup("split0" ); // group 0: split0
938 swap->grp[eGrpSplit1 ].molname = gmx_strdup("split1" ); // group 1: split1
939 swap->grp[eGrpSolvent].molname = gmx_strdup("solvent"); // group 2: solvent
940 swap->grp[3 ].molname = gmx_strdup("anions" ); // group 3: anions
941 swap->grp[4 ].molname = gmx_strdup("cations"); // group 4: cations
943 gmx_fio_do_int(fio, swap->grp[3].nat);
944 gmx_fio_do_int(fio, swap->grp[eGrpSolvent].nat);
945 gmx_fio_do_int(fio, swap->grp[eGrpSplit0].nat);
946 gmx_fio_do_gmx_bool(fio, swap->massw_split[eChannel0]);
947 gmx_fio_do_int(fio, swap->grp[eGrpSplit1].nat);
948 gmx_fio_do_gmx_bool(fio, swap->massw_split[eChannel1]);
949 gmx_fio_do_int(fio, swap->nstswap);
950 gmx_fio_do_int(fio, swap->nAverage);
951 gmx_fio_do_real(fio, swap->threshold);
952 gmx_fio_do_real(fio, swap->cyl0r);
953 gmx_fio_do_real(fio, swap->cyl0u);
954 gmx_fio_do_real(fio, swap->cyl0l);
955 gmx_fio_do_real(fio, swap->cyl1r);
956 gmx_fio_do_real(fio, swap->cyl1u);
957 gmx_fio_do_real(fio, swap->cyl1l);
959 // The order[] array keeps compatibility with older .tpr files
960 // by reading in the groups in the classic order
962 const int order[4] = {3, eGrpSolvent, eGrpSplit0, eGrpSplit1 };
964 for (int ig = 0; ig < 4; ig++)
967 snew(swap->grp[g].ind, swap->grp[g].nat);
968 gmx_fio_ndo_int(fio, swap->grp[g].ind, swap->grp[g].nat);
972 for (int j = eCompA; j <= eCompB; j++)
974 gmx_fio_do_int(fio, swap->grp[3].nmolReq[j]); // group 3 = anions
975 gmx_fio_do_int(fio, swap->grp[4].nmolReq[j]); // group 4 = cations
977 } /* End support reading older CompEl .tpr files */
979 if (file_version >= tpxv_CompElWithSwapLayerOffset)
981 gmx_fio_do_real(fio, swap->bulkOffset[eCompA]);
982 gmx_fio_do_real(fio, swap->bulkOffset[eCompB]);
987 static void do_legacy_efield(t_fileio *fio, gmx::KeyValueTreeObjectBuilder *root)
989 const char *const dimName[] = { "x", "y", "z" };
991 auto appliedForcesObj = root->addObject("applied-forces");
992 auto efieldObj = appliedForcesObj.addObject("electric-field");
993 // The content of the tpr file for this feature has
994 // been the same since gromacs 4.0 that was used for
996 for (int j = 0; j < DIM; ++j)
999 gmx_fio_do_int(fio, n);
1000 gmx_fio_do_int(fio, nt);
1001 std::vector<real> aa(n+1), phi(nt+1), at(nt+1), phit(nt+1);
1002 gmx_fio_ndo_real(fio, aa.data(), n);
1003 gmx_fio_ndo_real(fio, phi.data(), n);
1004 gmx_fio_ndo_real(fio, at.data(), nt);
1005 gmx_fio_ndo_real(fio, phit.data(), nt);
1008 if (n > 1 || nt > 1)
1010 gmx_fatal(FARGS, "Can not handle tpr files with more than one electric field term per direction.");
1012 auto dimObj = efieldObj.addObject(dimName[j]);
1013 dimObj.addValue<real>("E0", aa[0]);
1014 dimObj.addValue<real>("omega", at[0]);
1015 dimObj.addValue<real>("t0", phi[0]);
1016 dimObj.addValue<real>("sigma", phit[0]);
1022 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
1025 int i, j, k, idum = 0;
1027 gmx_bool bdum = false;
1029 if (file_version != tpx_version)
1031 /* Give a warning about features that are not accessible */
1032 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
1033 file_version, tpx_version);
1036 if (file_version == 0)
1041 gmx::KeyValueTreeBuilder paramsBuilder;
1042 gmx::KeyValueTreeObjectBuilder paramsObj = paramsBuilder.rootObject();
1044 /* Basic inputrec stuff */
1045 gmx_fio_do_int(fio, ir->eI);
1046 if (file_version >= 62)
1048 gmx_fio_do_int64(fio, ir->nsteps);
1052 gmx_fio_do_int(fio, idum);
1056 if (file_version >= 62)
1058 gmx_fio_do_int64(fio, ir->init_step);
1062 gmx_fio_do_int(fio, idum);
1063 ir->init_step = idum;
1066 gmx_fio_do_int(fio, ir->simulation_part);
1068 if (file_version >= 67)
1070 gmx_fio_do_int(fio, ir->nstcalcenergy);
1074 ir->nstcalcenergy = 1;
1076 if (file_version >= 81)
1078 gmx_fio_do_int(fio, ir->cutoff_scheme);
1079 if (file_version < 94)
1081 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
1086 ir->cutoff_scheme = ecutsGROUP;
1088 gmx_fio_do_int(fio, ir->ns_type);
1089 gmx_fio_do_int(fio, ir->nstlist);
1090 gmx_fio_do_int(fio, idum); /* used to be ndelta; not used anymore */
1092 gmx_fio_do_real(fio, ir->rtpi);
1094 gmx_fio_do_int(fio, ir->nstcomm);
1095 gmx_fio_do_int(fio, ir->comm_mode);
1097 /* ignore nstcheckpoint */
1098 if (file_version < tpxv_RemoveObsoleteParameters1)
1100 gmx_fio_do_int(fio, idum);
1103 gmx_fio_do_int(fio, ir->nstcgsteep);
1105 gmx_fio_do_int(fio, ir->nbfgscorr);
1107 gmx_fio_do_int(fio, ir->nstlog);
1108 gmx_fio_do_int(fio, ir->nstxout);
1109 gmx_fio_do_int(fio, ir->nstvout);
1110 gmx_fio_do_int(fio, ir->nstfout);
1111 gmx_fio_do_int(fio, ir->nstenergy);
1112 gmx_fio_do_int(fio, ir->nstxout_compressed);
1113 if (file_version >= 59)
1115 gmx_fio_do_double(fio, ir->init_t);
1116 gmx_fio_do_double(fio, ir->delta_t);
1120 gmx_fio_do_real(fio, rdum);
1122 gmx_fio_do_real(fio, rdum);
1125 gmx_fio_do_real(fio, ir->x_compression_precision);
1126 if (file_version >= 81)
1128 gmx_fio_do_real(fio, ir->verletbuf_tol);
1132 ir->verletbuf_tol = 0;
1134 gmx_fio_do_real(fio, ir->rlist);
1135 if (file_version >= 67 && file_version < tpxv_RemoveTwinRange)
1139 // Reading such a file version could invoke the twin-range
1140 // scheme, about which mdrun should give a fatal error.
1141 real dummy_rlistlong = -1;
1142 gmx_fio_do_real(fio, dummy_rlistlong);
1144 if (ir->rlist > 0 && (dummy_rlistlong == 0 || dummy_rlistlong > ir->rlist))
1146 // Get mdrun to issue an error (regardless of
1147 // ir->cutoff_scheme).
1148 ir->useTwinRange = true;
1152 // grompp used to set rlistlong actively. Users were
1153 // probably also confused and set rlistlong == rlist.
1154 // However, in all remaining cases, it is safe to let
1155 // mdrun proceed normally.
1156 ir->useTwinRange = false;
1162 // No need to read or write anything
1163 ir->useTwinRange = false;
1165 if (file_version >= 82 && file_version != 90)
1167 // Multiple time-stepping is no longer enabled, but the old
1168 // support required the twin-range scheme, for which mdrun
1169 // already emits a fatal error.
1170 int dummy_nstcalclr = -1;
1171 gmx_fio_do_int(fio, dummy_nstcalclr);
1173 gmx_fio_do_int(fio, ir->coulombtype);
1174 if (file_version >= 81)
1176 gmx_fio_do_int(fio, ir->coulomb_modifier);
1180 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1182 gmx_fio_do_real(fio, ir->rcoulomb_switch);
1183 gmx_fio_do_real(fio, ir->rcoulomb);
1184 gmx_fio_do_int(fio, ir->vdwtype);
1185 if (file_version >= 81)
1187 gmx_fio_do_int(fio, ir->vdw_modifier);
1191 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1193 gmx_fio_do_real(fio, ir->rvdw_switch);
1194 gmx_fio_do_real(fio, ir->rvdw);
1195 gmx_fio_do_int(fio, ir->eDispCorr);
1196 gmx_fio_do_real(fio, ir->epsilon_r);
1197 gmx_fio_do_real(fio, ir->epsilon_rf);
1198 gmx_fio_do_real(fio, ir->tabext);
1200 // This permits reading a .tpr file that used implicit solvent,
1201 // and later permitting mdrun to refuse to run it.
1204 if (file_version < tpxv_RemoveImplicitSolvation)
1206 gmx_fio_do_int(fio, idum);
1207 gmx_fio_do_int(fio, idum);
1208 gmx_fio_do_real(fio, rdum);
1209 gmx_fio_do_real(fio, rdum);
1210 gmx_fio_do_int(fio, idum);
1211 ir->implicit_solvent = (idum > 0);
1215 ir->implicit_solvent = false;
1217 if (file_version < tpxv_RemoveImplicitSolvation)
1219 gmx_fio_do_real(fio, rdum);
1220 gmx_fio_do_real(fio, rdum);
1221 gmx_fio_do_real(fio, rdum);
1222 gmx_fio_do_real(fio, rdum);
1223 if (file_version >= 60)
1225 gmx_fio_do_real(fio, rdum);
1226 gmx_fio_do_int(fio, idum);
1228 gmx_fio_do_real(fio, rdum);
1232 if (file_version >= 81)
1234 gmx_fio_do_real(fio, ir->fourier_spacing);
1238 ir->fourier_spacing = 0.0;
1240 gmx_fio_do_int(fio, ir->nkx);
1241 gmx_fio_do_int(fio, ir->nky);
1242 gmx_fio_do_int(fio, ir->nkz);
1243 gmx_fio_do_int(fio, ir->pme_order);
1244 gmx_fio_do_real(fio, ir->ewald_rtol);
1246 if (file_version >= 93)
1248 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1252 ir->ewald_rtol_lj = ir->ewald_rtol;
1254 gmx_fio_do_int(fio, ir->ewald_geometry);
1255 gmx_fio_do_real(fio, ir->epsilon_surface);
1257 /* ignore bOptFFT */
1258 if (file_version < tpxv_RemoveObsoleteParameters1)
1260 gmx_fio_do_gmx_bool(fio, bdum);
1263 if (file_version >= 93)
1265 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1267 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1268 gmx_fio_do_int(fio, ir->etc);
1269 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1270 * but the values 0 and 1 still mean no and
1271 * berendsen temperature coupling, respectively.
1273 if (file_version >= 79)
1275 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1277 if (file_version >= 71)
1279 gmx_fio_do_int(fio, ir->nsttcouple);
1283 ir->nsttcouple = ir->nstcalcenergy;
1285 gmx_fio_do_int(fio, ir->epc);
1286 gmx_fio_do_int(fio, ir->epct);
1287 if (file_version >= 71)
1289 gmx_fio_do_int(fio, ir->nstpcouple);
1293 ir->nstpcouple = ir->nstcalcenergy;
1295 gmx_fio_do_real(fio, ir->tau_p);
1296 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1297 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1298 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1299 gmx_fio_do_rvec(fio, ir->compress[XX]);
1300 gmx_fio_do_rvec(fio, ir->compress[YY]);
1301 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1302 gmx_fio_do_int(fio, ir->refcoord_scaling);
1303 gmx_fio_do_rvec(fio, ir->posres_com);
1304 gmx_fio_do_rvec(fio, ir->posres_comB);
1306 if (file_version < 79)
1308 gmx_fio_do_int(fio, ir->andersen_seed);
1312 ir->andersen_seed = 0;
1315 gmx_fio_do_real(fio, ir->shake_tol);
1317 gmx_fio_do_int(fio, ir->efep);
1318 do_fepvals(fio, ir->fepvals, bRead, file_version);
1320 if (file_version >= 79)
1322 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1325 ir->bSimTemp = TRUE;
1330 ir->bSimTemp = FALSE;
1334 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1337 if (file_version >= 79)
1339 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1342 ir->bExpanded = TRUE;
1346 ir->bExpanded = FALSE;
1351 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1354 gmx_fio_do_int(fio, ir->eDisre);
1355 gmx_fio_do_int(fio, ir->eDisreWeighting);
1356 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1357 gmx_fio_do_real(fio, ir->dr_fc);
1358 gmx_fio_do_real(fio, ir->dr_tau);
1359 gmx_fio_do_int(fio, ir->nstdisreout);
1360 gmx_fio_do_real(fio, ir->orires_fc);
1361 gmx_fio_do_real(fio, ir->orires_tau);
1362 gmx_fio_do_int(fio, ir->nstorireout);
1364 /* ignore dihre_fc */
1365 if (file_version < 79)
1367 gmx_fio_do_real(fio, rdum);
1370 gmx_fio_do_real(fio, ir->em_stepsize);
1371 gmx_fio_do_real(fio, ir->em_tol);
1372 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1373 gmx_fio_do_int(fio, ir->niter);
1374 gmx_fio_do_real(fio, ir->fc_stepsize);
1375 gmx_fio_do_int(fio, ir->eConstrAlg);
1376 gmx_fio_do_int(fio, ir->nProjOrder);
1377 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1378 gmx_fio_do_int(fio, ir->nLincsIter);
1379 gmx_fio_do_real(fio, ir->bd_fric);
1380 if (file_version >= tpxv_Use64BitRandomSeed)
1382 gmx_fio_do_int64(fio, ir->ld_seed);
1386 gmx_fio_do_int(fio, idum);
1390 for (i = 0; i < DIM; i++)
1392 gmx_fio_do_rvec(fio, ir->deform[i]);
1394 gmx_fio_do_real(fio, ir->cos_accel);
1396 gmx_fio_do_int(fio, ir->userint1);
1397 gmx_fio_do_int(fio, ir->userint2);
1398 gmx_fio_do_int(fio, ir->userint3);
1399 gmx_fio_do_int(fio, ir->userint4);
1400 gmx_fio_do_real(fio, ir->userreal1);
1401 gmx_fio_do_real(fio, ir->userreal2);
1402 gmx_fio_do_real(fio, ir->userreal3);
1403 gmx_fio_do_real(fio, ir->userreal4);
1405 /* AdResS is removed, but we need to be able to read old files,
1406 and let mdrun refuse to run them */
1407 if (file_version >= 77 && file_version < tpxv_RemoveAdress)
1409 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1412 int idum, numThermoForceGroups, numEnergyGroups;
1415 gmx_fio_do_int(fio, idum);
1416 gmx_fio_do_real(fio, rdum);
1417 gmx_fio_do_real(fio, rdum);
1418 gmx_fio_do_real(fio, rdum);
1419 gmx_fio_do_int(fio, idum);
1420 gmx_fio_do_int(fio, idum);
1421 gmx_fio_do_rvec(fio, rvecdum);
1422 gmx_fio_do_int(fio, numThermoForceGroups);
1423 gmx_fio_do_real(fio, rdum);
1424 gmx_fio_do_int(fio, numEnergyGroups);
1425 gmx_fio_do_int(fio, idum);
1427 if (numThermoForceGroups > 0)
1429 std::vector<int> idumn(numThermoForceGroups);
1430 gmx_fio_ndo_int(fio, idumn.data(), idumn.size());
1432 if (numEnergyGroups > 0)
1434 std::vector<int> idumn(numEnergyGroups);
1435 gmx_fio_ndo_int(fio, idumn.data(), idumn.size());
1441 ir->bAdress = FALSE;
1448 if (file_version >= tpxv_PullCoordTypeGeom)
1450 gmx_fio_do_gmx_bool(fio, ir->bPull);
1454 gmx_fio_do_int(fio, ePullOld);
1455 ir->bPull = (ePullOld > 0);
1456 /* We removed the first ePull=ePullNo for the enum */
1465 do_pull(fio, ir->pull, bRead, file_version, ePullOld);
1469 if (file_version >= tpxv_AcceleratedWeightHistogram)
1471 gmx_fio_do_gmx_bool(fio, ir->bDoAwh);
1477 snew(ir->awhParams, 1);
1479 do_awh(fio, ir->awhParams, bRead, file_version);
1487 /* Enforced rotation */
1488 if (file_version >= 74)
1490 gmx_fio_do_gmx_bool(fio, ir->bRot);
1497 do_rot(fio, ir->rot, bRead);
1505 /* Interactive molecular dynamics */
1506 if (file_version >= tpxv_InteractiveMolecularDynamics)
1508 gmx_fio_do_gmx_bool(fio, ir->bIMD);
1515 do_imd(fio, ir->imd, bRead);
1520 /* We don't support IMD sessions for old .tpr files */
1525 gmx_fio_do_int(fio, ir->opts.ngtc);
1526 if (file_version >= 69)
1528 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1532 ir->opts.nhchainlength = 1;
1534 gmx_fio_do_int(fio, ir->opts.ngacc);
1535 gmx_fio_do_int(fio, ir->opts.ngfrz);
1536 gmx_fio_do_int(fio, ir->opts.ngener);
1540 snew(ir->opts.nrdf, ir->opts.ngtc);
1541 snew(ir->opts.ref_t, ir->opts.ngtc);
1542 snew(ir->opts.annealing, ir->opts.ngtc);
1543 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1544 snew(ir->opts.anneal_time, ir->opts.ngtc);
1545 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1546 snew(ir->opts.tau_t, ir->opts.ngtc);
1547 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1548 snew(ir->opts.acc, ir->opts.ngacc);
1549 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1551 if (ir->opts.ngtc > 0)
1553 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1554 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1555 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1557 if (ir->opts.ngfrz > 0)
1559 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1561 if (ir->opts.ngacc > 0)
1563 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1565 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1566 ir->opts.ngener*ir->opts.ngener);
1568 /* First read the lists with annealing and npoints for each group */
1569 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1570 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1571 for (j = 0; j < (ir->opts.ngtc); j++)
1573 k = ir->opts.anneal_npoints[j];
1576 snew(ir->opts.anneal_time[j], k);
1577 snew(ir->opts.anneal_temp[j], k);
1579 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1580 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1584 gmx_fio_do_int(fio, ir->nwall);
1585 gmx_fio_do_int(fio, ir->wall_type);
1586 gmx_fio_do_real(fio, ir->wall_r_linpot);
1587 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1588 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1589 gmx_fio_do_real(fio, ir->wall_density[0]);
1590 gmx_fio_do_real(fio, ir->wall_density[1]);
1591 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1594 /* Cosine stuff for electric fields */
1595 if (file_version < tpxv_GenericParamsForElectricField)
1597 do_legacy_efield(fio, ¶msObj);
1601 if (file_version >= tpxv_ComputationalElectrophysiology)
1603 gmx_fio_do_int(fio, ir->eSwapCoords);
1604 if (ir->eSwapCoords != eswapNO)
1610 do_swapcoords_tpx(fio, ir->swap, bRead, file_version);
1616 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1617 gmx_fio_do_int(fio, ir->QMMMscheme);
1618 gmx_fio_do_real(fio, ir->scalefactor);
1619 gmx_fio_do_int(fio, ir->opts.ngQM);
1622 snew(ir->opts.QMmethod, ir->opts.ngQM);
1623 snew(ir->opts.QMbasis, ir->opts.ngQM);
1624 snew(ir->opts.QMcharge, ir->opts.ngQM);
1625 snew(ir->opts.QMmult, ir->opts.ngQM);
1626 snew(ir->opts.bSH, ir->opts.ngQM);
1627 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1628 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1629 snew(ir->opts.SAon, ir->opts.ngQM);
1630 snew(ir->opts.SAoff, ir->opts.ngQM);
1631 snew(ir->opts.SAsteps, ir->opts.ngQM);
1633 if (ir->opts.ngQM > 0 && ir->bQMMM)
1635 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1636 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1637 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1638 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1639 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1640 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1641 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1642 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1643 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1644 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1645 /* We leave in dummy i/o for removed parameters to avoid
1646 * changing the tpr format for every QMMM change.
1648 std::vector<int> dummy(ir->opts.ngQM, 0);
1649 gmx_fio_ndo_int(fio, dummy.data(), ir->opts.ngQM);
1650 gmx_fio_ndo_int(fio, dummy.data(), ir->opts.ngQM);
1652 /* end of QMMM stuff */
1655 if (file_version >= tpxv_GenericParamsForElectricField)
1657 gmx::FileIOXdrSerializer serializer(fio);
1660 paramsObj.mergeObject(
1661 gmx::deserializeKeyValueTree(&serializer));
1665 GMX_RELEASE_ASSERT(ir->params != nullptr,
1666 "Parameters should be present when writing inputrec");
1667 gmx::serializeKeyValueTree(*ir->params, &serializer);
1672 ir->params = new gmx::KeyValueTreeObject(paramsBuilder.build());
1677 static void do_harm(t_fileio *fio, t_iparams *iparams)
1679 gmx_fio_do_real(fio, iparams->harmonic.rA);
1680 gmx_fio_do_real(fio, iparams->harmonic.krA);
1681 gmx_fio_do_real(fio, iparams->harmonic.rB);
1682 gmx_fio_do_real(fio, iparams->harmonic.krB);
1685 static void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1686 gmx_bool bRead, int file_version)
1699 do_harm(fio, iparams);
1700 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1702 /* Correct incorrect storage of parameters */
1703 iparams->pdihs.phiB = iparams->pdihs.phiA;
1704 iparams->pdihs.cpB = iparams->pdihs.cpA;
1708 gmx_fio_do_real(fio, iparams->harmonic.rA);
1709 gmx_fio_do_real(fio, iparams->harmonic.krA);
1711 case F_LINEAR_ANGLES:
1712 gmx_fio_do_real(fio, iparams->linangle.klinA);
1713 gmx_fio_do_real(fio, iparams->linangle.aA);
1714 gmx_fio_do_real(fio, iparams->linangle.klinB);
1715 gmx_fio_do_real(fio, iparams->linangle.aB);
1718 gmx_fio_do_real(fio, iparams->fene.bm);
1719 gmx_fio_do_real(fio, iparams->fene.kb);
1723 gmx_fio_do_real(fio, iparams->restraint.lowA);
1724 gmx_fio_do_real(fio, iparams->restraint.up1A);
1725 gmx_fio_do_real(fio, iparams->restraint.up2A);
1726 gmx_fio_do_real(fio, iparams->restraint.kA);
1727 gmx_fio_do_real(fio, iparams->restraint.lowB);
1728 gmx_fio_do_real(fio, iparams->restraint.up1B);
1729 gmx_fio_do_real(fio, iparams->restraint.up2B);
1730 gmx_fio_do_real(fio, iparams->restraint.kB);
1736 gmx_fio_do_real(fio, iparams->tab.kA);
1737 gmx_fio_do_int(fio, iparams->tab.table);
1738 gmx_fio_do_real(fio, iparams->tab.kB);
1740 case F_CROSS_BOND_BONDS:
1741 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1742 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1743 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1745 case F_CROSS_BOND_ANGLES:
1746 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1747 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1748 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1749 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1751 case F_UREY_BRADLEY:
1752 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1753 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1754 gmx_fio_do_real(fio, iparams->u_b.r13A);
1755 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1756 if (file_version >= 79)
1758 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1759 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1760 gmx_fio_do_real(fio, iparams->u_b.r13B);
1761 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1765 iparams->u_b.thetaB = iparams->u_b.thetaA;
1766 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1767 iparams->u_b.r13B = iparams->u_b.r13A;
1768 iparams->u_b.kUBB = iparams->u_b.kUBA;
1771 case F_QUARTIC_ANGLES:
1772 gmx_fio_do_real(fio, iparams->qangle.theta);
1773 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1776 gmx_fio_do_real(fio, iparams->bham.a);
1777 gmx_fio_do_real(fio, iparams->bham.b);
1778 gmx_fio_do_real(fio, iparams->bham.c);
1781 gmx_fio_do_real(fio, iparams->morse.b0A);
1782 gmx_fio_do_real(fio, iparams->morse.cbA);
1783 gmx_fio_do_real(fio, iparams->morse.betaA);
1784 if (file_version >= 79)
1786 gmx_fio_do_real(fio, iparams->morse.b0B);
1787 gmx_fio_do_real(fio, iparams->morse.cbB);
1788 gmx_fio_do_real(fio, iparams->morse.betaB);
1792 iparams->morse.b0B = iparams->morse.b0A;
1793 iparams->morse.cbB = iparams->morse.cbA;
1794 iparams->morse.betaB = iparams->morse.betaA;
1798 gmx_fio_do_real(fio, iparams->cubic.b0);
1799 gmx_fio_do_real(fio, iparams->cubic.kb);
1800 gmx_fio_do_real(fio, iparams->cubic.kcub);
1804 case F_POLARIZATION:
1805 gmx_fio_do_real(fio, iparams->polarize.alpha);
1808 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1809 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1810 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1813 gmx_fio_do_real(fio, iparams->wpol.al_x);
1814 gmx_fio_do_real(fio, iparams->wpol.al_y);
1815 gmx_fio_do_real(fio, iparams->wpol.al_z);
1816 gmx_fio_do_real(fio, iparams->wpol.rOH);
1817 gmx_fio_do_real(fio, iparams->wpol.rHH);
1818 gmx_fio_do_real(fio, iparams->wpol.rOD);
1821 gmx_fio_do_real(fio, iparams->thole.a);
1822 gmx_fio_do_real(fio, iparams->thole.alpha1);
1823 gmx_fio_do_real(fio, iparams->thole.alpha2);
1824 gmx_fio_do_real(fio, iparams->thole.rfac);
1827 gmx_fio_do_real(fio, iparams->lj.c6);
1828 gmx_fio_do_real(fio, iparams->lj.c12);
1831 gmx_fio_do_real(fio, iparams->lj14.c6A);
1832 gmx_fio_do_real(fio, iparams->lj14.c12A);
1833 gmx_fio_do_real(fio, iparams->lj14.c6B);
1834 gmx_fio_do_real(fio, iparams->lj14.c12B);
1837 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1838 gmx_fio_do_real(fio, iparams->ljc14.qi);
1839 gmx_fio_do_real(fio, iparams->ljc14.qj);
1840 gmx_fio_do_real(fio, iparams->ljc14.c6);
1841 gmx_fio_do_real(fio, iparams->ljc14.c12);
1843 case F_LJC_PAIRS_NB:
1844 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1845 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1846 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1847 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1853 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1854 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1855 gmx_fio_do_real(fio, iparams->pdihs.phiB);
1856 gmx_fio_do_real(fio, iparams->pdihs.cpB);
1857 gmx_fio_do_int(fio, iparams->pdihs.mult);
1860 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1861 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1864 gmx_fio_do_int(fio, iparams->disres.label);
1865 gmx_fio_do_int(fio, iparams->disres.type);
1866 gmx_fio_do_real(fio, iparams->disres.low);
1867 gmx_fio_do_real(fio, iparams->disres.up1);
1868 gmx_fio_do_real(fio, iparams->disres.up2);
1869 gmx_fio_do_real(fio, iparams->disres.kfac);
1872 gmx_fio_do_int(fio, iparams->orires.ex);
1873 gmx_fio_do_int(fio, iparams->orires.label);
1874 gmx_fio_do_int(fio, iparams->orires.power);
1875 gmx_fio_do_real(fio, iparams->orires.c);
1876 gmx_fio_do_real(fio, iparams->orires.obs);
1877 gmx_fio_do_real(fio, iparams->orires.kfac);
1880 if (file_version < 82)
1882 gmx_fio_do_int(fio, idum);
1883 gmx_fio_do_int(fio, idum);
1885 gmx_fio_do_real(fio, iparams->dihres.phiA);
1886 gmx_fio_do_real(fio, iparams->dihres.dphiA);
1887 gmx_fio_do_real(fio, iparams->dihres.kfacA);
1888 if (file_version >= 82)
1890 gmx_fio_do_real(fio, iparams->dihres.phiB);
1891 gmx_fio_do_real(fio, iparams->dihres.dphiB);
1892 gmx_fio_do_real(fio, iparams->dihres.kfacB);
1896 iparams->dihres.phiB = iparams->dihres.phiA;
1897 iparams->dihres.dphiB = iparams->dihres.dphiA;
1898 iparams->dihres.kfacB = iparams->dihres.kfacA;
1902 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
1903 gmx_fio_do_rvec(fio, iparams->posres.fcA);
1904 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
1905 gmx_fio_do_rvec(fio, iparams->posres.fcB);
1908 gmx_fio_do_int(fio, iparams->fbposres.geom);
1909 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
1910 gmx_fio_do_real(fio, iparams->fbposres.r);
1911 gmx_fio_do_real(fio, iparams->fbposres.k);
1914 gmx_fio_ndo_real(fio, iparams->cbtdihs.cbtcA, NR_CBTDIHS);
1917 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1918 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1921 /* Fourier dihedrals are internally represented
1922 * as Ryckaert-Bellemans since those are faster to compute.
1924 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1925 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1929 gmx_fio_do_real(fio, iparams->constr.dA);
1930 gmx_fio_do_real(fio, iparams->constr.dB);
1933 gmx_fio_do_real(fio, iparams->settle.doh);
1934 gmx_fio_do_real(fio, iparams->settle.dhh);
1937 gmx_fio_do_real(fio, iparams->vsite.a);
1942 gmx_fio_do_real(fio, iparams->vsite.a);
1943 gmx_fio_do_real(fio, iparams->vsite.b);
1948 gmx_fio_do_real(fio, iparams->vsite.a);
1949 gmx_fio_do_real(fio, iparams->vsite.b);
1950 gmx_fio_do_real(fio, iparams->vsite.c);
1953 gmx_fio_do_int(fio, iparams->vsiten.n);
1954 gmx_fio_do_real(fio, iparams->vsiten.a);
1956 case F_GB12_NOLONGERUSED:
1957 case F_GB13_NOLONGERUSED:
1958 case F_GB14_NOLONGERUSED:
1959 // Implicit solvent parameters can still be read, but never used
1962 if (file_version < 68)
1964 gmx_fio_do_real(fio, rdum);
1965 gmx_fio_do_real(fio, rdum);
1966 gmx_fio_do_real(fio, rdum);
1967 gmx_fio_do_real(fio, rdum);
1969 if (file_version < tpxv_RemoveImplicitSolvation)
1971 gmx_fio_do_real(fio, rdum);
1972 gmx_fio_do_real(fio, rdum);
1973 gmx_fio_do_real(fio, rdum);
1974 gmx_fio_do_real(fio, rdum);
1975 gmx_fio_do_real(fio, rdum);
1980 gmx_fio_do_int(fio, iparams->cmap.cmapA);
1981 gmx_fio_do_int(fio, iparams->cmap.cmapB);
1984 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1985 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1989 static void do_ilist(t_fileio *fio, InteractionList *ilist, gmx_bool bRead)
1991 int nr = ilist->size();
1992 gmx_fio_do_int(fio, nr);
1995 ilist->iatoms.resize(nr);
1997 gmx_fio_ndo_int(fio, ilist->iatoms.data(), ilist->size());
2000 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2001 gmx_bool bRead, int file_version)
2003 gmx_fio_do_int(fio, ffparams->atnr);
2004 int numTypes = ffparams->numTypes();
2005 gmx_fio_do_int(fio, numTypes);
2008 ffparams->functype.resize(numTypes);
2009 ffparams->iparams.resize(numTypes);
2011 /* Read/write all the function types */
2012 gmx_fio_ndo_int(fio, ffparams->functype.data(), ffparams->functype.size());
2014 if (file_version >= 66)
2016 gmx_fio_do_double(fio, ffparams->reppow);
2020 ffparams->reppow = 12.0;
2023 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2025 /* Check whether all these function types are supported by the code.
2026 * In practice the code is backwards compatible, which means that the
2027 * numbering may have to be altered from old numbering to new numbering
2029 for (int i = 0; i < ffparams->numTypes(); i++)
2033 /* Loop over file versions */
2034 for (int k = 0; k < NFTUPD; k++)
2036 /* Compare the read file_version to the update table */
2037 if ((file_version < ftupd[k].fvnr) &&
2038 (ffparams->functype[i] >= ftupd[k].ftype))
2040 ffparams->functype[i] += 1;
2045 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2050 static void add_settle_atoms(InteractionList *ilist)
2054 /* Settle used to only store the first atom: add the other two */
2055 ilist->iatoms.resize(2*ilist->size());
2056 for (i = ilist->size()/4 - 1; i >= 0; i--)
2058 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2059 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2060 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2061 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2065 static void do_ilists(t_fileio *fio, InteractionLists *ilists, gmx_bool bRead,
2068 GMX_RELEASE_ASSERT(ilists, "Need a valid ilists object");
2069 GMX_RELEASE_ASSERT(ilists->size() == F_NRE, "The code needs to be in sync with InteractionLists");
2071 for (int j = 0; j < F_NRE; j++)
2073 InteractionList &ilist = (*ilists)[j];
2074 gmx_bool bClear = FALSE;
2077 for (int k = 0; k < NFTUPD; k++)
2079 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2087 ilist.iatoms.clear();
2091 do_ilist(fio, &ilist, bRead);
2092 if (file_version < 78 && j == F_SETTLE && ilist.size() > 0)
2094 add_settle_atoms(&ilist);
2100 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead)
2102 gmx_fio_do_int(fio, block->nr);
2105 if ((block->nalloc_index > 0) && (nullptr != block->index))
2107 sfree(block->index);
2109 block->nalloc_index = block->nr+1;
2110 snew(block->index, block->nalloc_index);
2112 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2115 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead)
2117 gmx_fio_do_int(fio, block->nr);
2118 gmx_fio_do_int(fio, block->nra);
2121 block->nalloc_index = block->nr+1;
2122 snew(block->index, block->nalloc_index);
2123 block->nalloc_a = block->nra;
2124 snew(block->a, block->nalloc_a);
2126 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2127 gmx_fio_ndo_int(fio, block->a, block->nra);
2130 /* This is a primitive routine to make it possible to translate atomic numbers
2131 * to element names when reading TPR files, without making the Gromacs library
2132 * directory a dependency on mdrun (which is the case if we need elements.dat).
2135 atomicnumber_to_element(int atomicnumber)
2139 /* This does not have to be complete, so we only include elements likely
2140 * to occur in PDB files.
2142 switch (atomicnumber)
2144 case 1: p = "H"; break;
2145 case 5: p = "B"; break;
2146 case 6: p = "C"; break;
2147 case 7: p = "N"; break;
2148 case 8: p = "O"; break;
2149 case 9: p = "F"; break;
2150 case 11: p = "Na"; break;
2151 case 12: p = "Mg"; break;
2152 case 15: p = "P"; break;
2153 case 16: p = "S"; break;
2154 case 17: p = "Cl"; break;
2155 case 18: p = "Ar"; break;
2156 case 19: p = "K"; break;
2157 case 20: p = "Ca"; break;
2158 case 25: p = "Mn"; break;
2159 case 26: p = "Fe"; break;
2160 case 28: p = "Ni"; break;
2161 case 29: p = "Cu"; break;
2162 case 30: p = "Zn"; break;
2163 case 35: p = "Br"; break;
2164 case 47: p = "Ag"; break;
2165 default: p = ""; break;
2171 static void do_atom(t_fileio *fio, t_atom *atom, gmx_bool bRead)
2173 gmx_fio_do_real(fio, atom->m);
2174 gmx_fio_do_real(fio, atom->q);
2175 gmx_fio_do_real(fio, atom->mB);
2176 gmx_fio_do_real(fio, atom->qB);
2177 gmx_fio_do_ushort(fio, atom->type);
2178 gmx_fio_do_ushort(fio, atom->typeB);
2179 gmx_fio_do_int(fio, atom->ptype);
2180 gmx_fio_do_int(fio, atom->resind);
2181 gmx_fio_do_int(fio, atom->atomnumber);
2184 /* Set element string from atomic number if present.
2185 * This routine returns an empty string if the name is not found.
2187 std::strncpy(atom->elem, atomicnumber_to_element(atom->atomnumber), 4);
2188 /* avoid warnings about potentially unterminated string */
2189 atom->elem[3] = '\0';
2193 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead)
2195 for (int j = 0; j < ngrp; j++)
2197 gmx_fio_do_int(fio, grps[j].nr);
2200 snew(grps[j].nm_ind, grps[j].nr);
2202 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2206 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2212 gmx_fio_do_int(fio, ls);
2213 *nm = get_symtab_handle(symtab, ls);
2217 ls = lookup_symtab(symtab, *nm);
2218 gmx_fio_do_int(fio, ls);
2222 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2227 for (j = 0; (j < nstr); j++)
2229 do_symstr(fio, &(nm[j]), bRead, symtab);
2233 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2234 t_symtab *symtab, int file_version)
2238 for (j = 0; (j < n); j++)
2240 do_symstr(fio, &(ri[j].name), bRead, symtab);
2241 if (file_version >= 63)
2243 gmx_fio_do_int(fio, ri[j].nr);
2244 gmx_fio_do_uchar(fio, ri[j].ic);
2254 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2259 gmx_fio_do_int(fio, atoms->nr);
2260 gmx_fio_do_int(fio, atoms->nres);
2263 /* Since we have always written all t_atom properties in the tpr file
2264 * (at least for all backward compatible versions), we don't store
2265 * but simple set the booleans here.
2267 atoms->haveMass = TRUE;
2268 atoms->haveCharge = TRUE;
2269 atoms->haveType = TRUE;
2270 atoms->haveBState = TRUE;
2271 atoms->havePdbInfo = FALSE;
2273 snew(atoms->atom, atoms->nr);
2274 snew(atoms->atomname, atoms->nr);
2275 snew(atoms->atomtype, atoms->nr);
2276 snew(atoms->atomtypeB, atoms->nr);
2277 snew(atoms->resinfo, atoms->nres);
2278 atoms->pdbinfo = nullptr;
2282 GMX_RELEASE_ASSERT(atoms->haveMass && atoms->haveCharge && atoms->haveType && atoms->haveBState, "Mass, charge, atomtype and B-state parameters should be present in t_atoms when writing a tpr file");
2284 for (i = 0; (i < atoms->nr); i++)
2286 do_atom(fio, &atoms->atom[i], bRead);
2288 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2289 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2290 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2292 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2295 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2296 gmx_bool bRead, t_symtab *symtab)
2300 do_grps(fio, egcNR, groups->grps, bRead);
2301 gmx_fio_do_int(fio, groups->ngrpname);
2304 snew(groups->grpname, groups->ngrpname);
2306 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2307 for (g = 0; g < egcNR; g++)
2309 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2310 if (groups->ngrpnr[g] == 0)
2314 groups->grpnr[g] = nullptr;
2321 snew(groups->grpnr[g], groups->ngrpnr[g]);
2323 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2328 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2333 gmx_fio_do_int(fio, atomtypes->nr);
2337 snew(atomtypes->atomnumber, j);
2339 if (bRead && file_version < tpxv_RemoveImplicitSolvation)
2341 std::vector<real> dummy(atomtypes->nr, 0);
2342 gmx_fio_ndo_real(fio, dummy.data(), dummy.size());
2343 gmx_fio_ndo_real(fio, dummy.data(), dummy.size());
2344 gmx_fio_ndo_real(fio, dummy.data(), dummy.size());
2346 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2348 if (bRead && file_version >= 60 && file_version < tpxv_RemoveImplicitSolvation)
2350 std::vector<real> dummy(atomtypes->nr, 0);
2351 gmx_fio_ndo_real(fio, dummy.data(), dummy.size());
2352 gmx_fio_ndo_real(fio, dummy.data(), dummy.size());
2356 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2362 gmx_fio_do_int(fio, symtab->nr);
2366 snew(symtab->symbuf, 1);
2367 symbuf = symtab->symbuf;
2368 symbuf->bufsize = nr;
2369 snew(symbuf->buf, nr);
2370 for (i = 0; (i < nr); i++)
2372 gmx_fio_do_string(fio, buf);
2373 symbuf->buf[i] = gmx_strdup(buf);
2378 symbuf = symtab->symbuf;
2379 while (symbuf != nullptr)
2381 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2383 gmx_fio_do_string(fio, symbuf->buf[i]);
2386 symbuf = symbuf->next;
2390 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2395 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2398 int ngrid = cmap_grid->cmapdata.size();
2399 gmx_fio_do_int(fio, ngrid);
2400 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2402 int gs = cmap_grid->grid_spacing;
2403 int nelem = gs * gs;
2407 cmap_grid->cmapdata.resize(ngrid);
2409 for (int i = 0; i < ngrid; i++)
2411 cmap_grid->cmapdata[i].cmap.resize(4*nelem);
2415 for (int i = 0; i < ngrid; i++)
2417 for (int j = 0; j < nelem; j++)
2419 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2420 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2421 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2422 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2428 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2429 t_symtab *symtab, int file_version)
2431 do_symstr(fio, &(molt->name), bRead, symtab);
2433 do_atoms(fio, &molt->atoms, bRead, symtab, file_version);
2435 do_ilists(fio, &molt->ilist, bRead, file_version);
2437 do_block(fio, &molt->cgs, bRead);
2439 /* This used to be in the atoms struct */
2440 do_blocka(fio, &molt->excls, bRead);
2443 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb,
2444 int numAtomsPerMolecule,
2447 gmx_fio_do_int(fio, molb->type);
2448 gmx_fio_do_int(fio, molb->nmol);
2449 /* To maintain forward topology reading compatibility, we store #atoms.
2450 * TODO: Change this to conditional reading of a dummy int when we
2451 * increase tpx_generation.
2453 gmx_fio_do_int(fio, numAtomsPerMolecule);
2454 /* Position restraint coordinates */
2455 int numPosres_xA = molb->posres_xA.size();
2456 gmx_fio_do_int(fio, numPosres_xA);
2457 if (numPosres_xA > 0)
2461 molb->posres_xA.resize(numPosres_xA);
2463 gmx_fio_ndo_rvec(fio, as_rvec_array(molb->posres_xA.data()), numPosres_xA);
2465 int numPosres_xB = molb->posres_xB.size();
2466 gmx_fio_do_int(fio, numPosres_xB);
2467 if (numPosres_xB > 0)
2471 molb->posres_xB.resize(numPosres_xB);
2473 gmx_fio_ndo_rvec(fio, as_rvec_array(molb->posres_xB.data()), numPosres_xB);
2478 static void set_disres_npair(gmx_mtop_t *mtop)
2480 gmx_mtop_ilistloop_t iloop;
2483 gmx::ArrayRef<t_iparams> ip = mtop->ffparams.iparams;
2485 iloop = gmx_mtop_ilistloop_init(mtop);
2486 while (const InteractionLists *ilist = gmx_mtop_ilistloop_next(iloop, &nmol))
2488 const InteractionList &il = (*ilist)[F_DISRES];
2492 gmx::ArrayRef<const int> a = il.iatoms;
2494 for (int i = 0; i < il.size(); i += 3)
2497 if (i+3 == il.size() || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2499 ip[a[i]].disres.npair = npair;
2507 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2510 do_symtab(fio, &(mtop->symtab), bRead);
2512 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2514 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
2516 int nmoltype = mtop->moltype.size();
2517 gmx_fio_do_int(fio, nmoltype);
2520 mtop->moltype.resize(nmoltype);
2522 for (gmx_moltype_t &moltype : mtop->moltype)
2524 do_moltype(fio, &moltype, bRead, &mtop->symtab, file_version);
2527 int nmolblock = mtop->molblock.size();
2528 gmx_fio_do_int(fio, nmolblock);
2531 mtop->molblock.resize(nmolblock);
2533 for (gmx_molblock_t &molblock : mtop->molblock)
2535 int numAtomsPerMolecule = (bRead ? 0 : mtop->moltype[molblock.type].atoms.nr);
2536 do_molblock(fio, &molblock, numAtomsPerMolecule, bRead);
2538 gmx_fio_do_int(fio, mtop->natoms);
2540 if (file_version >= tpxv_IntermolecularBondeds)
2542 gmx_fio_do_gmx_bool(fio, mtop->bIntermolecularInteractions);
2543 if (mtop->bIntermolecularInteractions)
2547 mtop->intermolecular_ilist = gmx::compat::make_unique<InteractionLists>();
2549 do_ilists(fio, mtop->intermolecular_ilist.get(), bRead, file_version);
2554 mtop->bIntermolecularInteractions = FALSE;
2557 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
2559 if (file_version >= 65)
2561 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
2565 mtop->ffparams.cmap_grid.grid_spacing = 0;
2566 mtop->ffparams.cmap_grid.cmapdata.clear();
2569 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab));
2571 mtop->haveMoleculeIndices = true;
2575 close_symtab(&(mtop->symtab));
2579 /* If TopOnlyOK is TRUE then we can read even future versions
2580 * of tpx files, provided the fileGeneration hasn't changed.
2581 * If it is FALSE, we need the inputrecord too, and bail out
2582 * if the file is newer than the program.
2584 * The version and generation of the topology (see top of this file)
2585 * are returned in the two last arguments, if those arguments are non-NULL.
2587 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2589 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
2590 gmx_bool TopOnlyOK, int *fileVersionPointer, int *fileGenerationPointer)
2593 char file_tag[STRLEN];
2598 int fileVersion; /* Version number of the code that wrote the file */
2599 int fileGeneration; /* Generation version number of the code that wrote the file */
2601 /* XDR binary topology file */
2602 precision = sizeof(real);
2605 gmx_fio_do_string(fio, buf);
2606 if (std::strncmp(buf, "VERSION", 7) != 0)
2608 gmx_fatal(FARGS, "Can not read file %s,\n"
2609 " this file is from a GROMACS version which is older than 2.0\n"
2610 " Make a new one with grompp or use a gro or pdb file, if possible",
2611 gmx_fio_getname(fio));
2613 gmx_fio_do_int(fio, precision);
2614 bDouble = (precision == sizeof(double));
2615 if ((precision != sizeof(float)) && !bDouble)
2617 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
2618 "instead of %zu or %zu",
2619 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
2621 gmx_fio_setprecision(fio, bDouble);
2622 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
2623 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
2627 snprintf(buf, STRLEN, "VERSION %s", gmx_version());
2628 gmx_fio_write_string(fio, buf);
2629 bDouble = (precision == sizeof(double));
2630 gmx_fio_setprecision(fio, bDouble);
2631 gmx_fio_do_int(fio, precision);
2632 fileVersion = tpx_version;
2633 sprintf(file_tag, "%s", tpx_tag);
2634 fileGeneration = tpx_generation;
2637 /* Check versions! */
2638 gmx_fio_do_int(fio, fileVersion);
2640 /* This is for backward compatibility with development versions 77-79
2641 * where the tag was, mistakenly, placed before the generation,
2642 * which would cause a segv instead of a proper error message
2643 * when reading the topology only from tpx with <77 code.
2645 if (fileVersion >= 77 && fileVersion <= 79)
2647 gmx_fio_do_string(fio, file_tag);
2650 gmx_fio_do_int(fio, fileGeneration);
2652 if (fileVersion >= 81)
2654 gmx_fio_do_string(fio, file_tag);
2658 if (fileVersion < 77)
2660 /* Versions before 77 don't have the tag, set it to release */
2661 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
2664 if (std::strcmp(file_tag, tpx_tag) != 0)
2666 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
2669 /* We only support reading tpx files with the same tag as the code
2670 * or tpx files with the release tag and with lower version number.
2672 if (std::strcmp(file_tag, TPX_TAG_RELEASE) != 0 && fileVersion < tpx_version)
2674 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2675 gmx_fio_getname(fio), fileVersion, file_tag,
2676 tpx_version, tpx_tag);
2681 if (fileVersionPointer)
2683 *fileVersionPointer = fileVersion;
2685 if (fileGenerationPointer)
2687 *fileGenerationPointer = fileGeneration;
2690 if ((fileVersion <= tpx_incompatible_version) ||
2691 ((fileVersion > tpx_version) && !TopOnlyOK) ||
2692 (fileGeneration > tpx_generation) ||
2693 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
2695 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
2696 gmx_fio_getname(fio), fileVersion, tpx_version);
2699 gmx_fio_do_int(fio, tpx->natoms);
2700 gmx_fio_do_int(fio, tpx->ngtc);
2702 if (fileVersion < 62)
2704 gmx_fio_do_int(fio, idum);
2705 gmx_fio_do_real(fio, rdum);
2707 if (fileVersion >= 79)
2709 gmx_fio_do_int(fio, tpx->fep_state);
2711 gmx_fio_do_real(fio, tpx->lambda);
2712 gmx_fio_do_gmx_bool(fio, tpx->bIr);
2713 gmx_fio_do_gmx_bool(fio, tpx->bTop);
2714 gmx_fio_do_gmx_bool(fio, tpx->bX);
2715 gmx_fio_do_gmx_bool(fio, tpx->bV);
2716 gmx_fio_do_gmx_bool(fio, tpx->bF);
2717 gmx_fio_do_gmx_bool(fio, tpx->bBox);
2719 if ((fileGeneration > tpx_generation))
2721 /* This can only happen if TopOnlyOK=TRUE */
2726 static int do_tpx(t_fileio *fio, gmx_bool bRead,
2727 t_inputrec *ir, t_state *state, rvec *x, rvec *v,
2733 gmx_bool bPeriodicMols;
2737 GMX_RELEASE_ASSERT(state != nullptr, "Cannot write a null state");
2738 GMX_RELEASE_ASSERT(x == nullptr && v == nullptr, "Passing separate x and v pointers to do_tpx() is not supported when writing");
2740 tpx.natoms = state->natoms;
2741 tpx.ngtc = state->ngtc;
2742 tpx.fep_state = state->fep_state;
2743 tpx.lambda = state->lambda[efptFEP];
2744 tpx.bIr = ir != nullptr;
2745 tpx.bTop = mtop != nullptr;
2746 tpx.bX = (state->flags & (1 << estX)) != 0;
2747 tpx.bV = (state->flags & (1 << estV)) != 0;
2753 GMX_RELEASE_ASSERT(!(x == nullptr && v != nullptr), "Passing x==NULL and v!=NULL is not supported");
2756 TopOnlyOK = (ir == nullptr);
2758 int fileVersion; /* Version number of the code that wrote the file */
2759 int fileGeneration; /* Generation version number of the code that wrote the file */
2760 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &fileVersion, &fileGeneration);
2765 init_gtc_state(state, tpx.ngtc, 0, 0);
2768 // v is also nullptr by the above assertion, so we may
2769 // need to make memory in state for storing the contents
2773 state->flags |= (1 << estX);
2777 state->flags |= (1 << estV);
2779 state_change_natoms(state, tpx.natoms);
2785 x = state->x.rvec_array();
2786 v = state->v.rvec_array();
2789 #define do_test(fio, b, p) if (bRead && ((p) != NULL) && !(b)) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
2791 do_test(fio, tpx.bBox, state->box);
2794 gmx_fio_ndo_rvec(fio, state->box, DIM);
2795 if (fileVersion >= 51)
2797 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
2801 /* We initialize box_rel after reading the inputrec */
2802 clear_mat(state->box_rel);
2804 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
2805 if (fileVersion < 56)
2808 gmx_fio_ndo_rvec(fio, mdum, DIM);
2812 if (state->ngtc > 0)
2815 snew(dumv, state->ngtc);
2816 if (fileVersion < 69)
2818 gmx_fio_ndo_real(fio, dumv, state->ngtc);
2820 /* These used to be the Berendsen tcoupl_lambda's */
2821 gmx_fio_ndo_real(fio, dumv, state->ngtc);
2825 do_test(fio, tpx.bTop, mtop);
2830 do_mtop(fio, mtop, bRead, fileVersion);
2835 do_mtop(fio, &dum_top, bRead, fileVersion);
2838 do_test(fio, tpx.bX, x);
2843 state->flags |= (1<<estX);
2845 gmx_fio_ndo_rvec(fio, x, tpx.natoms);
2848 do_test(fio, tpx.bV, v);
2853 state->flags |= (1<<estV);
2855 gmx_fio_ndo_rvec(fio, v, tpx.natoms);
2858 // No need to run do_test when the last argument is NULL
2862 snew(dummyForces, state->natoms);
2863 gmx_fio_ndo_rvec(fio, dummyForces, tpx.natoms);
2867 /* Starting with tpx version 26, we have the inputrec
2868 * at the end of the file, so we can ignore it
2869 * if the file is never than the software (but still the
2870 * same generation - see comments at the top of this file.
2875 bPeriodicMols = FALSE;
2877 do_test(fio, tpx.bIr, ir);
2880 if (fileVersion >= 53)
2882 /* Removed the pbc info from do_inputrec, since we always want it */
2886 bPeriodicMols = ir->bPeriodicMols;
2888 gmx_fio_do_int(fio, ePBC);
2889 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
2891 if (fileGeneration <= tpx_generation && ir)
2893 do_inputrec(fio, ir, bRead, fileVersion);
2894 if (fileVersion < 51)
2896 set_box_rel(ir, state);
2898 if (fileVersion < 53)
2901 bPeriodicMols = ir->bPeriodicMols;
2904 if (bRead && ir && fileVersion >= 53)
2906 /* We need to do this after do_inputrec, since that initializes ir */
2908 ir->bPeriodicMols = bPeriodicMols;
2916 if (state->ngtc == 0)
2918 /* Reading old version without tcoupl state data: set it */
2919 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
2921 if (tpx.bTop && mtop)
2923 if (fileVersion < 57)
2925 if (mtop->moltype[0].ilist[F_DISRES].size() > 0)
2927 ir->eDisre = edrSimple;
2931 ir->eDisre = edrNone;
2934 set_disres_npair(mtop);
2938 if (tpx.bTop && mtop)
2940 gmx_mtop_finalize(mtop);
2947 static t_fileio *open_tpx(const char *fn, const char *mode)
2949 return gmx_fio_open(fn, mode);
2952 static void close_tpx(t_fileio *fio)
2957 /************************************************************
2959 * The following routines are the exported ones
2961 ************************************************************/
2963 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK)
2967 fio = open_tpx(fn, "r");
2968 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, nullptr, nullptr);
2972 void write_tpx_state(const char *fn,
2973 const t_inputrec *ir, const t_state *state, const gmx_mtop_t *mtop)
2977 fio = open_tpx(fn, "w");
2979 const_cast<t_inputrec *>(ir),
2980 const_cast<t_state *>(state), nullptr, nullptr,
2981 const_cast<gmx_mtop_t *>(mtop));
2985 void read_tpx_state(const char *fn,
2986 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
2990 fio = open_tpx(fn, "r");
2991 do_tpx(fio, TRUE, ir, state, nullptr, nullptr, mtop);
2995 int read_tpx(const char *fn,
2996 t_inputrec *ir, matrix box, int *natoms,
2997 rvec *x, rvec *v, gmx_mtop_t *mtop)
3003 fio = open_tpx(fn, "r");
3004 ePBC = do_tpx(fio, TRUE, ir, &state, x, v, mtop);
3006 if (mtop != nullptr && natoms != nullptr)
3008 *natoms = mtop->natoms;
3012 copy_mat(state.box, box);
3018 int read_tpx_top(const char *fn,
3019 t_inputrec *ir, matrix box, int *natoms,
3020 rvec *x, rvec *v, t_topology *top)
3025 ePBC = read_tpx(fn, ir, box, natoms, x, v, &mtop);
3027 *top = gmx_mtop_t_to_t_topology(&mtop, true);
3032 gmx_bool fn2bTPX(const char *file)
3034 return (efTPR == fn2ftp(file));
3037 void pr_tpxheader(FILE *fp, int indent, const char *title, const t_tpxheader *sh)
3039 if (available(fp, sh, indent, title))
3041 indent = pr_title(fp, indent, title);
3042 pr_indent(fp, indent);
3043 fprintf(fp, "bIr = %spresent\n", sh->bIr ? "" : "not ");
3044 pr_indent(fp, indent);
3045 fprintf(fp, "bBox = %spresent\n", sh->bBox ? "" : "not ");
3046 pr_indent(fp, indent);
3047 fprintf(fp, "bTop = %spresent\n", sh->bTop ? "" : "not ");
3048 pr_indent(fp, indent);
3049 fprintf(fp, "bX = %spresent\n", sh->bX ? "" : "not ");
3050 pr_indent(fp, indent);
3051 fprintf(fp, "bV = %spresent\n", sh->bV ? "" : "not ");
3052 pr_indent(fp, indent);
3053 fprintf(fp, "bF = %spresent\n", sh->bF ? "" : "not ");
3055 pr_indent(fp, indent);
3056 fprintf(fp, "natoms = %d\n", sh->natoms);
3057 pr_indent(fp, indent);
3058 fprintf(fp, "lambda = %e\n", sh->lambda);