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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 /* This file is completely threadsafe - keep it that way! */
42 #include "gromacs/fileio/tpxio.h"
52 #include "gromacs/applied_forces/awh/read_params.h"
53 #include "gromacs/fileio/filetypes.h"
54 #include "gromacs/fileio/gmxfio.h"
55 #include "gromacs/fileio/gmxfio_xdr.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdtypes/awh_history.h"
59 #include "gromacs/mdtypes/awh_params.h"
60 #include "gromacs/mdtypes/inputrec.h"
61 #include "gromacs/mdtypes/md_enums.h"
62 #include "gromacs/mdtypes/multipletimestepping.h"
63 #include "gromacs/mdtypes/pull_params.h"
64 #include "gromacs/mdtypes/state.h"
65 #include "gromacs/pbcutil/boxutilities.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/topology/block.h"
68 #include "gromacs/topology/ifunc.h"
69 #include "gromacs/topology/mtop_util.h"
70 #include "gromacs/topology/symtab.h"
71 #include "gromacs/topology/topology.h"
72 #include "gromacs/utility/arraysize.h"
73 #include "gromacs/utility/baseversion.h"
74 #include "gromacs/utility/cstringutil.h"
75 #include "gromacs/utility/fatalerror.h"
76 #include "gromacs/utility/futil.h"
77 #include "gromacs/utility/gmxassert.h"
78 #include "gromacs/utility/inmemoryserializer.h"
79 #include "gromacs/utility/iserializer.h"
80 #include "gromacs/utility/keyvaluetreebuilder.h"
81 #include "gromacs/utility/keyvaluetreeserializer.h"
82 #include "gromacs/utility/smalloc.h"
83 #include "gromacs/utility/snprintf.h"
84 #include "gromacs/utility/txtdump.h"
86 #define TPX_TAG_RELEASE "release"
88 /*! \brief Tag string for the file format written to run input files
89 * written by this version of the code.
91 * Change this if you want to change the run input format in a feature
92 * branch. This ensures that there will not be different run input
93 * formats around which cannot be distinguished, while not causing
94 * problems rebasing the feature branch onto upstream changes. When
95 * merging with mainstream GROMACS, set this tag string back to
96 * TPX_TAG_RELEASE, and instead add an element to tpxv.
98 static const std::string tpx_tag = TPX_TAG_RELEASE;
100 /*! \brief Enum of values that describe the contents of a tpr file
101 * whose format matches a version number
103 * The enum helps the code be more self-documenting and ensure merges
104 * do not silently resolve when two patches make the same bump. When
105 * adding new functionality, add a new element just above tpxv_Count
106 * in this enumeration, and write code below that does the right thing
107 * according to the value of file_version.
111 tpxv_ComputationalElectrophysiology =
112 96, /**< support for ion/water position swaps (computational electrophysiology) */
113 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to int64_t */
114 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
115 tpxv_InteractiveMolecularDynamics, /**< interactive molecular dynamics (IMD) */
116 tpxv_RemoveObsoleteParameters1, /**< remove optimize_fft, dihre_fc, nstcheckpoint */
117 tpxv_PullCoordTypeGeom, /**< add pull type and geometry per group and flat-bottom */
118 tpxv_PullGeomDirRel, /**< add pull geometry direction-relative */
119 tpxv_IntermolecularBondeds, /**< permit inter-molecular bonded interactions in the topology */
120 tpxv_CompElWithSwapLayerOffset, /**< added parameters for improved CompEl setups */
121 tpxv_CompElPolyatomicIonsAndMultipleIonTypes, /**< CompEl now can handle polyatomic ions and more than two types of ions */
122 tpxv_RemoveAdress, /**< removed support for AdResS */
123 tpxv_PullCoordNGroup, /**< add ngroup to pull coord */
124 tpxv_RemoveTwinRange, /**< removed support for twin-range interactions */
125 tpxv_ReplacePullPrintCOM12, /**< Replaced print-com-1, 2 with pull-print-com */
126 tpxv_PullExternalPotential, /**< Added pull type external potential */
127 tpxv_GenericParamsForElectricField, /**< Introduced KeyValueTree and moved electric field parameters */
128 tpxv_AcceleratedWeightHistogram, /**< sampling with accelerated weight histogram method (AWH) */
129 tpxv_RemoveImplicitSolvation, /**< removed support for implicit solvation */
130 tpxv_PullPrevStepCOMAsReference, /**< Enabled using the COM of the pull group of the last frame as reference for PBC */
131 tpxv_MimicQMMM, /**< Introduced support for MiMiC QM/MM interface */
132 tpxv_PullAverage, /**< Added possibility to output average pull force and position */
133 tpxv_GenericInternalParameters, /**< Added internal parameters for mdrun modules*/
134 tpxv_VSite2FD, /**< Added 2FD type virtual site */
135 tpxv_AddSizeField, /**< Added field with information about the size of the serialized tpr file in bytes, excluding the header */
136 tpxv_StoreNonBondedInteractionExclusionGroup, /**< Store the non bonded interaction exclusion group in the topology */
137 tpxv_VSite1, /**< Added 1 type virtual site */
138 tpxv_MTS, /**< Added multiple time stepping */
139 tpxv_RemovedConstantAcceleration, /**< Removed support for constant acceleration NEMD. */
140 tpxv_TransformationPullCoord, /**< Support for transformation pull coordinates */
141 tpxv_SoftcoreGapsys, /**< Added gapsys softcore function */
142 tpxv_ReaddedConstantAcceleration, /**< Re-added support for constant acceleration NEMD. */
143 tpxv_RemoveTholeRfac, /**< Remove unused rfac parameter from thole listed force */
144 tpxv_Count /**< the total number of tpxv versions */
147 /*! \brief Version number of the file format written to run input
148 * files by this version of the code.
150 * The tpx_version increases whenever the file format in the main
151 * development branch changes, due to an extension of the tpxv enum above.
152 * Backward compatibility for reading old run input files is maintained
153 * by checking this version number against that of the file and then using
154 * the correct code path.
156 * When developing a feature branch that needs to change the run input
157 * file format, change tpx_tag instead. */
158 static const int tpx_version = tpxv_Count - 1;
162 * Enum keeping track of incompatible changes for older TPR versions.
164 * The enum should be updated with a new field when editing the TOPOLOGY
165 * or HEADER of the tpx format. In particular, updating ftupd or
166 * changing the fields of TprHeaderVersion often trigger such needs.
168 * This way we can maintain forward compatibility too for all analysis tools
169 * and/or external programs that only need to know the atom/residue names,
170 * charges, and bond connectivity.
172 * It first appeared in tpx version 26, when I also moved the inputrecord
173 * to the end of the tpx file, so we can just skip it if we only
176 * In particular, it must be increased when adding new elements to
177 * ftupd, so that old code can read new .tpr files.
179 enum class TpxGeneration : int
181 Initial = 26, //! First version is 26
182 AddSizeField, //! TPR header modified for writing as a block.
183 AddVSite1, //! ftupd changed to include VSite1 type.
184 Count //! Number of entries.
187 //! Value of Current TPR generation.
188 static const int tpx_generation = static_cast<int>(TpxGeneration::Count) - 1;
190 /* This number should be the most recent backwards incompatible version
191 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
193 static const int tpx_incompatible_version = 57; // GMX4.0 has version 58
196 /* Struct used to maintain tpx compatibility when function types are added */
199 int fvnr; /* file version number in which the function type first appeared */
200 int ftype; /* function type */
204 * TODO The following three lines make little sense, please clarify if
205 * you've had to work out how ftupd works.
207 * The entries should be ordered in:
208 * 1. ascending function type number
209 * 2. ascending file version number
211 * Because we support reading of old .tpr file versions (even when
212 * mdrun can no longer run the simulation), we need to be able to read
213 * obsolete t_interaction_function types. Any data read from such
214 * fields is discarded. Their names have _NOLONGERUSED appended to
215 * them to make things clear.
217 * When adding to or making breaking changes to reading this struct,
218 * update TpxGeneration.
220 static const t_ftupd ftupd[] = {
221 { 70, F_RESTRBONDS },
222 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
223 { 76, F_LINEAR_ANGLES },
224 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
225 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
227 { 60, F_GB12_NOLONGERUSED },
228 { 61, F_GB13_NOLONGERUSED },
229 { 61, F_GB14_NOLONGERUSED },
230 { 72, F_GBPOL_NOLONGERUSED },
231 { 72, F_NPSOLVATION_NOLONGERUSED },
233 { 76, F_ANHARM_POL },
235 { tpxv_VSite1, F_VSITE1 },
236 { tpxv_VSite2FD, F_VSITE2FD },
237 { tpxv_GenericInternalParameters, F_DENSITYFITTING },
238 { 69, F_VTEMP_NOLONGERUSED },
249 { 79, F_DVDL_RESTRAINT },
250 { 79, F_DVDL_TEMPERATURE },
252 #define NFTUPD asize(ftupd)
254 /**************************************************************
256 * Now the higer level routines that do io of the structures and arrays
258 **************************************************************/
259 static void do_pullgrp_tpx_pre95(gmx::ISerializer* serializer, t_pull_group* pgrp, t_pull_coord* pcrd)
263 int numAtoms = pgrp->ind.size();
264 serializer->doInt(&numAtoms);
265 pgrp->ind.resize(numAtoms);
266 serializer->doIntArray(pgrp->ind.data(), numAtoms);
267 int numWeights = pgrp->weight.size();
268 serializer->doInt(&numWeights);
269 pgrp->weight.resize(numWeights);
270 serializer->doRealArray(pgrp->weight.data(), numWeights);
271 serializer->doInt(&pgrp->pbcatom);
272 serializer->doRvec(&pcrd->vec.as_vec());
273 clear_rvec(pcrd->origin);
274 serializer->doRvec(&tmp);
276 serializer->doReal(&pcrd->rate);
277 serializer->doReal(&pcrd->k);
278 serializer->doReal(&pcrd->kB);
281 static void do_pull_group(gmx::ISerializer* serializer, t_pull_group* pgrp)
283 int numAtoms = pgrp->ind.size();
284 serializer->doInt(&numAtoms);
285 pgrp->ind.resize(numAtoms);
286 serializer->doIntArray(pgrp->ind.data(), numAtoms);
287 int numWeights = pgrp->weight.size();
288 serializer->doInt(&numWeights);
289 pgrp->weight.resize(numWeights);
290 serializer->doRealArray(pgrp->weight.data(), numWeights);
291 serializer->doInt(&pgrp->pbcatom);
294 static void do_pull_coord(gmx::ISerializer* serializer,
297 PullingAlgorithm ePullOld,
298 PullGroupGeometry eGeomOld,
301 if (file_version >= tpxv_PullCoordNGroup)
303 serializer->doEnumAsInt(&pcrd->eType);
304 if (file_version >= tpxv_PullExternalPotential)
306 if (pcrd->eType == PullingAlgorithm::External)
308 serializer->doString(&pcrd->externalPotentialProvider);
312 pcrd->externalPotentialProvider.clear();
317 if (serializer->reading())
319 pcrd->externalPotentialProvider.clear();
322 /* Note that we try to support adding new geometries without
323 * changing the tpx version. This requires checks when printing the
324 * geometry string and a check and fatal_error in init_pull.
326 serializer->doEnumAsInt(&pcrd->eGeom);
327 serializer->doInt(&pcrd->ngroup);
328 if (pcrd->ngroup <= c_pullCoordNgroupMax)
330 serializer->doIntArray(pcrd->group.data(), pcrd->ngroup);
334 /* More groups in file than supported, this must be a new geometry
335 * that is not supported by our current code. Since we will not
336 * use the groups for this coord (checks in the pull and WHAM code
337 * ensure this), we can ignore the groups and set ngroup=0.
340 snew(dum, pcrd->ngroup);
341 serializer->doIntArray(dum, pcrd->ngroup);
346 serializer->doIvec(&pcrd->dim.as_vec());
347 if (file_version >= tpxv_TransformationPullCoord)
349 serializer->doString(&pcrd->expression);
353 if (serializer->reading())
355 pcrd->expression.clear();
362 serializer->doInt(&pcrd->group[0]);
363 serializer->doInt(&pcrd->group[1]);
364 if (file_version >= tpxv_PullCoordTypeGeom)
366 pcrd->ngroup = (pcrd->eGeom == PullGroupGeometry::DirectionRelative ? 4 : 2);
367 serializer->doEnumAsInt(&pcrd->eType);
368 serializer->doEnumAsInt(&pcrd->eGeom);
369 if (pcrd->ngroup == 4)
371 serializer->doInt(&pcrd->group[2]);
372 serializer->doInt(&pcrd->group[3]);
374 serializer->doIvec(&pcrd->dim.as_vec());
378 pcrd->eType = ePullOld;
379 pcrd->eGeom = eGeomOld;
380 copy_ivec(dimOld, pcrd->dim);
383 serializer->doRvec(&pcrd->origin.as_vec());
384 serializer->doRvec(&pcrd->vec.as_vec());
385 if (file_version >= tpxv_PullCoordTypeGeom)
387 serializer->doBool(&pcrd->bStart);
391 /* This parameter is only printed, but not actually used by mdrun */
392 pcrd->bStart = FALSE;
394 serializer->doReal(&pcrd->init);
395 serializer->doReal(&pcrd->rate);
396 serializer->doReal(&pcrd->k);
397 serializer->doReal(&pcrd->kB);
400 static void do_expandedvals(gmx::ISerializer* serializer, t_expanded* expand, t_lambda* fepvals, int file_version)
402 int n_lambda = fepvals->n_lambda;
404 /* reset the lambda calculation window */
405 fepvals->lambda_start_n = 0;
406 fepvals->lambda_stop_n = n_lambda;
407 if (file_version >= 79)
411 expand->init_lambda_weights.resize(n_lambda);
412 serializer->doRealArray(expand->init_lambda_weights.data(), n_lambda);
413 serializer->doBool(&expand->bInit_weights);
416 serializer->doInt(&expand->nstexpanded);
417 serializer->doEnumAsInt(&expand->elmcmove);
418 serializer->doEnumAsInt(&expand->elamstats);
419 serializer->doInt(&expand->lmc_repeats);
420 serializer->doInt(&expand->gibbsdeltalam);
421 serializer->doInt(&expand->lmc_forced_nstart);
422 serializer->doInt(&expand->lmc_seed);
423 serializer->doReal(&expand->mc_temp);
424 serializer->doBool(&expand->bSymmetrizedTMatrix);
425 serializer->doInt(&expand->nstTij);
426 serializer->doInt(&expand->minvarmin);
427 serializer->doInt(&expand->c_range);
428 serializer->doReal(&expand->wl_scale);
429 serializer->doReal(&expand->wl_ratio);
430 serializer->doReal(&expand->init_wl_delta);
431 serializer->doBool(&expand->bWLoneovert);
432 serializer->doEnumAsInt(&expand->elmceq);
433 serializer->doInt(&expand->equil_steps);
434 serializer->doInt(&expand->equil_samples);
435 serializer->doInt(&expand->equil_n_at_lam);
436 serializer->doReal(&expand->equil_wl_delta);
437 serializer->doReal(&expand->equil_ratio);
441 static void do_simtempvals(gmx::ISerializer* serializer, t_simtemp* simtemp, int n_lambda, int file_version)
443 if (file_version >= 79)
445 serializer->doEnumAsInt(&simtemp->eSimTempScale);
446 serializer->doReal(&simtemp->simtemp_high);
447 serializer->doReal(&simtemp->simtemp_low);
450 if (serializer->reading())
452 simtemp->temperatures.resize(n_lambda);
454 serializer->doRealArray(simtemp->temperatures.data(), n_lambda);
459 static void do_imd(gmx::ISerializer* serializer, t_IMD* imd)
461 serializer->doInt(&imd->nat);
462 if (serializer->reading())
464 snew(imd->ind, imd->nat);
466 serializer->doIntArray(imd->ind, imd->nat);
469 static void do_fepvals(gmx::ISerializer* serializer, t_lambda* fepvals, int file_version)
471 /* i is defined in the ndo_double macro; use g to iterate. */
474 /* free energy values */
476 if (file_version >= 79)
478 serializer->doInt(&fepvals->init_fep_state);
479 serializer->doDouble(&fepvals->init_lambda);
480 serializer->doDouble(&fepvals->delta_lambda);
482 else if (file_version >= 59)
484 serializer->doDouble(&fepvals->init_lambda);
485 serializer->doDouble(&fepvals->delta_lambda);
489 serializer->doReal(&rdum);
490 fepvals->init_lambda = rdum;
491 serializer->doReal(&rdum);
492 fepvals->delta_lambda = rdum;
494 if (file_version >= 79)
496 serializer->doInt(&fepvals->n_lambda);
497 for (auto g : keysOf(fepvals->all_lambda))
499 if (fepvals->n_lambda > 0)
501 fepvals->all_lambda[g].resize(fepvals->n_lambda);
502 serializer->doDoubleArray(fepvals->all_lambda[g].data(), fepvals->n_lambda);
503 serializer->doBoolArray(fepvals->separate_dvdl.begin(), fepvals->separate_dvdl.size());
505 else if (fepvals->init_lambda >= 0)
507 fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Fep] = TRUE;
511 else if (file_version >= 64)
513 serializer->doInt(&fepvals->n_lambda);
514 if (serializer->reading())
516 /* still allocate the all_lambda array's contents. */
517 for (auto g : keysOf(fepvals->all_lambda))
519 if (fepvals->n_lambda > 0)
521 fepvals->all_lambda[g].resize(fepvals->n_lambda);
525 serializer->doDoubleArray(fepvals->all_lambda[FreeEnergyPerturbationCouplingType::Fep].data(),
527 if (fepvals->init_lambda >= 0)
529 fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Fep] = TRUE;
531 if (serializer->reading())
533 /* copy the contents of the efptFEP lambda component to all
534 the other components */
535 for (auto g : keysOf(fepvals->all_lambda))
537 for (int h = 0; h < fepvals->n_lambda; h++)
539 if (g != FreeEnergyPerturbationCouplingType::Fep)
541 fepvals->all_lambda[g][h] =
542 fepvals->all_lambda[FreeEnergyPerturbationCouplingType::Fep][h];
551 fepvals->n_lambda = 0;
552 if (fepvals->init_lambda >= 0)
554 fepvals->separate_dvdl[FreeEnergyPerturbationCouplingType::Fep] = TRUE;
557 serializer->doReal(&fepvals->sc_alpha);
558 serializer->doInt(&fepvals->sc_power);
559 if (file_version >= 79)
561 serializer->doReal(&fepvals->sc_r_power);
565 fepvals->sc_r_power = 6.0;
567 if (fepvals->sc_r_power != 6.0)
569 gmx_fatal(FARGS, "sc-r-power=48 is no longer supported");
571 serializer->doReal(&fepvals->sc_sigma);
572 if (serializer->reading())
574 if (file_version >= 71)
576 fepvals->sc_sigma_min = fepvals->sc_sigma;
580 fepvals->sc_sigma_min = 0;
583 if (file_version >= 79)
585 serializer->doBool(&fepvals->bScCoul);
589 fepvals->bScCoul = TRUE;
591 if (file_version >= 64)
593 serializer->doInt(&fepvals->nstdhdl);
597 fepvals->nstdhdl = 1;
600 if (file_version >= 73)
602 serializer->doEnumAsInt(&fepvals->separate_dhdl_file);
603 serializer->doEnumAsInt(&fepvals->dhdl_derivatives);
607 fepvals->separate_dhdl_file = SeparateDhdlFile::Yes;
608 fepvals->dhdl_derivatives = DhDlDerivativeCalculation::Yes;
610 if (file_version >= 71)
612 serializer->doInt(&fepvals->dh_hist_size);
613 serializer->doDouble(&fepvals->dh_hist_spacing);
617 fepvals->dh_hist_size = 0;
618 fepvals->dh_hist_spacing = 0.1;
620 if (file_version >= 79)
622 serializer->doEnumAsInt(&fepvals->edHdLPrintEnergy);
626 fepvals->edHdLPrintEnergy = FreeEnergyPrintEnergy::No;
628 if (file_version >= tpxv_SoftcoreGapsys)
630 serializer->doInt(reinterpret_cast<int*>(&fepvals->softcoreFunction));
631 serializer->doReal(&fepvals->scGapsysScaleLinpointLJ);
632 serializer->doReal(&fepvals->scGapsysScaleLinpointQ);
633 serializer->doReal(&fepvals->scGapsysSigmaLJ);
637 fepvals->softcoreFunction = SoftcoreType::Beutler;
638 fepvals->scGapsysScaleLinpointLJ = 0.85;
639 fepvals->scGapsysScaleLinpointQ = 0.3;
640 fepvals->scGapsysSigmaLJ = 0.3;
643 /* handle lambda_neighbors */
644 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
646 serializer->doInt(&fepvals->lambda_neighbors);
647 if ((fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0)
648 && (fepvals->init_lambda < 0))
650 fepvals->lambda_start_n = (fepvals->init_fep_state - fepvals->lambda_neighbors);
651 fepvals->lambda_stop_n = (fepvals->init_fep_state + fepvals->lambda_neighbors + 1);
652 if (fepvals->lambda_start_n < 0)
654 fepvals->lambda_start_n = 0;
656 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
658 fepvals->lambda_stop_n = fepvals->n_lambda;
663 fepvals->lambda_start_n = 0;
664 fepvals->lambda_stop_n = fepvals->n_lambda;
669 fepvals->lambda_start_n = 0;
670 fepvals->lambda_stop_n = fepvals->n_lambda;
674 static void do_pull(gmx::ISerializer* serializer, pull_params_t* pull, int file_version, PullingAlgorithm ePullOld)
676 PullGroupGeometry eGeomOld = PullGroupGeometry::Count;
680 if (file_version >= 95)
682 serializer->doInt(&pull->ngroup);
684 serializer->doInt(&pull->ncoord);
685 if (file_version < 95)
687 pull->ngroup = pull->ncoord + 1;
689 if (file_version < tpxv_PullCoordTypeGeom)
693 serializer->doEnumAsInt(&eGeomOld);
694 serializer->doIvec(&dimOld);
695 /* The inner cylinder radius, now removed */
696 serializer->doReal(&dum);
698 serializer->doReal(&pull->cylinder_r);
699 serializer->doReal(&pull->constr_tol);
700 if (file_version >= 95)
702 serializer->doBool(&pull->bPrintCOM);
703 /* With file_version < 95 this value is set below */
705 if (file_version >= tpxv_ReplacePullPrintCOM12)
707 serializer->doBool(&pull->bPrintRefValue);
708 serializer->doBool(&pull->bPrintComp);
710 else if (file_version >= tpxv_PullCoordTypeGeom)
713 serializer->doInt(&idum); /* used to be bPrintCOM2 */
714 serializer->doBool(&pull->bPrintRefValue);
715 serializer->doBool(&pull->bPrintComp);
719 pull->bPrintRefValue = FALSE;
720 pull->bPrintComp = TRUE;
722 serializer->doInt(&pull->nstxout);
723 serializer->doInt(&pull->nstfout);
724 if (file_version >= tpxv_PullPrevStepCOMAsReference)
726 serializer->doBool(&pull->bSetPbcRefToPrevStepCOM);
730 pull->bSetPbcRefToPrevStepCOM = FALSE;
732 pull->group.resize(pull->ngroup);
733 pull->coord.resize(pull->ncoord);
734 if (file_version < 95)
736 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
737 if (eGeomOld == PullGroupGeometry::DirectionPBC)
739 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
741 if (eGeomOld > PullGroupGeometry::DirectionPBC)
745 case (PullGroupGeometry::DirectionRelative):
746 eGeomOld = PullGroupGeometry::DirectionPBC;
748 case (PullGroupGeometry::Angle):
749 eGeomOld = PullGroupGeometry::DirectionRelative;
751 case (PullGroupGeometry::Dihedral): eGeomOld = PullGroupGeometry::Angle; break;
752 case (PullGroupGeometry::AngleAxis): eGeomOld = PullGroupGeometry::Dihedral; break;
753 case (PullGroupGeometry::Count): eGeomOld = PullGroupGeometry::AngleAxis; break;
754 default: GMX_RELEASE_ASSERT(false, "Unhandled old pull type");
758 for (g = 0; g < pull->ngroup; g++)
760 /* We read and ignore a pull coordinate for group 0 */
761 do_pullgrp_tpx_pre95(serializer, &pull->group[g], &pull->coord[std::max(g - 1, 0)]);
764 pull->coord[g - 1].group[0] = 0;
765 pull->coord[g - 1].group[1] = g;
769 pull->bPrintCOM = (!pull->group[0].ind.empty());
773 for (g = 0; g < pull->ngroup; g++)
775 do_pull_group(serializer, &pull->group[g]);
777 for (g = 0; g < pull->ncoord; g++)
779 do_pull_coord(serializer, &pull->coord[g], file_version, ePullOld, eGeomOld, dimOld);
780 if (serializer->reading())
782 pull->coord[g].coordIndex = g;
786 if (file_version >= tpxv_PullAverage)
790 v = pull->bXOutAverage;
791 serializer->doBool(&v);
792 pull->bXOutAverage = v;
793 v = pull->bFOutAverage;
794 serializer->doBool(&v);
795 pull->bFOutAverage = v;
800 static void do_rotgrp(gmx::ISerializer* serializer, t_rotgrp* rotg)
802 serializer->doEnumAsInt(&rotg->eType);
803 if (serializer->reading())
806 serializer->doInt(&temp);
807 rotg->bMassW = static_cast<bool>(temp);
811 int temp = static_cast<int>(rotg->bMassW);
812 serializer->doInt(&temp);
814 serializer->doInt(&rotg->nat);
815 if (serializer->reading())
817 snew(rotg->ind, rotg->nat);
819 serializer->doIntArray(rotg->ind, rotg->nat);
820 if (serializer->reading())
822 snew(rotg->x_ref, rotg->nat);
824 serializer->doRvecArray(rotg->x_ref, rotg->nat);
825 serializer->doRvec(&rotg->inputVec);
826 serializer->doRvec(&rotg->pivot);
827 serializer->doReal(&rotg->rate);
828 serializer->doReal(&rotg->k);
829 serializer->doReal(&rotg->slab_dist);
830 serializer->doReal(&rotg->min_gaussian);
831 serializer->doReal(&rotg->eps);
832 serializer->doEnumAsInt(&rotg->eFittype);
833 serializer->doInt(&rotg->PotAngle_nstep);
834 serializer->doReal(&rotg->PotAngle_step);
837 static void do_rot(gmx::ISerializer* serializer, t_rot* rot)
841 serializer->doInt(&rot->ngrp);
842 serializer->doInt(&rot->nstrout);
843 serializer->doInt(&rot->nstsout);
844 if (serializer->reading())
846 snew(rot->grp, rot->ngrp);
848 for (g = 0; g < rot->ngrp; g++)
850 do_rotgrp(serializer, &rot->grp[g]);
855 static void do_swapgroup(gmx::ISerializer* serializer, t_swapGroup* g)
858 /* Name of the group or molecule */
860 if (serializer->reading())
862 serializer->doString(&buf);
863 g->molname = gmx_strdup(buf.c_str());
868 serializer->doString(&buf);
871 /* Number of atoms in the group */
872 serializer->doInt(&g->nat);
874 /* The group's atom indices */
875 if (serializer->reading())
877 snew(g->ind, g->nat);
879 serializer->doIntArray(g->ind, g->nat);
881 /* Requested counts for compartments A and B */
882 serializer->doIntArray(g->nmolReq.data(), static_cast<int>(Compartment::Count));
885 static void do_swapcoords_tpx(gmx::ISerializer* serializer, t_swapcoords* swap, int file_version)
887 /* Enums for better readability of the code */
900 if (file_version >= tpxv_CompElPolyatomicIonsAndMultipleIonTypes)
902 /* The total number of swap groups is the sum of the fixed groups
903 * (split0, split1, solvent), and the user-defined groups (2+ types of ions)
905 serializer->doInt(&swap->ngrp);
906 if (serializer->reading())
908 snew(swap->grp, swap->ngrp);
910 for (int ig = 0; ig < swap->ngrp; ig++)
912 do_swapgroup(serializer, &swap->grp[ig]);
914 serializer->doBool(&swap->massw_split[eChannel0]);
915 serializer->doBool(&swap->massw_split[eChannel1]);
916 serializer->doInt(&swap->nstswap);
917 serializer->doInt(&swap->nAverage);
918 serializer->doReal(&swap->threshold);
919 serializer->doReal(&swap->cyl0r);
920 serializer->doReal(&swap->cyl0u);
921 serializer->doReal(&swap->cyl0l);
922 serializer->doReal(&swap->cyl1r);
923 serializer->doReal(&swap->cyl1u);
924 serializer->doReal(&swap->cyl1l);
928 /*** Support reading older CompEl .tpr files ***/
930 /* In the original CompEl .tpr files, we always have 5 groups: */
932 snew(swap->grp, swap->ngrp);
934 swap->grp[static_cast<int>(SwapGroupSplittingType::Split0)].molname = gmx_strdup("split0"); // group 0: split0
935 swap->grp[static_cast<int>(SwapGroupSplittingType::Split1)].molname = gmx_strdup("split1"); // group 1: split1
936 swap->grp[static_cast<int>(SwapGroupSplittingType::Solvent)].molname =
937 gmx_strdup("solvent"); // group 2: solvent
938 swap->grp[3].molname = gmx_strdup("anions"); // group 3: anions
939 swap->grp[4].molname = gmx_strdup("cations"); // group 4: cations
941 serializer->doInt(&swap->grp[3].nat);
942 serializer->doInt(&swap->grp[static_cast<int>(SwapGroupSplittingType::Solvent)].nat);
943 serializer->doInt(&swap->grp[static_cast<int>(SwapGroupSplittingType::Split0)].nat);
944 serializer->doBool(&swap->massw_split[eChannel0]);
945 serializer->doInt(&swap->grp[static_cast<int>(SwapGroupSplittingType::Split1)].nat);
946 serializer->doBool(&swap->massw_split[eChannel1]);
947 serializer->doInt(&swap->nstswap);
948 serializer->doInt(&swap->nAverage);
949 serializer->doReal(&swap->threshold);
950 serializer->doReal(&swap->cyl0r);
951 serializer->doReal(&swap->cyl0u);
952 serializer->doReal(&swap->cyl0l);
953 serializer->doReal(&swap->cyl1r);
954 serializer->doReal(&swap->cyl1u);
955 serializer->doReal(&swap->cyl1l);
957 // The order[] array keeps compatibility with older .tpr files
958 // by reading in the groups in the classic order
960 const int order[4] = { 3,
961 static_cast<int>(SwapGroupSplittingType::Solvent),
962 static_cast<int>(SwapGroupSplittingType::Split0),
963 static_cast<int>(SwapGroupSplittingType::Split1) };
965 for (int ig = 0; ig < 4; ig++)
968 snew(swap->grp[g].ind, swap->grp[g].nat);
969 serializer->doIntArray(swap->grp[g].ind, swap->grp[g].nat);
973 for (int j = eCompA; j <= eCompB; j++)
975 serializer->doInt(&swap->grp[3].nmolReq[j]); // group 3 = anions
976 serializer->doInt(&swap->grp[4].nmolReq[j]); // group 4 = cations
978 } /* End support reading older CompEl .tpr files */
980 if (file_version >= tpxv_CompElWithSwapLayerOffset)
982 serializer->doReal(&swap->bulkOffset[eCompA]);
983 serializer->doReal(&swap->bulkOffset[eCompB]);
987 static void do_legacy_efield(gmx::ISerializer* serializer, 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 serializer->doInt(&n);
1000 serializer->doInt(&nt);
1001 std::vector<real> aa(n + 1), phi(nt + 1), at(nt + 1), phit(nt + 1);
1002 serializer->doRealArray(aa.data(), n);
1003 serializer->doRealArray(phi.data(), n);
1004 serializer->doRealArray(at.data(), nt);
1005 serializer->doRealArray(phit.data(), nt);
1008 if (n > 1 || nt > 1)
1011 "Can not handle tpr files with more than one electric field term per "
1014 auto dimObj = efieldObj.addObject(dimName[j]);
1015 dimObj.addValue<real>("E0", aa[0]);
1016 dimObj.addValue<real>("omega", at[0]);
1017 dimObj.addValue<real>("t0", phi[0]);
1018 dimObj.addValue<real>("sigma", phit[0]);
1024 static void do_inputrec(gmx::ISerializer* serializer, t_inputrec* ir, int file_version)
1026 int i, j, k, idum = 0;
1028 gmx_bool bdum = false;
1030 if (file_version != tpx_version)
1032 /* Give a warning about features that are not accessible */
1033 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n", file_version, tpx_version);
1036 if (file_version == 0)
1041 gmx::KeyValueTreeBuilder paramsBuilder;
1042 gmx::KeyValueTreeObjectBuilder paramsObj = paramsBuilder.rootObject();
1044 /* Basic inputrec stuff */
1045 serializer->doEnumAsInt(&ir->eI);
1046 if (file_version >= 62)
1048 serializer->doInt64(&ir->nsteps);
1052 serializer->doInt(&idum);
1056 if (file_version >= 62)
1058 serializer->doInt64(&ir->init_step);
1062 serializer->doInt(&idum);
1063 ir->init_step = idum;
1066 serializer->doInt(&ir->simulation_part);
1068 if (file_version >= tpxv_MTS)
1070 serializer->doBool(&ir->useMts);
1071 int numLevels = ir->mtsLevels.size();
1074 serializer->doInt(&numLevels);
1076 ir->mtsLevels.resize(numLevels);
1077 for (auto& mtsLevel : ir->mtsLevels)
1079 int forceGroups = mtsLevel.forceGroups.to_ulong();
1080 serializer->doInt(&forceGroups);
1081 mtsLevel.forceGroups = std::bitset<static_cast<int>(gmx::MtsForceGroups::Count)>(forceGroups);
1082 serializer->doInt(&mtsLevel.stepFactor);
1088 ir->mtsLevels.clear();
1091 if (file_version >= 67)
1093 serializer->doInt(&ir->nstcalcenergy);
1097 ir->nstcalcenergy = 1;
1099 if (file_version >= 81)
1101 serializer->doEnumAsInt(&ir->cutoff_scheme);
1102 if (file_version < 94)
1104 // Need to invert the scheme order
1105 switch (ir->cutoff_scheme)
1107 case (CutoffScheme::Group): ir->cutoff_scheme = CutoffScheme::Verlet; break;
1108 case (CutoffScheme::Verlet): ir->cutoff_scheme = CutoffScheme::Group; break;
1109 default: GMX_RELEASE_ASSERT(false, "Unhandled cutoff scheme type");
1115 ir->cutoff_scheme = CutoffScheme::Group;
1117 serializer->doInt(&idum); /* used to be ns_type; not used anymore */
1118 serializer->doInt(&ir->nstlist);
1119 serializer->doInt(&idum); /* used to be ndelta; not used anymore */
1121 serializer->doReal(&ir->rtpi);
1123 serializer->doInt(&ir->nstcomm);
1124 serializer->doEnumAsInt(&ir->comm_mode);
1126 /* ignore nstcheckpoint */
1127 if (file_version < tpxv_RemoveObsoleteParameters1)
1129 serializer->doInt(&idum);
1132 serializer->doInt(&ir->nstcgsteep);
1134 serializer->doInt(&ir->nbfgscorr);
1136 serializer->doInt(&ir->nstlog);
1137 serializer->doInt(&ir->nstxout);
1138 serializer->doInt(&ir->nstvout);
1139 serializer->doInt(&ir->nstfout);
1140 serializer->doInt(&ir->nstenergy);
1141 serializer->doInt(&ir->nstxout_compressed);
1142 if (file_version >= 59)
1144 serializer->doDouble(&ir->init_t);
1145 serializer->doDouble(&ir->delta_t);
1149 serializer->doReal(&rdum);
1151 serializer->doReal(&rdum);
1154 serializer->doReal(&ir->x_compression_precision);
1155 if (file_version >= 81)
1157 serializer->doReal(&ir->verletbuf_tol);
1161 ir->verletbuf_tol = 0;
1163 serializer->doReal(&ir->rlist);
1164 if (file_version >= 67 && file_version < tpxv_RemoveTwinRange)
1166 if (serializer->reading())
1168 // Reading such a file version could invoke the twin-range
1169 // scheme, about which mdrun should give a fatal error.
1170 real dummy_rlistlong = -1;
1171 serializer->doReal(&dummy_rlistlong);
1173 ir->useTwinRange = (ir->rlist > 0 && (dummy_rlistlong == 0 || dummy_rlistlong > ir->rlist));
1174 // When true, this forces mdrun to issue an error (regardless of
1175 // ir->cutoff_scheme).
1177 // Otherwise, grompp used to set rlistlong actively. Users
1178 // were probably also confused and set rlistlong == rlist.
1179 // However, in all remaining cases, it is safe to let
1180 // mdrun proceed normally.
1185 // No need to read or write anything
1186 ir->useTwinRange = false;
1188 if (file_version >= 82 && file_version != 90)
1190 // Multiple time-stepping is no longer enabled, but the old
1191 // support required the twin-range scheme, for which mdrun
1192 // already emits a fatal error.
1193 int dummy_nstcalclr = -1;
1194 serializer->doInt(&dummy_nstcalclr);
1196 serializer->doEnumAsInt(&ir->coulombtype);
1197 if (file_version >= 81)
1199 serializer->doEnumAsInt(&ir->coulomb_modifier);
1203 ir->coulomb_modifier = (ir->cutoff_scheme == CutoffScheme::Verlet ? InteractionModifiers::PotShift
1204 : InteractionModifiers::None);
1206 serializer->doReal(&ir->rcoulomb_switch);
1207 serializer->doReal(&ir->rcoulomb);
1208 serializer->doEnumAsInt(&ir->vdwtype);
1209 if (file_version >= 81)
1211 serializer->doEnumAsInt(&ir->vdw_modifier);
1215 ir->vdw_modifier = (ir->cutoff_scheme == CutoffScheme::Verlet ? InteractionModifiers::PotShift
1216 : InteractionModifiers::None);
1218 serializer->doReal(&ir->rvdw_switch);
1219 serializer->doReal(&ir->rvdw);
1220 serializer->doEnumAsInt(&ir->eDispCorr);
1221 serializer->doReal(&ir->epsilon_r);
1222 serializer->doReal(&ir->epsilon_rf);
1223 serializer->doReal(&ir->tabext);
1225 // This permits reading a .tpr file that used implicit solvent,
1226 // and later permitting mdrun to refuse to run it.
1227 if (serializer->reading())
1229 if (file_version < tpxv_RemoveImplicitSolvation)
1231 serializer->doInt(&idum);
1232 serializer->doInt(&idum);
1233 serializer->doReal(&rdum);
1234 serializer->doReal(&rdum);
1235 serializer->doInt(&idum);
1236 ir->implicit_solvent = (idum > 0);
1240 ir->implicit_solvent = false;
1242 if (file_version < tpxv_RemoveImplicitSolvation)
1244 serializer->doReal(&rdum);
1245 serializer->doReal(&rdum);
1246 serializer->doReal(&rdum);
1247 serializer->doReal(&rdum);
1248 if (file_version >= 60)
1250 serializer->doReal(&rdum);
1251 serializer->doInt(&idum);
1253 serializer->doReal(&rdum);
1257 if (file_version >= 81)
1259 serializer->doReal(&ir->fourier_spacing);
1263 ir->fourier_spacing = 0.0;
1265 serializer->doInt(&ir->nkx);
1266 serializer->doInt(&ir->nky);
1267 serializer->doInt(&ir->nkz);
1268 serializer->doInt(&ir->pme_order);
1269 serializer->doReal(&ir->ewald_rtol);
1271 if (file_version >= 93)
1273 serializer->doReal(&ir->ewald_rtol_lj);
1277 ir->ewald_rtol_lj = ir->ewald_rtol;
1279 serializer->doEnumAsInt(&ir->ewald_geometry);
1280 serializer->doReal(&ir->epsilon_surface);
1282 /* ignore bOptFFT */
1283 if (file_version < tpxv_RemoveObsoleteParameters1)
1285 serializer->doBool(&bdum);
1288 if (file_version >= 93)
1290 serializer->doEnumAsInt(&ir->ljpme_combination_rule);
1292 serializer->doBool(&ir->bContinuation);
1293 serializer->doEnumAsInt(&ir->etc);
1294 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1295 * but the values 0 and 1 still mean no and
1296 * berendsen temperature coupling, respectively.
1298 if (file_version >= 79)
1300 serializer->doBool(&ir->bPrintNHChains);
1302 if (file_version >= 71)
1304 serializer->doInt(&ir->nsttcouple);
1308 ir->nsttcouple = ir->nstcalcenergy;
1310 serializer->doEnumAsInt(&ir->epc);
1311 serializer->doEnumAsInt(&ir->epct);
1312 if (file_version >= 71)
1314 serializer->doInt(&ir->nstpcouple);
1318 ir->nstpcouple = ir->nstcalcenergy;
1320 serializer->doReal(&ir->tau_p);
1321 serializer->doRvec(&ir->ref_p[XX]);
1322 serializer->doRvec(&ir->ref_p[YY]);
1323 serializer->doRvec(&ir->ref_p[ZZ]);
1324 serializer->doRvec(&ir->compress[XX]);
1325 serializer->doRvec(&ir->compress[YY]);
1326 serializer->doRvec(&ir->compress[ZZ]);
1327 serializer->doEnumAsInt(&ir->refcoord_scaling);
1328 serializer->doRvec(&ir->posres_com);
1329 serializer->doRvec(&ir->posres_comB);
1331 if (file_version < 79)
1333 serializer->doInt(&ir->andersen_seed);
1337 ir->andersen_seed = 0;
1340 serializer->doReal(&ir->shake_tol);
1342 serializer->doEnumAsInt(&ir->efep);
1343 if (serializer->reading())
1345 ir->fepvals = std::make_unique<t_lambda>();
1347 do_fepvals(serializer, ir->fepvals.get(), file_version);
1349 if (file_version >= 79)
1351 serializer->doBool(&ir->bSimTemp);
1354 ir->bSimTemp = TRUE;
1359 ir->bSimTemp = FALSE;
1363 if (serializer->reading())
1365 ir->simtempvals = std::make_unique<t_simtemp>();
1367 do_simtempvals(serializer, ir->simtempvals.get(), ir->fepvals->n_lambda, file_version);
1370 if (file_version >= 79)
1372 serializer->doBool(&ir->bExpanded);
1376 if (serializer->reading())
1378 ir->expandedvals = std::make_unique<t_expanded>();
1380 do_expandedvals(serializer, ir->expandedvals.get(), ir->fepvals.get(), file_version);
1383 serializer->doEnumAsInt(&ir->eDisre);
1384 serializer->doEnumAsInt(&ir->eDisreWeighting);
1385 serializer->doBool(&ir->bDisreMixed);
1386 serializer->doReal(&ir->dr_fc);
1387 serializer->doReal(&ir->dr_tau);
1388 serializer->doInt(&ir->nstdisreout);
1389 serializer->doReal(&ir->orires_fc);
1390 serializer->doReal(&ir->orires_tau);
1391 serializer->doInt(&ir->nstorireout);
1393 /* ignore dihre_fc */
1394 if (file_version < 79)
1396 serializer->doReal(&rdum);
1399 serializer->doReal(&ir->em_stepsize);
1400 serializer->doReal(&ir->em_tol);
1401 serializer->doBool(&ir->bShakeSOR);
1402 serializer->doInt(&ir->niter);
1403 serializer->doReal(&ir->fc_stepsize);
1404 serializer->doEnumAsInt(&ir->eConstrAlg);
1405 serializer->doInt(&ir->nProjOrder);
1406 serializer->doReal(&ir->LincsWarnAngle);
1407 serializer->doInt(&ir->nLincsIter);
1408 serializer->doReal(&ir->bd_fric);
1409 if (file_version >= tpxv_Use64BitRandomSeed)
1411 serializer->doInt64(&ir->ld_seed);
1415 serializer->doInt(&idum);
1419 for (i = 0; i < DIM; i++)
1421 serializer->doRvec(&ir->deform[i]);
1423 serializer->doReal(&ir->cos_accel);
1425 serializer->doInt(&ir->userint1);
1426 serializer->doInt(&ir->userint2);
1427 serializer->doInt(&ir->userint3);
1428 serializer->doInt(&ir->userint4);
1429 serializer->doReal(&ir->userreal1);
1430 serializer->doReal(&ir->userreal2);
1431 serializer->doReal(&ir->userreal3);
1432 serializer->doReal(&ir->userreal4);
1434 /* AdResS is removed, but we need to be able to read old files,
1435 and let mdrun refuse to run them */
1436 if (file_version >= 77 && file_version < tpxv_RemoveAdress)
1438 serializer->doBool(&ir->bAdress);
1441 int idum, numThermoForceGroups, numEnergyGroups;
1444 serializer->doInt(&idum);
1445 serializer->doReal(&rdum);
1446 serializer->doReal(&rdum);
1447 serializer->doReal(&rdum);
1448 serializer->doInt(&idum);
1449 serializer->doInt(&idum);
1450 serializer->doRvec(&rvecdum);
1451 serializer->doInt(&numThermoForceGroups);
1452 serializer->doReal(&rdum);
1453 serializer->doInt(&numEnergyGroups);
1454 serializer->doInt(&idum);
1456 if (numThermoForceGroups > 0)
1458 std::vector<int> idumn(numThermoForceGroups);
1459 serializer->doIntArray(idumn.data(), idumn.size());
1461 if (numEnergyGroups > 0)
1463 std::vector<int> idumn(numEnergyGroups);
1464 serializer->doIntArray(idumn.data(), idumn.size());
1470 ir->bAdress = FALSE;
1475 PullingAlgorithm ePullOld = PullingAlgorithm::Umbrella;
1477 if (file_version >= tpxv_PullCoordTypeGeom)
1479 serializer->doBool(&ir->bPull);
1483 serializer->doEnumAsInt(&ePullOld);
1484 ir->bPull = (ePullOld != PullingAlgorithm::Umbrella);
1485 /* We removed the first ePull=ePullNo for the enum */
1488 case (PullingAlgorithm::Umbrella): break; // this is equal to not using pulling
1489 case (PullingAlgorithm::Constraint): ePullOld = PullingAlgorithm::Umbrella; break;
1490 case (PullingAlgorithm::ConstantForce):
1491 ePullOld = PullingAlgorithm::Constraint;
1493 case (PullingAlgorithm::FlatBottom):
1494 ePullOld = PullingAlgorithm::ConstantForce;
1496 case (PullingAlgorithm::FlatBottomHigh):
1497 ePullOld = PullingAlgorithm::FlatBottom;
1499 case (PullingAlgorithm::External):
1500 ePullOld = PullingAlgorithm::FlatBottomHigh;
1502 case (PullingAlgorithm::Count): ePullOld = PullingAlgorithm::External; break;
1503 default: GMX_RELEASE_ASSERT(false, "Unhandled old pull algorithm");
1508 if (serializer->reading())
1510 ir->pull = std::make_unique<pull_params_t>();
1512 do_pull(serializer, ir->pull.get(), file_version, ePullOld);
1516 if (file_version >= tpxv_AcceleratedWeightHistogram)
1518 serializer->doBool(&ir->bDoAwh);
1522 if (serializer->reading())
1524 ir->awhParams = std::make_unique<gmx::AwhParams>(serializer);
1528 ir->awhParams->serialize(serializer);
1537 /* Enforced rotation */
1538 if (file_version >= 74)
1540 serializer->doBool(&ir->bRot);
1543 if (serializer->reading())
1547 do_rot(serializer, ir->rot);
1555 /* Interactive molecular dynamics */
1556 if (file_version >= tpxv_InteractiveMolecularDynamics)
1558 serializer->doBool(&ir->bIMD);
1561 if (serializer->reading())
1565 do_imd(serializer, ir->imd);
1570 /* We don't support IMD sessions for old .tpr files */
1575 serializer->doInt(&ir->opts.ngtc);
1576 if (file_version >= 69)
1578 serializer->doInt(&ir->opts.nhchainlength);
1582 ir->opts.nhchainlength = 1;
1584 if (serializer->reading() && file_version >= tpxv_RemovedConstantAcceleration
1585 && file_version < tpxv_ReaddedConstantAcceleration)
1591 serializer->doInt(&ir->opts.ngacc);
1593 serializer->doInt(&ir->opts.ngfrz);
1594 serializer->doInt(&ir->opts.ngener);
1596 if (serializer->reading())
1598 snew(ir->opts.nrdf, ir->opts.ngtc);
1599 snew(ir->opts.ref_t, ir->opts.ngtc);
1600 snew(ir->opts.annealing, ir->opts.ngtc);
1601 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1602 snew(ir->opts.anneal_time, ir->opts.ngtc);
1603 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1604 snew(ir->opts.tau_t, ir->opts.ngtc);
1605 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1606 snew(ir->opts.acceleration, ir->opts.ngacc);
1607 snew(ir->opts.egp_flags, ir->opts.ngener * ir->opts.ngener);
1609 if (ir->opts.ngtc > 0)
1611 serializer->doRealArray(ir->opts.nrdf, ir->opts.ngtc);
1612 serializer->doRealArray(ir->opts.ref_t, ir->opts.ngtc);
1613 serializer->doRealArray(ir->opts.tau_t, ir->opts.ngtc);
1615 if (ir->opts.ngfrz > 0)
1617 serializer->doIvecArray(ir->opts.nFreeze, ir->opts.ngfrz);
1619 if (ir->opts.ngacc > 0)
1621 serializer->doRvecArray(ir->opts.acceleration, ir->opts.ngacc);
1623 if (serializer->reading())
1625 ir->useConstantAcceleration = false;
1626 for (int g = 0; g < ir->opts.ngacc; g++)
1628 if (norm2(ir->opts.acceleration[g]) != 0)
1630 ir->useConstantAcceleration = true;
1634 serializer->doIntArray(ir->opts.egp_flags, ir->opts.ngener * ir->opts.ngener);
1636 /* First read the lists with annealing and npoints for each group */
1637 serializer->doEnumArrayAsInt(ir->opts.annealing, ir->opts.ngtc);
1638 serializer->doIntArray(ir->opts.anneal_npoints, ir->opts.ngtc);
1639 for (j = 0; j < (ir->opts.ngtc); j++)
1641 k = ir->opts.anneal_npoints[j];
1642 if (serializer->reading())
1644 snew(ir->opts.anneal_time[j], k);
1645 snew(ir->opts.anneal_temp[j], k);
1647 serializer->doRealArray(ir->opts.anneal_time[j], k);
1648 serializer->doRealArray(ir->opts.anneal_temp[j], k);
1652 serializer->doInt(&ir->nwall);
1653 serializer->doEnumAsInt(&ir->wall_type);
1654 serializer->doReal(&ir->wall_r_linpot);
1655 serializer->doInt(&ir->wall_atomtype[0]);
1656 serializer->doInt(&ir->wall_atomtype[1]);
1657 serializer->doReal(&ir->wall_density[0]);
1658 serializer->doReal(&ir->wall_density[1]);
1659 serializer->doReal(&ir->wall_ewald_zfac);
1662 /* Cosine stuff for electric fields */
1663 if (file_version < tpxv_GenericParamsForElectricField)
1665 do_legacy_efield(serializer, ¶msObj);
1669 if (file_version >= tpxv_ComputationalElectrophysiology)
1671 serializer->doEnumAsInt(&ir->eSwapCoords);
1672 if (ir->eSwapCoords != SwapType::No)
1674 if (serializer->reading())
1678 do_swapcoords_tpx(serializer, ir->swap, file_version);
1682 /* QMMM reading - despite defunct we require reading for backwards
1683 * compability and correct serialization
1686 serializer->doBool(&ir->bQMMM);
1688 serializer->doInt(&qmmmScheme);
1689 real unusedScalefactor;
1690 serializer->doReal(&unusedScalefactor);
1692 // this is still used in Mimic
1693 serializer->doInt(&ir->opts.ngQM);
1694 if (ir->opts.ngQM > 0 && ir->bQMMM)
1696 /* We leave in dummy i/o for removed parameters to avoid
1697 * changing the tpr format.
1699 std::vector<int> dummyIntVec(4 * ir->opts.ngQM, 0);
1700 serializer->doIntArray(dummyIntVec.data(), dummyIntVec.size());
1701 dummyIntVec.clear();
1703 // std::vector<bool> has no data()
1704 std::vector<char> dummyBoolVec(ir->opts.ngQM * sizeof(bool) / sizeof(char));
1705 serializer->doBoolArray(reinterpret_cast<bool*>(dummyBoolVec.data()), dummyBoolVec.size());
1706 dummyBoolVec.clear();
1708 dummyIntVec.resize(2 * ir->opts.ngQM, 0);
1709 serializer->doIntArray(dummyIntVec.data(), dummyIntVec.size());
1710 dummyIntVec.clear();
1712 std::vector<real> dummyRealVec(2 * ir->opts.ngQM, 0);
1713 serializer->doRealArray(dummyRealVec.data(), dummyRealVec.size());
1714 dummyRealVec.clear();
1716 dummyIntVec.resize(3 * ir->opts.ngQM, 0);
1717 serializer->doIntArray(dummyIntVec.data(), dummyIntVec.size());
1718 dummyIntVec.clear();
1720 /* end of QMMM stuff */
1723 if (file_version >= tpxv_GenericParamsForElectricField)
1725 if (serializer->reading())
1727 paramsObj.mergeObject(gmx::deserializeKeyValueTree(serializer));
1731 GMX_RELEASE_ASSERT(ir->params != nullptr,
1732 "Parameters should be present when writing inputrec");
1733 gmx::serializeKeyValueTree(*ir->params, serializer);
1736 if (serializer->reading())
1738 ir->params = new gmx::KeyValueTreeObject(paramsBuilder.build());
1739 // Initialize internal parameters to an empty kvt for all tpr versions
1740 ir->internalParameters = std::make_unique<gmx::KeyValueTreeObject>();
1743 if (file_version >= tpxv_GenericInternalParameters)
1745 if (serializer->reading())
1747 ir->internalParameters =
1748 std::make_unique<gmx::KeyValueTreeObject>(gmx::deserializeKeyValueTree(serializer));
1752 GMX_RELEASE_ASSERT(ir->internalParameters != nullptr,
1753 "Parameters should be present when writing inputrec");
1754 gmx::serializeKeyValueTree(*ir->internalParameters, serializer);
1760 static void do_harm(gmx::ISerializer* serializer, t_iparams* iparams)
1762 serializer->doReal(&iparams->harmonic.rA);
1763 serializer->doReal(&iparams->harmonic.krA);
1764 serializer->doReal(&iparams->harmonic.rB);
1765 serializer->doReal(&iparams->harmonic.krB);
1768 static void do_iparams(gmx::ISerializer* serializer, t_functype ftype, t_iparams* iparams, int file_version)
1781 do_harm(serializer, iparams);
1782 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && serializer->reading())
1784 /* Correct incorrect storage of parameters */
1785 iparams->pdihs.phiB = iparams->pdihs.phiA;
1786 iparams->pdihs.cpB = iparams->pdihs.cpA;
1790 serializer->doReal(&iparams->harmonic.rA);
1791 serializer->doReal(&iparams->harmonic.krA);
1793 case F_LINEAR_ANGLES:
1794 serializer->doReal(&iparams->linangle.klinA);
1795 serializer->doReal(&iparams->linangle.aA);
1796 serializer->doReal(&iparams->linangle.klinB);
1797 serializer->doReal(&iparams->linangle.aB);
1800 serializer->doReal(&iparams->fene.bm);
1801 serializer->doReal(&iparams->fene.kb);
1805 serializer->doReal(&iparams->restraint.lowA);
1806 serializer->doReal(&iparams->restraint.up1A);
1807 serializer->doReal(&iparams->restraint.up2A);
1808 serializer->doReal(&iparams->restraint.kA);
1809 serializer->doReal(&iparams->restraint.lowB);
1810 serializer->doReal(&iparams->restraint.up1B);
1811 serializer->doReal(&iparams->restraint.up2B);
1812 serializer->doReal(&iparams->restraint.kB);
1818 serializer->doReal(&iparams->tab.kA);
1819 serializer->doInt(&iparams->tab.table);
1820 serializer->doReal(&iparams->tab.kB);
1822 case F_CROSS_BOND_BONDS:
1823 serializer->doReal(&iparams->cross_bb.r1e);
1824 serializer->doReal(&iparams->cross_bb.r2e);
1825 serializer->doReal(&iparams->cross_bb.krr);
1827 case F_CROSS_BOND_ANGLES:
1828 serializer->doReal(&iparams->cross_ba.r1e);
1829 serializer->doReal(&iparams->cross_ba.r2e);
1830 serializer->doReal(&iparams->cross_ba.r3e);
1831 serializer->doReal(&iparams->cross_ba.krt);
1833 case F_UREY_BRADLEY:
1834 serializer->doReal(&iparams->u_b.thetaA);
1835 serializer->doReal(&iparams->u_b.kthetaA);
1836 serializer->doReal(&iparams->u_b.r13A);
1837 serializer->doReal(&iparams->u_b.kUBA);
1838 if (file_version >= 79)
1840 serializer->doReal(&iparams->u_b.thetaB);
1841 serializer->doReal(&iparams->u_b.kthetaB);
1842 serializer->doReal(&iparams->u_b.r13B);
1843 serializer->doReal(&iparams->u_b.kUBB);
1847 iparams->u_b.thetaB = iparams->u_b.thetaA;
1848 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1849 iparams->u_b.r13B = iparams->u_b.r13A;
1850 iparams->u_b.kUBB = iparams->u_b.kUBA;
1853 case F_QUARTIC_ANGLES:
1854 serializer->doReal(&iparams->qangle.theta);
1855 serializer->doRealArray(iparams->qangle.c, 5);
1858 serializer->doReal(&iparams->bham.a);
1859 serializer->doReal(&iparams->bham.b);
1860 serializer->doReal(&iparams->bham.c);
1863 serializer->doReal(&iparams->morse.b0A);
1864 serializer->doReal(&iparams->morse.cbA);
1865 serializer->doReal(&iparams->morse.betaA);
1866 if (file_version >= 79)
1868 serializer->doReal(&iparams->morse.b0B);
1869 serializer->doReal(&iparams->morse.cbB);
1870 serializer->doReal(&iparams->morse.betaB);
1874 iparams->morse.b0B = iparams->morse.b0A;
1875 iparams->morse.cbB = iparams->morse.cbA;
1876 iparams->morse.betaB = iparams->morse.betaA;
1880 serializer->doReal(&iparams->cubic.b0);
1881 serializer->doReal(&iparams->cubic.kb);
1882 serializer->doReal(&iparams->cubic.kcub);
1884 case F_CONNBONDS: break;
1885 case F_POLARIZATION: serializer->doReal(&iparams->polarize.alpha); break;
1887 serializer->doReal(&iparams->anharm_polarize.alpha);
1888 serializer->doReal(&iparams->anharm_polarize.drcut);
1889 serializer->doReal(&iparams->anharm_polarize.khyp);
1892 serializer->doReal(&iparams->wpol.al_x);
1893 serializer->doReal(&iparams->wpol.al_y);
1894 serializer->doReal(&iparams->wpol.al_z);
1895 serializer->doReal(&iparams->wpol.rOH);
1896 serializer->doReal(&iparams->wpol.rHH);
1897 serializer->doReal(&iparams->wpol.rOD);
1900 serializer->doReal(&iparams->thole.a);
1901 serializer->doReal(&iparams->thole.alpha1);
1902 serializer->doReal(&iparams->thole.alpha2);
1903 if (file_version < tpxv_RemoveTholeRfac)
1906 serializer->doReal(&noRfac);
1911 serializer->doReal(&iparams->lj.c6);
1912 serializer->doReal(&iparams->lj.c12);
1915 serializer->doReal(&iparams->lj14.c6A);
1916 serializer->doReal(&iparams->lj14.c12A);
1917 serializer->doReal(&iparams->lj14.c6B);
1918 serializer->doReal(&iparams->lj14.c12B);
1921 serializer->doReal(&iparams->ljc14.fqq);
1922 serializer->doReal(&iparams->ljc14.qi);
1923 serializer->doReal(&iparams->ljc14.qj);
1924 serializer->doReal(&iparams->ljc14.c6);
1925 serializer->doReal(&iparams->ljc14.c12);
1927 case F_LJC_PAIRS_NB:
1928 serializer->doReal(&iparams->ljcnb.qi);
1929 serializer->doReal(&iparams->ljcnb.qj);
1930 serializer->doReal(&iparams->ljcnb.c6);
1931 serializer->doReal(&iparams->ljcnb.c12);
1937 serializer->doReal(&iparams->pdihs.phiA);
1938 serializer->doReal(&iparams->pdihs.cpA);
1939 serializer->doReal(&iparams->pdihs.phiB);
1940 serializer->doReal(&iparams->pdihs.cpB);
1941 serializer->doInt(&iparams->pdihs.mult);
1944 serializer->doReal(&iparams->pdihs.phiA);
1945 serializer->doReal(&iparams->pdihs.cpA);
1948 serializer->doInt(&iparams->disres.label);
1949 serializer->doInt(&iparams->disres.type);
1950 serializer->doReal(&iparams->disres.low);
1951 serializer->doReal(&iparams->disres.up1);
1952 serializer->doReal(&iparams->disres.up2);
1953 serializer->doReal(&iparams->disres.kfac);
1956 serializer->doInt(&iparams->orires.ex);
1957 serializer->doInt(&iparams->orires.label);
1958 serializer->doInt(&iparams->orires.power);
1959 serializer->doReal(&iparams->orires.c);
1960 serializer->doReal(&iparams->orires.obs);
1961 serializer->doReal(&iparams->orires.kfac);
1964 if (file_version < 82)
1966 serializer->doInt(&idum);
1967 serializer->doInt(&idum);
1969 serializer->doReal(&iparams->dihres.phiA);
1970 serializer->doReal(&iparams->dihres.dphiA);
1971 serializer->doReal(&iparams->dihres.kfacA);
1972 if (file_version >= 82)
1974 serializer->doReal(&iparams->dihres.phiB);
1975 serializer->doReal(&iparams->dihres.dphiB);
1976 serializer->doReal(&iparams->dihres.kfacB);
1980 iparams->dihres.phiB = iparams->dihres.phiA;
1981 iparams->dihres.dphiB = iparams->dihres.dphiA;
1982 iparams->dihres.kfacB = iparams->dihres.kfacA;
1986 serializer->doRvec(&iparams->posres.pos0A);
1987 serializer->doRvec(&iparams->posres.fcA);
1988 serializer->doRvec(&iparams->posres.pos0B);
1989 serializer->doRvec(&iparams->posres.fcB);
1992 serializer->doInt(&iparams->fbposres.geom);
1993 serializer->doRvec(&iparams->fbposres.pos0);
1994 serializer->doReal(&iparams->fbposres.r);
1995 serializer->doReal(&iparams->fbposres.k);
1997 case F_CBTDIHS: serializer->doRealArray(iparams->cbtdihs.cbtcA, NR_CBTDIHS); break;
1999 // Fall-through intended
2001 /* Fourier dihedrals are internally represented
2002 * as Ryckaert-Bellemans since those are faster to compute.
2004 serializer->doRealArray(iparams->rbdihs.rbcA, NR_RBDIHS);
2005 serializer->doRealArray(iparams->rbdihs.rbcB, NR_RBDIHS);
2009 serializer->doReal(&iparams->constr.dA);
2010 serializer->doReal(&iparams->constr.dB);
2013 serializer->doReal(&iparams->settle.doh);
2014 serializer->doReal(&iparams->settle.dhh);
2016 case F_VSITE1: break; // VSite1 has 0 parameters
2018 case F_VSITE2FD: serializer->doReal(&iparams->vsite.a); break;
2022 serializer->doReal(&iparams->vsite.a);
2023 serializer->doReal(&iparams->vsite.b);
2028 serializer->doReal(&iparams->vsite.a);
2029 serializer->doReal(&iparams->vsite.b);
2030 serializer->doReal(&iparams->vsite.c);
2033 serializer->doInt(&iparams->vsiten.n);
2034 serializer->doReal(&iparams->vsiten.a);
2036 case F_GB12_NOLONGERUSED:
2037 case F_GB13_NOLONGERUSED:
2038 case F_GB14_NOLONGERUSED:
2039 // Implicit solvent parameters can still be read, but never used
2040 if (serializer->reading())
2042 if (file_version < 68)
2044 serializer->doReal(&rdum);
2045 serializer->doReal(&rdum);
2046 serializer->doReal(&rdum);
2047 serializer->doReal(&rdum);
2049 if (file_version < tpxv_RemoveImplicitSolvation)
2051 serializer->doReal(&rdum);
2052 serializer->doReal(&rdum);
2053 serializer->doReal(&rdum);
2054 serializer->doReal(&rdum);
2055 serializer->doReal(&rdum);
2060 serializer->doInt(&iparams->cmap.cmapA);
2061 serializer->doInt(&iparams->cmap.cmapB);
2065 "unknown function type %d (%s) in %s line %d",
2067 interaction_function[ftype].name,
2073 static void do_ilist(gmx::ISerializer* serializer, InteractionList* ilist)
2075 int nr = ilist->size();
2076 serializer->doInt(&nr);
2077 if (serializer->reading())
2079 ilist->iatoms.resize(nr);
2081 serializer->doIntArray(ilist->iatoms.data(), ilist->size());
2084 static void do_ffparams(gmx::ISerializer* serializer, gmx_ffparams_t* ffparams, int file_version)
2086 serializer->doInt(&ffparams->atnr);
2087 int numTypes = ffparams->numTypes();
2088 serializer->doInt(&numTypes);
2089 if (serializer->reading())
2091 ffparams->functype.resize(numTypes);
2092 ffparams->iparams.resize(numTypes);
2094 /* Read/write all the function types */
2095 serializer->doIntArray(ffparams->functype.data(), ffparams->functype.size());
2097 if (file_version >= 66)
2099 serializer->doDouble(&ffparams->reppow);
2103 ffparams->reppow = 12.0;
2106 serializer->doReal(&ffparams->fudgeQQ);
2108 /* Check whether all these function types are supported by the code.
2109 * In practice the code is backwards compatible, which means that the
2110 * numbering may have to be altered from old numbering to new numbering
2112 for (int i = 0; i < ffparams->numTypes(); i++)
2114 if (serializer->reading())
2116 /* Loop over file versions */
2117 for (int k = 0; k < NFTUPD; k++)
2119 /* Compare the read file_version to the update table */
2120 if ((file_version < ftupd[k].fvnr) && (ffparams->functype[i] >= ftupd[k].ftype))
2122 ffparams->functype[i] += 1;
2127 do_iparams(serializer, ffparams->functype[i], &ffparams->iparams[i], file_version);
2131 static void add_settle_atoms(InteractionList* ilist)
2135 /* Settle used to only store the first atom: add the other two */
2136 ilist->iatoms.resize(2 * ilist->size());
2137 for (i = ilist->size() / 4 - 1; i >= 0; i--)
2139 ilist->iatoms[4 * i + 0] = ilist->iatoms[2 * i + 0];
2140 ilist->iatoms[4 * i + 1] = ilist->iatoms[2 * i + 1];
2141 ilist->iatoms[4 * i + 2] = ilist->iatoms[2 * i + 1] + 1;
2142 ilist->iatoms[4 * i + 3] = ilist->iatoms[2 * i + 1] + 2;
2146 static void do_ilists(gmx::ISerializer* serializer, InteractionLists* ilists, int file_version)
2148 GMX_RELEASE_ASSERT(ilists, "Need a valid ilists object");
2149 GMX_RELEASE_ASSERT(ilists->size() == F_NRE,
2150 "The code needs to be in sync with InteractionLists");
2152 for (int j = 0; j < F_NRE; j++)
2154 InteractionList& ilist = (*ilists)[j];
2155 gmx_bool bClear = FALSE;
2156 if (serializer->reading())
2158 for (int k = 0; k < NFTUPD; k++)
2160 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2168 ilist.iatoms.clear();
2172 do_ilist(serializer, &ilist);
2173 if (file_version < 78 && j == F_SETTLE && !ilist.empty())
2175 add_settle_atoms(&ilist);
2181 static void do_block(gmx::ISerializer* serializer, t_block* block)
2183 serializer->doInt(&block->nr);
2184 if (serializer->reading())
2186 if ((block->nalloc_index > 0) && (nullptr != block->index))
2188 sfree(block->index);
2190 block->nalloc_index = block->nr + 1;
2191 snew(block->index, block->nalloc_index);
2193 serializer->doIntArray(block->index, block->nr + 1);
2196 static void doListOfLists(gmx::ISerializer* serializer, gmx::ListOfLists<int>* listOfLists)
2198 int numLists = listOfLists->ssize();
2199 serializer->doInt(&numLists);
2200 int numElements = listOfLists->elementsView().ssize();
2201 serializer->doInt(&numElements);
2202 if (serializer->reading())
2204 std::vector<int> listRanges(numLists + 1);
2205 serializer->doIntArray(listRanges.data(), numLists + 1);
2206 std::vector<int> elements(numElements);
2207 serializer->doIntArray(elements.data(), numElements);
2208 *listOfLists = gmx::ListOfLists<int>(std::move(listRanges), std::move(elements));
2212 serializer->doIntArray(const_cast<int*>(listOfLists->listRangesView().data()), numLists + 1);
2213 serializer->doIntArray(const_cast<int*>(listOfLists->elementsView().data()), numElements);
2217 /* This is a primitive routine to make it possible to translate atomic numbers
2218 * to element names when reading TPR files, without making the Gromacs library
2219 * directory a dependency on mdrun (which is the case if we need elements.dat).
2221 static const char* atomicnumber_to_element(int atomicnumber)
2225 /* This does not have to be complete, so we only include elements likely
2226 * to occur in PDB files.
2228 switch (atomicnumber)
2230 case 1: p = "H"; break;
2231 case 5: p = "B"; break;
2232 case 6: p = "C"; break;
2233 case 7: p = "N"; break;
2234 case 8: p = "O"; break;
2235 case 9: p = "F"; break;
2236 case 11: p = "Na"; break;
2237 case 12: p = "Mg"; break;
2238 case 15: p = "P"; break;
2239 case 16: p = "S"; break;
2240 case 17: p = "Cl"; break;
2241 case 18: p = "Ar"; break;
2242 case 19: p = "K"; break;
2243 case 20: p = "Ca"; break;
2244 case 25: p = "Mn"; break;
2245 case 26: p = "Fe"; break;
2246 case 28: p = "Ni"; break;
2247 case 29: p = "Cu"; break;
2248 case 30: p = "Zn"; break;
2249 case 35: p = "Br"; break;
2250 case 47: p = "Ag"; break;
2251 default: p = ""; break;
2257 static void do_atom(gmx::ISerializer* serializer, t_atom* atom)
2259 serializer->doReal(&atom->m);
2260 serializer->doReal(&atom->q);
2261 serializer->doReal(&atom->mB);
2262 serializer->doReal(&atom->qB);
2263 serializer->doUShort(&atom->type);
2264 serializer->doUShort(&atom->typeB);
2265 serializer->doEnumAsInt(&atom->ptype);
2266 serializer->doInt(&atom->resind);
2267 serializer->doInt(&atom->atomnumber);
2268 if (serializer->reading())
2270 /* Set element string from atomic number if present.
2271 * This routine returns an empty string if the name is not found.
2273 std::strncpy(atom->elem, atomicnumber_to_element(atom->atomnumber), 4);
2274 /* avoid warnings about potentially unterminated string */
2275 atom->elem[3] = '\0';
2279 static void do_grps(gmx::ISerializer* serializer, gmx::ArrayRef<AtomGroupIndices> grps)
2281 for (auto& group : grps)
2283 int size = group.size();
2284 serializer->doInt(&size);
2285 if (serializer->reading())
2289 serializer->doIntArray(group.data(), size);
2293 static void do_symstr(gmx::ISerializer* serializer, char*** nm, t_symtab* symtab)
2297 if (serializer->reading())
2299 serializer->doInt(&ls);
2300 *nm = get_symtab_handle(symtab, ls);
2304 ls = lookup_symtab(symtab, *nm);
2305 serializer->doInt(&ls);
2309 static void do_strstr(gmx::ISerializer* serializer, int nstr, char*** nm, t_symtab* symtab)
2313 for (j = 0; (j < nstr); j++)
2315 do_symstr(serializer, &(nm[j]), symtab);
2319 static void do_resinfo(gmx::ISerializer* serializer, int n, t_resinfo* ri, t_symtab* symtab, int file_version)
2323 for (j = 0; (j < n); j++)
2325 do_symstr(serializer, &(ri[j].name), symtab);
2326 if (file_version >= 63)
2328 serializer->doInt(&ri[j].nr);
2329 serializer->doUChar(&ri[j].ic);
2339 static void do_atoms(gmx::ISerializer* serializer, t_atoms* atoms, t_symtab* symtab, int file_version)
2343 serializer->doInt(&atoms->nr);
2344 serializer->doInt(&atoms->nres);
2345 if (serializer->reading())
2347 /* Since we have always written all t_atom properties in the tpr file
2348 * (at least for all backward compatible versions), we don't store
2349 * but simple set the booleans here.
2351 atoms->haveMass = TRUE;
2352 atoms->haveCharge = TRUE;
2353 atoms->haveType = TRUE;
2354 atoms->haveBState = TRUE;
2355 atoms->havePdbInfo = FALSE;
2357 snew(atoms->atom, atoms->nr);
2358 snew(atoms->atomname, atoms->nr);
2359 snew(atoms->atomtype, atoms->nr);
2360 snew(atoms->atomtypeB, atoms->nr);
2361 snew(atoms->resinfo, atoms->nres);
2362 atoms->pdbinfo = nullptr;
2366 GMX_RELEASE_ASSERT(atoms->haveMass && atoms->haveCharge && atoms->haveType && atoms->haveBState,
2367 "Mass, charge, atomtype and B-state parameters should be present in "
2368 "t_atoms when writing a tpr file");
2370 for (i = 0; (i < atoms->nr); i++)
2372 do_atom(serializer, &atoms->atom[i]);
2374 do_strstr(serializer, atoms->nr, atoms->atomname, symtab);
2375 do_strstr(serializer, atoms->nr, atoms->atomtype, symtab);
2376 do_strstr(serializer, atoms->nr, atoms->atomtypeB, symtab);
2378 do_resinfo(serializer, atoms->nres, atoms->resinfo, symtab, file_version);
2381 static void do_groups(gmx::ISerializer* serializer, SimulationGroups* groups, t_symtab* symtab)
2383 do_grps(serializer, groups->groups);
2384 int numberOfGroupNames = groups->groupNames.size();
2385 serializer->doInt(&numberOfGroupNames);
2386 if (serializer->reading())
2388 groups->groupNames.resize(numberOfGroupNames);
2390 do_strstr(serializer, numberOfGroupNames, groups->groupNames.data(), symtab);
2391 for (auto group : gmx::keysOf(groups->groupNumbers))
2393 int numberOfGroupNumbers = groups->numberOfGroupNumbers(group);
2394 serializer->doInt(&numberOfGroupNumbers);
2395 if (numberOfGroupNumbers != 0)
2397 if (serializer->reading())
2399 groups->groupNumbers[group].resize(numberOfGroupNumbers);
2401 serializer->doUCharArray(groups->groupNumbers[group].data(), numberOfGroupNumbers);
2406 static void do_atomtypes(gmx::ISerializer* serializer, t_atomtypes* atomtypes, int file_version)
2410 serializer->doInt(&atomtypes->nr);
2412 if (serializer->reading())
2414 snew(atomtypes->atomnumber, j);
2416 if (serializer->reading() && file_version < tpxv_RemoveImplicitSolvation)
2418 std::vector<real> dummy(atomtypes->nr, 0);
2419 serializer->doRealArray(dummy.data(), dummy.size());
2420 serializer->doRealArray(dummy.data(), dummy.size());
2421 serializer->doRealArray(dummy.data(), dummy.size());
2423 serializer->doIntArray(atomtypes->atomnumber, j);
2425 if (serializer->reading() && file_version >= 60 && file_version < tpxv_RemoveImplicitSolvation)
2427 std::vector<real> dummy(atomtypes->nr, 0);
2428 serializer->doRealArray(dummy.data(), dummy.size());
2429 serializer->doRealArray(dummy.data(), dummy.size());
2433 static void do_symtab(gmx::ISerializer* serializer, t_symtab* symtab)
2438 serializer->doInt(&symtab->nr);
2440 if (serializer->reading())
2442 snew(symtab->symbuf, 1);
2443 symbuf = symtab->symbuf;
2444 symbuf->bufsize = nr;
2445 snew(symbuf->buf, nr);
2446 for (i = 0; (i < nr); i++)
2449 serializer->doString(&buf);
2450 symbuf->buf[i] = gmx_strdup(buf.c_str());
2455 symbuf = symtab->symbuf;
2456 while (symbuf != nullptr)
2458 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2460 std::string buf = symbuf->buf[i];
2461 serializer->doString(&buf);
2464 symbuf = symbuf->next;
2468 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2473 static void do_cmap(gmx::ISerializer* serializer, gmx_cmap_t* cmap_grid)
2476 int ngrid = cmap_grid->cmapdata.size();
2477 serializer->doInt(&ngrid);
2478 serializer->doInt(&cmap_grid->grid_spacing);
2480 int gs = cmap_grid->grid_spacing;
2481 int nelem = gs * gs;
2483 if (serializer->reading())
2485 cmap_grid->cmapdata.resize(ngrid);
2487 for (int i = 0; i < ngrid; i++)
2489 cmap_grid->cmapdata[i].cmap.resize(4 * nelem);
2493 for (int i = 0; i < ngrid; i++)
2495 for (int j = 0; j < nelem; j++)
2497 serializer->doReal(&cmap_grid->cmapdata[i].cmap[j * 4]);
2498 serializer->doReal(&cmap_grid->cmapdata[i].cmap[j * 4 + 1]);
2499 serializer->doReal(&cmap_grid->cmapdata[i].cmap[j * 4 + 2]);
2500 serializer->doReal(&cmap_grid->cmapdata[i].cmap[j * 4 + 3]);
2506 static void do_moltype(gmx::ISerializer* serializer, gmx_moltype_t* molt, t_symtab* symtab, int file_version)
2508 do_symstr(serializer, &(molt->name), symtab);
2510 do_atoms(serializer, &molt->atoms, symtab, file_version);
2512 do_ilists(serializer, &molt->ilist, file_version);
2514 /* TODO: Remove the obsolete charge group index from the file */
2516 cgs.nr = molt->atoms.nr;
2517 cgs.nalloc_index = cgs.nr + 1;
2518 snew(cgs.index, cgs.nalloc_index);
2519 for (int i = 0; i < cgs.nr + 1; i++)
2523 do_block(serializer, &cgs);
2526 /* This used to be in the atoms struct */
2527 doListOfLists(serializer, &molt->excls);
2530 static void do_molblock(gmx::ISerializer* serializer, gmx_molblock_t* molb, int numAtomsPerMolecule)
2532 serializer->doInt(&molb->type);
2533 serializer->doInt(&molb->nmol);
2534 /* To maintain forward topology reading compatibility, we store #atoms.
2535 * TODO: Change this to conditional reading of a dummy int when we
2536 * increase tpx_generation.
2538 serializer->doInt(&numAtomsPerMolecule);
2539 /* Position restraint coordinates */
2540 int numPosres_xA = molb->posres_xA.size();
2541 serializer->doInt(&numPosres_xA);
2542 if (numPosres_xA > 0)
2544 if (serializer->reading())
2546 molb->posres_xA.resize(numPosres_xA);
2548 serializer->doRvecArray(as_rvec_array(molb->posres_xA.data()), numPosres_xA);
2550 int numPosres_xB = molb->posres_xB.size();
2551 serializer->doInt(&numPosres_xB);
2552 if (numPosres_xB > 0)
2554 if (serializer->reading())
2556 molb->posres_xB.resize(numPosres_xB);
2558 serializer->doRvecArray(as_rvec_array(molb->posres_xB.data()), numPosres_xB);
2562 static void set_disres_npair(gmx_mtop_t* mtop)
2564 gmx::ArrayRef<t_iparams> ip = mtop->ffparams.iparams;
2566 for (const auto ilist : IListRange(*mtop))
2568 const InteractionList& il = ilist.list()[F_DISRES];
2572 gmx::ArrayRef<const int> a = il.iatoms;
2574 for (int i = 0; i < il.size(); i += 3)
2577 if (i + 3 == il.size() || ip[a[i]].disres.label != ip[a[i + 3]].disres.label)
2579 ip[a[i]].disres.npair = npair;
2587 static void do_mtop(gmx::ISerializer* serializer, gmx_mtop_t* mtop, int file_version)
2589 do_symtab(serializer, &(mtop->symtab));
2591 do_symstr(serializer, &(mtop->name), &(mtop->symtab));
2593 do_ffparams(serializer, &mtop->ffparams, file_version);
2595 int nmoltype = mtop->moltype.size();
2596 serializer->doInt(&nmoltype);
2597 if (serializer->reading())
2599 mtop->moltype.resize(nmoltype);
2601 for (gmx_moltype_t& moltype : mtop->moltype)
2603 do_moltype(serializer, &moltype, &mtop->symtab, file_version);
2606 int nmolblock = mtop->molblock.size();
2607 serializer->doInt(&nmolblock);
2608 if (serializer->reading())
2610 mtop->molblock.resize(nmolblock);
2612 for (gmx_molblock_t& molblock : mtop->molblock)
2614 int numAtomsPerMolecule = (serializer->reading() ? 0 : mtop->moltype[molblock.type].atoms.nr);
2615 do_molblock(serializer, &molblock, numAtomsPerMolecule);
2617 serializer->doInt(&mtop->natoms);
2619 if (file_version >= tpxv_IntermolecularBondeds)
2621 serializer->doBool(&mtop->bIntermolecularInteractions);
2622 if (mtop->bIntermolecularInteractions)
2624 if (serializer->reading())
2626 mtop->intermolecular_ilist = std::make_unique<InteractionLists>();
2628 do_ilists(serializer, mtop->intermolecular_ilist.get(), file_version);
2633 mtop->bIntermolecularInteractions = FALSE;
2636 do_atomtypes(serializer, &(mtop->atomtypes), file_version);
2638 if (file_version >= 65)
2640 do_cmap(serializer, &mtop->ffparams.cmap_grid);
2644 mtop->ffparams.cmap_grid.grid_spacing = 0;
2645 mtop->ffparams.cmap_grid.cmapdata.clear();
2648 do_groups(serializer, &mtop->groups, &(mtop->symtab));
2650 mtop->haveMoleculeIndices = true;
2652 if (file_version >= tpxv_StoreNonBondedInteractionExclusionGroup)
2654 std::int64_t intermolecularExclusionGroupSize = gmx::ssize(mtop->intermolecularExclusionGroup);
2655 serializer->doInt64(&intermolecularExclusionGroupSize);
2656 GMX_RELEASE_ASSERT(intermolecularExclusionGroupSize >= 0,
2657 "Number of atoms with excluded intermolecular non-bonded interactions "
2659 mtop->intermolecularExclusionGroup.resize(intermolecularExclusionGroupSize); // no effect when writing
2660 serializer->doIntArray(mtop->intermolecularExclusionGroup.data(),
2661 mtop->intermolecularExclusionGroup.size());
2664 if (serializer->reading())
2666 close_symtab(&(mtop->symtab));
2671 * Read the first part of the TPR file to find general system information.
2673 * If \p TopOnlyOK is true then we can read even future versions
2674 * of tpx files, provided the \p fileGeneration hasn't changed.
2675 * If it is false, we need the \p ir too, and bail out
2676 * if the file is newer than the program.
2678 * The version and generation of the topology (see top of this file)
2679 * are returned in the two last arguments, if those arguments are non-nullptr.
2681 * If possible, we will read the \p ir even when \p TopOnlyOK is true.
2683 * \param[in,out] serializer The serializer used to handle header processing.
2684 * \param[in,out] tpx File header datastructure.
2685 * \param[in] filename The name of the file being read/written
2686 * \param[in,out] fio File handle.
2687 * \param[in] TopOnlyOK If not reading \p ir is fine or not.
2689 static void do_tpxheader(gmx::FileIOXdrSerializer* serializer,
2691 const char* filename,
2699 /* XDR binary topology file */
2700 precision = sizeof(real);
2702 std::string fileTag;
2703 if (serializer->reading())
2705 serializer->doString(&buf);
2706 if (std::strncmp(buf.c_str(), "VERSION", 7) != 0)
2710 "Can not read file %s,\n"
2711 " this file is from a GROMACS version which is older than 2.0\n"
2712 " Make a new one with grompp or use a gro or pdb file, if possible",
2715 // We need to know the precision used to write the TPR file, to match it
2716 // to the precision of the currently running binary. If the precisions match
2717 // there is no problem, but mismatching precision needs to be accounted for
2718 // by reading into temporary variables of the correct precision instead
2719 // of the desired target datastructures.
2720 serializer->doInt(&precision);
2721 tpx->isDouble = (precision == sizeof(double));
2722 if ((precision != sizeof(float)) && !tpx->isDouble)
2725 "Unknown precision in file %s: real is %d bytes "
2726 "instead of %zu or %zu",
2732 gmx_fio_setprecision(fio, tpx->isDouble);
2734 "Reading file %s, %s (%s precision)\n",
2737 tpx->isDouble ? "double" : "single");
2741 buf = gmx::formatString("VERSION %s", gmx_version());
2742 serializer->doString(&buf);
2743 gmx_fio_setprecision(fio, tpx->isDouble);
2744 serializer->doInt(&precision);
2748 /* Check versions! */
2749 serializer->doInt(&tpx->fileVersion);
2751 /* This is for backward compatibility with development versions 77-79
2752 * where the tag was, mistakenly, placed before the generation,
2753 * which would cause a segv instead of a proper error message
2754 * when reading the topology only from tpx with <77 code.
2756 if (tpx->fileVersion >= 77 && tpx->fileVersion <= 79)
2758 serializer->doString(&fileTag);
2761 serializer->doInt(&tpx->fileGeneration);
2763 if (tpx->fileVersion >= 81)
2765 serializer->doString(&fileTag);
2767 if (serializer->reading())
2769 if (tpx->fileVersion < 77)
2771 /* Versions before 77 don't have the tag, set it to release */
2772 fileTag = gmx::formatString("%s", TPX_TAG_RELEASE);
2775 if (fileTag != tpx_tag)
2777 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n", fileTag.c_str(), tpx_tag.c_str());
2779 /* We only support reading tpx files with the same tag as the code
2780 * or tpx files with the release tag and with lower version number.
2782 if (fileTag != TPX_TAG_RELEASE && tpx->fileVersion < tpx_version)
2785 "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' "
2786 "with program for tpx version %d, tag '%s'",
2796 if ((tpx->fileVersion <= tpx_incompatible_version)
2797 || ((tpx->fileVersion > tpx_version) && !TopOnlyOK) || (tpx->fileGeneration > tpx_generation)
2798 || tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
2801 "reading tpx file (%s) version %d with version %d program",
2807 serializer->doInt(&tpx->natoms);
2808 serializer->doInt(&tpx->ngtc);
2810 if (tpx->fileVersion < 62)
2812 serializer->doInt(&idum);
2813 serializer->doReal(&rdum);
2815 if (tpx->fileVersion >= 79)
2817 serializer->doInt(&tpx->fep_state);
2819 serializer->doReal(&tpx->lambda);
2820 serializer->doBool(&tpx->bIr);
2821 serializer->doBool(&tpx->bTop);
2822 serializer->doBool(&tpx->bX);
2823 serializer->doBool(&tpx->bV);
2824 serializer->doBool(&tpx->bF);
2825 serializer->doBool(&tpx->bBox);
2827 if (tpx->fileVersion >= tpxv_AddSizeField && tpx->fileGeneration >= 27)
2829 if (!serializer->reading())
2831 GMX_RELEASE_ASSERT(tpx->sizeOfTprBody != 0,
2832 "Not possible to write new file with zero TPR body size");
2834 serializer->doInt64(&tpx->sizeOfTprBody);
2837 if ((tpx->fileGeneration > tpx_generation))
2839 /* This can only happen if TopOnlyOK=TRUE */
2844 #define do_test(serializer, b, p) \
2845 if ((serializer)->reading() && ((p) != nullptr) && !(b)) \
2846 gmx_fatal(FARGS, "No %s in input file", #p)
2849 * Process the first part of the TPR into the state datastructure.
2851 * Due to the structure of the legacy code, it is necessary
2852 * to split up the state reading into two parts, with the
2853 * box and legacy temperature coupling processed before the
2854 * topology datastructures.
2856 * See the documentation for do_tpx_body for the correct order of
2857 * the operations for reading a tpr file.
2859 * \param[in] serializer Abstract serializer used to read/write data.
2860 * \param[in] tpx The file header data.
2861 * \param[in, out] state Global state data.
2863 static void do_tpx_state_first(gmx::ISerializer* serializer, TpxFileHeader* tpx, t_state* state)
2865 if (serializer->reading())
2868 init_gtc_state(state, tpx->ngtc, 0, 0);
2870 do_test(serializer, tpx->bBox, state->box);
2873 serializer->doRvecArray(state->box, DIM);
2874 if (tpx->fileVersion >= 51)
2876 serializer->doRvecArray(state->box_rel, DIM);
2880 /* We initialize box_rel after reading the inputrec */
2881 clear_mat(state->box_rel);
2883 serializer->doRvecArray(state->boxv, DIM);
2884 if (tpx->fileVersion < 56)
2887 serializer->doRvecArray(mdum, DIM);
2891 if (state->ngtc > 0)
2894 snew(dumv, state->ngtc);
2895 if (tpx->fileVersion < 69)
2897 serializer->doRealArray(dumv, state->ngtc);
2899 /* These used to be the Berendsen tcoupl_lambda's */
2900 serializer->doRealArray(dumv, state->ngtc);
2906 * Process global topology data.
2908 * See the documentation for do_tpx_body for the correct order of
2909 * the operations for reading a tpr file.
2911 * \param[in] serializer Abstract serializer used to read/write data.
2912 * \param[in] tpx The file header data.
2913 * \param[in,out] mtop Global topology.
2915 static void do_tpx_mtop(gmx::ISerializer* serializer, TpxFileHeader* tpx, gmx_mtop_t* mtop)
2917 do_test(serializer, tpx->bTop, mtop);
2922 do_mtop(serializer, mtop, tpx->fileVersion);
2923 set_disres_npair(mtop);
2929 do_mtop(serializer, &dum_top, tpx->fileVersion);
2934 * Process coordinate vectors for state data.
2936 * Main part of state gets processed here.
2938 * See the documentation for do_tpx_body for the correct order of
2939 * the operations for reading a tpr file.
2941 * \param[in] serializer Abstract serializer used to read/write data.
2942 * \param[in] tpx The file header data.
2943 * \param[in,out] state Global state data.
2944 * \param[in,out] x Individual coordinates for processing, deprecated.
2945 * \param[in,out] v Individual velocities for processing, deprecated.
2947 static void do_tpx_state_second(gmx::ISerializer* serializer, TpxFileHeader* tpx, t_state* state, rvec* x, rvec* v)
2949 if (!serializer->reading())
2952 x == nullptr && v == nullptr,
2953 "Passing separate x and v pointers to do_tpx() is not supported when writing");
2957 GMX_RELEASE_ASSERT(!(x == nullptr && v != nullptr),
2958 "Passing x==NULL and v!=NULL is not supported");
2961 if (serializer->reading())
2965 // v is also nullptr by the above assertion, so we may
2966 // need to make memory in state for storing the contents
2970 state->flags |= enumValueToBitMask(StateEntry::X);
2974 state->flags |= enumValueToBitMask(StateEntry::V);
2976 state_change_natoms(state, tpx->natoms);
2982 x = state->x.rvec_array();
2983 v = state->v.rvec_array();
2985 do_test(serializer, tpx->bX, x);
2988 if (serializer->reading())
2990 state->flags |= enumValueToBitMask(StateEntry::X);
2992 serializer->doRvecArray(x, tpx->natoms);
2995 do_test(serializer, tpx->bV, v);
2998 if (serializer->reading())
3000 state->flags |= enumValueToBitMask(StateEntry::V);
3004 std::vector<gmx::RVec> dummyVelocities(tpx->natoms);
3005 serializer->doRvecArray(as_rvec_array(dummyVelocities.data()), tpx->natoms);
3009 serializer->doRvecArray(v, tpx->natoms);
3013 // No need to run do_test when the last argument is NULL
3016 std::vector<gmx::RVec> dummyForces(state->natoms);
3017 serializer->doRvecArray(as_rvec_array(dummyForces.data()), tpx->natoms);
3021 * Process simulation parameters.
3023 * See the documentation for do_tpx_body for the correct order of
3024 * the operations for reading a tpr file.
3026 * \param[in] serializer Abstract serializer used to read/write data.
3027 * \param[in] tpx The file header data.
3028 * \param[in,out] ir Datastructure with simulation parameters.
3030 static PbcType do_tpx_ir(gmx::ISerializer* serializer, TpxFileHeader* tpx, t_inputrec* ir)
3033 gmx_bool bPeriodicMols;
3035 /* Starting with tpx version 26, we have the inputrec
3036 * at the end of the file, so we can ignore it
3037 * if the file is never than the software (but still the
3038 * same generation - see comments at the top of this file.
3042 pbcType = PbcType::Unset;
3043 bPeriodicMols = FALSE;
3045 do_test(serializer, tpx->bIr, ir);
3048 if (tpx->fileVersion >= 53)
3050 /* Removed the pbc info from do_inputrec, since we always want it */
3051 if (!serializer->reading())
3053 pbcType = ir->pbcType;
3054 bPeriodicMols = ir->bPeriodicMols;
3056 serializer->doInt(reinterpret_cast<int*>(&pbcType));
3057 serializer->doBool(&bPeriodicMols);
3059 if (tpx->fileGeneration <= tpx_generation && ir)
3061 do_inputrec(serializer, ir, tpx->fileVersion);
3062 if (tpx->fileVersion < 53)
3064 pbcType = ir->pbcType;
3065 bPeriodicMols = ir->bPeriodicMols;
3068 if (serializer->reading() && ir && tpx->fileVersion >= 53)
3070 /* We need to do this after do_inputrec, since that initializes ir */
3071 ir->pbcType = pbcType;
3072 ir->bPeriodicMols = bPeriodicMols;
3079 * Correct and finalize read information.
3081 * If \p state is nullptr, skip the parts dependent on it.
3083 * See the documentation for do_tpx_body for the correct order of
3084 * the operations for reading a tpr file.
3086 * \param[in] tpx The file header used to check version numbers.
3087 * \param[out] ir Input rec that needs correction.
3088 * \param[out] state State needing correction.
3089 * \param[out] mtop Topology to finalize.
3091 static void do_tpx_finalize(TpxFileHeader* tpx, t_inputrec* ir, t_state* state, gmx_mtop_t* mtop)
3093 if (tpx->fileVersion < 51 && state)
3095 set_box_rel(ir, state);
3099 if (state && state->ngtc == 0)
3101 /* Reading old version without tcoupl state data: set it */
3102 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3104 if (tpx->bTop && mtop)
3106 if (tpx->fileVersion < 57)
3108 ir->eDisre = !mtop->moltype[0].ilist[F_DISRES].empty()
3109 ? DistanceRestraintRefinement::Simple
3110 : DistanceRestraintRefinement::None;
3117 * Process TPR data for file reading/writing.
3119 * The TPR file gets processed in in four stages due to the organization
3120 * of the data within it.
3122 * First, state data for the box is processed in do_tpx_state_first.
3123 * This is followed by processing the topology in do_tpx_mtop.
3124 * Coordinate and velocity vectors are handled next in do_tpx_state_second.
3125 * The last file information processed is the collection of simulation parameters in do_tpx_ir.
3126 * When reading, a final processing step is undertaken at the end.
3128 * \param[in] serializer Abstract serializer used to read/write data.
3129 * \param[in] tpx The file header data.
3130 * \param[in,out] ir Datastructures with simulation parameters.
3131 * \param[in,out] state Global state data.
3132 * \param[in,out] x Individual coordinates for processing, deprecated.
3133 * \param[in,out] v Individual velocities for processing, deprecated.
3134 * \param[in,out] mtop Global topology.
3136 static PbcType do_tpx_body(gmx::ISerializer* serializer,
3146 do_tpx_state_first(serializer, tpx, state);
3148 do_tpx_mtop(serializer, tpx, mtop);
3151 do_tpx_state_second(serializer, tpx, state, x, v);
3153 PbcType pbcType = do_tpx_ir(serializer, tpx, ir);
3154 if (serializer->reading())
3156 do_tpx_finalize(tpx, ir, state, mtop);
3162 * Overload for do_tpx_body that defaults to state vectors being nullptr.
3164 * \param[in] serializer Abstract serializer used to read/write data.
3165 * \param[in] tpx The file header data.
3166 * \param[in,out] ir Datastructures with simulation parameters.
3167 * \param[in,out] mtop Global topology.
3169 static PbcType do_tpx_body(gmx::ISerializer* serializer, TpxFileHeader* tpx, t_inputrec* ir, gmx_mtop_t* mtop)
3171 return do_tpx_body(serializer, tpx, ir, nullptr, nullptr, nullptr, mtop);
3174 static t_fileio* open_tpx(const char* fn, const char* mode)
3176 return gmx_fio_open(fn, mode);
3179 static void close_tpx(t_fileio* fio)
3185 * Fill information into the header only from state before writing.
3187 * Populating the header needs to be independent from writing the information
3188 * to file to allow things like writing the raw byte stream.
3190 * \param[in] state The current simulation state. Can't write without it.
3191 * \param[in] ir Parameter and system information.
3192 * \param[in] mtop Global topology.
3193 * \returns Fully populated header.
3195 static TpxFileHeader populateTpxHeader(const t_state& state, const t_inputrec* ir, const gmx_mtop_t* mtop)
3197 TpxFileHeader header;
3198 header.natoms = state.natoms;
3199 header.ngtc = state.ngtc;
3200 header.fep_state = state.fep_state;
3201 header.lambda = state.lambda[FreeEnergyPerturbationCouplingType::Fep];
3202 header.bIr = ir != nullptr;
3203 header.bTop = mtop != nullptr;
3204 header.bX = (state.flags & enumValueToBitMask(StateEntry::X)) != 0;
3205 header.bV = (state.flags & enumValueToBitMask(StateEntry::V)) != 0;
3208 header.fileVersion = tpx_version;
3209 header.fileGeneration = tpx_generation;
3210 header.isDouble = (sizeof(real) == sizeof(double));
3216 * Process the body of a TPR file as an opaque data buffer.
3218 * Reads/writes the information in \p buffer from/to the \p serializer
3219 * provided to the function. Does not interact with the actual
3220 * TPR datastructures but with an in memory representation of the
3221 * data, so that this data can be efficiently read or written from/to
3222 * an original source.
3224 * \param[in] serializer The abstract serializer used for reading or writing
3225 * the information in \p buffer.
3226 * \param[in,out] buffer Information from TPR file as char buffer.
3228 static void doTpxBodyBuffer(gmx::ISerializer* serializer, gmx::ArrayRef<char> buffer)
3230 serializer->doOpaque(buffer.data(), buffer.size());
3234 * Populates simulation datastructures.
3236 * Here the information from the serialization interface \p serializer
3237 * is used to first populate the datastructures containing the simulation
3238 * information. Depending on the version found in the header \p tpx,
3239 * this is done using the new reading of the data as one block from disk,
3240 * followed by complete deserialization of the information read from there.
3241 * Otherwise, the datastructures are populated as before one by one from disk.
3242 * The second version is the default for the legacy tools that read the
3243 * coordinates and velocities separate from the state.
3245 * After reading in the data, a separate buffer is populated from them
3246 * containing only \p ir and \p mtop that can be communicated directly
3247 * to nodes needing the information to set up a simulation.
3249 * \param[in] tpx The file header.
3250 * \param[in] serializer The Serialization interface used to read the TPR.
3251 * \param[out] ir Input rec to populate.
3252 * \param[out] state State vectors to populate.
3253 * \param[out] x Coordinates to populate if needed.
3254 * \param[out] v Velocities to populate if needed.
3255 * \param[out] mtop Global topology to populate.
3257 * \returns Partial de-serialized TPR used for communication to nodes.
3259 static PartialDeserializedTprFile readTpxBody(TpxFileHeader* tpx,
3260 gmx::ISerializer* serializer,
3267 PartialDeserializedTprFile partialDeserializedTpr;
3268 if (tpx->fileVersion >= tpxv_AddSizeField && tpx->fileGeneration >= 27)
3270 partialDeserializedTpr.body.resize(tpx->sizeOfTprBody);
3271 partialDeserializedTpr.header = *tpx;
3272 doTpxBodyBuffer(serializer, partialDeserializedTpr.body);
3274 partialDeserializedTpr.pbcType =
3275 completeTprDeserialization(&partialDeserializedTpr, ir, state, x, v, mtop);
3279 partialDeserializedTpr.pbcType = do_tpx_body(serializer, tpx, ir, state, x, v, mtop);
3281 // Update header to system info for communication to nodes.
3282 // As we only need to communicate the inputrec and mtop to other nodes,
3283 // we prepare a new char buffer with the information we have already read
3285 partialDeserializedTpr.header = populateTpxHeader(*state, ir, mtop);
3286 // Long-term we should move to use little endian in files to avoid extra byte swapping,
3287 // but since we just used the default XDR format (which is big endian) for the TPR
3288 // header it would cause third-party libraries reading our raw data to tear their hair
3289 // if we swap the endian in the middle of the file, so we stick to big endian in the
3290 // TPR file for now - and thus we ask the serializer to swap if this host is little endian.
3291 gmx::InMemorySerializer tprBodySerializer(gmx::EndianSwapBehavior::SwapIfHostIsLittleEndian);
3292 do_tpx_body(&tprBodySerializer, &partialDeserializedTpr.header, ir, mtop);
3293 partialDeserializedTpr.body = tprBodySerializer.finishAndGetBuffer();
3295 return partialDeserializedTpr;
3298 /************************************************************
3300 * The following routines are the exported ones
3302 ************************************************************/
3304 TpxFileHeader readTpxHeader(const char* fileName, bool canReadTopologyOnly)
3308 fio = open_tpx(fileName, "r");
3309 gmx::FileIOXdrSerializer serializer(fio);
3312 do_tpxheader(&serializer, &tpx, fileName, fio, canReadTopologyOnly);
3317 void write_tpx_state(const char* fn, const t_inputrec* ir, const t_state* state, const gmx_mtop_t& mtop)
3319 /* To write a state, we first need to write the state information to a buffer before
3320 * we append the raw bytes to the file. For this, the header information needs to be
3321 * populated before we write the main body because it has some information that is
3322 * otherwise not available.
3327 TpxFileHeader tpx = populateTpxHeader(*state, ir, &mtop);
3328 // Long-term we should move to use little endian in files to avoid extra byte swapping,
3329 // but since we just used the default XDR format (which is big endian) for the TPR
3330 // header it would cause third-party libraries reading our raw data to tear their hair
3331 // if we swap the endian in the middle of the file, so we stick to big endian in the
3332 // TPR file for now - and thus we ask the serializer to swap if this host is little endian.
3333 gmx::InMemorySerializer tprBodySerializer(gmx::EndianSwapBehavior::SwapIfHostIsLittleEndian);
3335 do_tpx_body(&tprBodySerializer,
3337 const_cast<t_inputrec*>(ir),
3338 const_cast<t_state*>(state),
3341 const_cast<gmx_mtop_t*>(&mtop));
3343 std::vector<char> tprBody = tprBodySerializer.finishAndGetBuffer();
3344 tpx.sizeOfTprBody = tprBody.size();
3346 fio = open_tpx(fn, "w");
3347 gmx::FileIOXdrSerializer serializer(fio);
3348 do_tpxheader(&serializer, &tpx, fn, fio, ir == nullptr);
3349 doTpxBodyBuffer(&serializer, tprBody);
3354 PbcType completeTprDeserialization(PartialDeserializedTprFile* partialDeserializedTpr,
3361 // Long-term we should move to use little endian in files to avoid extra byte swapping,
3362 // but since we just used the default XDR format (which is big endian) for the TPR
3363 // header it would cause third-party libraries reading our raw data to tear their hair
3364 // if we swap the endian in the middle of the file, so we stick to big endian in the
3365 // TPR file for now - and thus we ask the serializer to swap if this host is little endian.
3366 gmx::InMemoryDeserializer tprBodyDeserializer(partialDeserializedTpr->body,
3367 partialDeserializedTpr->header.isDouble,
3368 gmx::EndianSwapBehavior::SwapIfHostIsLittleEndian);
3369 return do_tpx_body(&tprBodyDeserializer, &partialDeserializedTpr->header, ir, state, x, v, mtop);
3372 PbcType completeTprDeserialization(PartialDeserializedTprFile* partialDeserializedTpr,
3376 return completeTprDeserialization(partialDeserializedTpr, ir, nullptr, nullptr, nullptr, mtop);
3379 PartialDeserializedTprFile read_tpx_state(const char* fn, t_inputrec* ir, t_state* state, gmx_mtop_t* mtop)
3382 fio = open_tpx(fn, "r");
3383 gmx::FileIOXdrSerializer serializer(fio);
3384 PartialDeserializedTprFile partialDeserializedTpr;
3385 do_tpxheader(&serializer, &partialDeserializedTpr.header, fn, fio, ir == nullptr);
3386 partialDeserializedTpr =
3387 readTpxBody(&partialDeserializedTpr.header, &serializer, ir, state, nullptr, nullptr, mtop);
3389 return partialDeserializedTpr;
3392 PbcType read_tpx(const char* fn, t_inputrec* ir, matrix box, int* natoms, rvec* x, rvec* v, gmx_mtop_t* mtop)
3398 fio = open_tpx(fn, "r");
3399 gmx::FileIOXdrSerializer serializer(fio);
3400 do_tpxheader(&serializer, &tpx, fn, fio, ir == nullptr);
3401 PartialDeserializedTprFile partialDeserializedTpr =
3402 readTpxBody(&tpx, &serializer, ir, &state, x, v, mtop);
3404 if (mtop != nullptr && natoms != nullptr)
3406 *natoms = mtop->natoms;
3410 copy_mat(state.box, box);
3412 return partialDeserializedTpr.pbcType;
3415 PbcType read_tpx_top(const char* fn, t_inputrec* ir, matrix box, int* natoms, rvec* x, rvec* v, t_topology* top)
3420 pbcType = read_tpx(fn, ir, box, natoms, x, v, &mtop);
3422 *top = gmx_mtop_t_to_t_topology(&mtop, true);
3427 gmx_bool fn2bTPX(const char* file)
3429 return (efTPR == fn2ftp(file));
3432 void pr_tpxheader(FILE* fp, int indent, const char* title, const TpxFileHeader* sh)
3434 if (available(fp, sh, indent, title))
3436 indent = pr_title(fp, indent, title);
3437 pr_indent(fp, indent);
3438 fprintf(fp, "bIr = %spresent\n", sh->bIr ? "" : "not ");
3439 pr_indent(fp, indent);
3440 fprintf(fp, "bBox = %spresent\n", sh->bBox ? "" : "not ");
3441 pr_indent(fp, indent);
3442 fprintf(fp, "bTop = %spresent\n", sh->bTop ? "" : "not ");
3443 pr_indent(fp, indent);
3444 fprintf(fp, "bX = %spresent\n", sh->bX ? "" : "not ");
3445 pr_indent(fp, indent);
3446 fprintf(fp, "bV = %spresent\n", sh->bV ? "" : "not ");
3447 pr_indent(fp, indent);
3448 fprintf(fp, "bF = %spresent\n", sh->bF ? "" : "not ");
3450 pr_indent(fp, indent);
3451 fprintf(fp, "natoms = %d\n", sh->natoms);
3452 pr_indent(fp, indent);
3453 fprintf(fp, "lambda = %e\n", sh->lambda);
3454 pr_indent(fp, indent);
3455 fprintf(fp, "buffer size = %" PRId64 "\n", sh->sizeOfTprBody);