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, 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/fileio/filetypes.h"
53 #include "gromacs/fileio/gmxfio.h"
54 #include "gromacs/fileio/gmxfio_xdr.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdtypes/awh_history.h"
58 #include "gromacs/mdtypes/awh_params.h"
59 #include "gromacs/mdtypes/inputrec.h"
60 #include "gromacs/mdtypes/md_enums.h"
61 #include "gromacs/mdtypes/multipletimestepping.h"
62 #include "gromacs/mdtypes/pull_params.h"
63 #include "gromacs/mdtypes/state.h"
64 #include "gromacs/pbcutil/boxutilities.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/topology/block.h"
67 #include "gromacs/topology/ifunc.h"
68 #include "gromacs/topology/mtop_util.h"
69 #include "gromacs/topology/symtab.h"
70 #include "gromacs/topology/topology.h"
71 #include "gromacs/utility/arraysize.h"
72 #include "gromacs/utility/baseversion.h"
73 #include "gromacs/utility/cstringutil.h"
74 #include "gromacs/utility/fatalerror.h"
75 #include "gromacs/utility/futil.h"
76 #include "gromacs/utility/gmxassert.h"
77 #include "gromacs/utility/inmemoryserializer.h"
78 #include "gromacs/utility/keyvaluetreebuilder.h"
79 #include "gromacs/utility/keyvaluetreeserializer.h"
80 #include "gromacs/utility/smalloc.h"
81 #include "gromacs/utility/snprintf.h"
82 #include "gromacs/utility/txtdump.h"
84 #define TPX_TAG_RELEASE "release"
86 /*! \brief Tag string for the file format written to run input files
87 * written by this version of the code.
89 * Change this if you want to change the run input format in a feature
90 * branch. This ensures that there will not be different run input
91 * formats around which cannot be distinguished, while not causing
92 * problems rebasing the feature branch onto upstream changes. When
93 * merging with mainstream GROMACS, set this tag string back to
94 * TPX_TAG_RELEASE, and instead add an element to tpxv.
96 static const char* tpx_tag = TPX_TAG_RELEASE;
98 /*! \brief Enum of values that describe the contents of a tpr file
99 * whose format matches a version number
101 * The enum helps the code be more self-documenting and ensure merges
102 * do not silently resolve when two patches make the same bump. When
103 * adding new functionality, add a new element just above tpxv_Count
104 * in this enumeration, and write code below that does the right thing
105 * according to the value of file_version.
109 tpxv_ComputationalElectrophysiology =
110 96, /**< support for ion/water position swaps (computational electrophysiology) */
111 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to int64_t */
112 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
113 tpxv_InteractiveMolecularDynamics, /**< interactive molecular dynamics (IMD) */
114 tpxv_RemoveObsoleteParameters1, /**< remove optimize_fft, dihre_fc, nstcheckpoint */
115 tpxv_PullCoordTypeGeom, /**< add pull type and geometry per group and flat-bottom */
116 tpxv_PullGeomDirRel, /**< add pull geometry direction-relative */
117 tpxv_IntermolecularBondeds, /**< permit inter-molecular bonded interactions in the topology */
118 tpxv_CompElWithSwapLayerOffset, /**< added parameters for improved CompEl setups */
119 tpxv_CompElPolyatomicIonsAndMultipleIonTypes, /**< CompEl now can handle polyatomic ions and more than two types of ions */
120 tpxv_RemoveAdress, /**< removed support for AdResS */
121 tpxv_PullCoordNGroup, /**< add ngroup to pull coord */
122 tpxv_RemoveTwinRange, /**< removed support for twin-range interactions */
123 tpxv_ReplacePullPrintCOM12, /**< Replaced print-com-1, 2 with pull-print-com */
124 tpxv_PullExternalPotential, /**< Added pull type external potential */
125 tpxv_GenericParamsForElectricField, /**< Introduced KeyValueTree and moved electric field parameters */
126 tpxv_AcceleratedWeightHistogram, /**< sampling with accelerated weight histogram method (AWH) */
127 tpxv_RemoveImplicitSolvation, /**< removed support for implicit solvation */
128 tpxv_PullPrevStepCOMAsReference, /**< Enabled using the COM of the pull group of the last frame as reference for PBC */
129 tpxv_MimicQMMM, /**< Introduced support for MiMiC QM/MM interface */
130 tpxv_PullAverage, /**< Added possibility to output average pull force and position */
131 tpxv_GenericInternalParameters, /**< Added internal parameters for mdrun modules*/
132 tpxv_VSite2FD, /**< Added 2FD type virtual site */
133 tpxv_AddSizeField, /**< Added field with information about the size of the serialized tpr file in bytes, excluding the header */
134 tpxv_StoreNonBondedInteractionExclusionGroup, /**< Store the non bonded interaction exclusion group in the topology */
135 tpxv_VSite1, /**< Added 1 type virtual site */
136 tpxv_MTS, /**< Added multiple time stepping */
137 tpxv_Count /**< the total number of tpxv versions */
140 /*! \brief Version number of the file format written to run input
141 * files by this version of the code.
143 * The tpx_version increases whenever the file format in the main
144 * development branch changes, due to an extension of the tpxv enum above.
145 * Backward compatibility for reading old run input files is maintained
146 * by checking this version number against that of the file and then using
147 * the correct code path.
149 * When developing a feature branch that needs to change the run input
150 * file format, change tpx_tag instead. */
151 static const int tpx_version = tpxv_Count - 1;
154 /* This number should only be increased when you edit the TOPOLOGY section
155 * or the HEADER of the tpx format.
156 * This way we can maintain forward compatibility too for all analysis tools
157 * and/or external programs that only need to know the atom/residue names,
158 * charges, and bond connectivity.
160 * It first appeared in tpx version 26, when I also moved the inputrecord
161 * to the end of the tpx file, so we can just skip it if we only
164 * In particular, it must be increased when adding new elements to
165 * ftupd, so that old code can read new .tpr files.
167 * Updated for added field that contains the number of bytes of the tpr body, excluding the header.
169 static const int tpx_generation = 27;
171 /* This number should be the most recent backwards incompatible version
172 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
174 static const int tpx_incompatible_version = 57; // GMX4.0 has version 58
177 /* Struct used to maintain tpx compatibility when function types are added */
180 int fvnr; /* file version number in which the function type first appeared */
181 int ftype; /* function type */
185 * TODO The following three lines make little sense, please clarify if
186 * you've had to work out how ftupd works.
188 * The entries should be ordered in:
189 * 1. ascending function type number
190 * 2. ascending file version number
192 * Because we support reading of old .tpr file versions (even when
193 * mdrun can no longer run the simulation), we need to be able to read
194 * obsolete t_interaction_function types. Any data read from such
195 * fields is discarded. Their names have _NOLONGERUSED appended to
196 * them to make things clear.
198 static const t_ftupd ftupd[] = {
199 { 70, F_RESTRBONDS },
200 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
201 { 76, F_LINEAR_ANGLES },
202 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
203 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
205 { 60, F_GB12_NOLONGERUSED },
206 { 61, F_GB13_NOLONGERUSED },
207 { 61, F_GB14_NOLONGERUSED },
208 { 72, F_GBPOL_NOLONGERUSED },
209 { 72, F_NPSOLVATION_NOLONGERUSED },
211 { 76, F_ANHARM_POL },
213 { tpxv_VSite1, F_VSITE1 },
214 { tpxv_VSite2FD, F_VSITE2FD },
215 { tpxv_GenericInternalParameters, F_DENSITYFITTING },
216 { 69, F_VTEMP_NOLONGERUSED },
227 { 79, F_DVDL_RESTRAINT },
228 { 79, F_DVDL_TEMPERATURE },
230 #define NFTUPD asize(ftupd)
232 /* Needed for backward compatibility */
235 /**************************************************************
237 * Now the higer level routines that do io of the structures and arrays
239 **************************************************************/
240 static void do_pullgrp_tpx_pre95(gmx::ISerializer* serializer, t_pull_group* pgrp, t_pull_coord* pcrd)
244 int numAtoms = pgrp->ind.size();
245 serializer->doInt(&numAtoms);
246 pgrp->ind.resize(numAtoms);
247 serializer->doIntArray(pgrp->ind.data(), numAtoms);
248 int numWeights = pgrp->weight.size();
249 serializer->doInt(&numWeights);
250 pgrp->weight.resize(numWeights);
251 serializer->doRealArray(pgrp->weight.data(), numWeights);
252 serializer->doInt(&pgrp->pbcatom);
253 serializer->doRvec(&pcrd->vec.as_vec());
254 clear_rvec(pcrd->origin);
255 serializer->doRvec(&tmp);
257 serializer->doReal(&pcrd->rate);
258 serializer->doReal(&pcrd->k);
259 serializer->doReal(&pcrd->kB);
262 static void do_pull_group(gmx::ISerializer* serializer, t_pull_group* pgrp)
264 int numAtoms = pgrp->ind.size();
265 serializer->doInt(&numAtoms);
266 pgrp->ind.resize(numAtoms);
267 serializer->doIntArray(pgrp->ind.data(), numAtoms);
268 int numWeights = pgrp->weight.size();
269 serializer->doInt(&numWeights);
270 pgrp->weight.resize(numWeights);
271 serializer->doRealArray(pgrp->weight.data(), numWeights);
272 serializer->doInt(&pgrp->pbcatom);
275 static void do_pull_coord(gmx::ISerializer* serializer,
282 if (file_version >= tpxv_PullCoordNGroup)
284 serializer->doInt(&pcrd->eType);
285 if (file_version >= tpxv_PullExternalPotential)
287 if (pcrd->eType == epullEXTERNAL)
290 if (serializer->reading())
292 serializer->doString(&buf);
293 pcrd->externalPotentialProvider = gmx_strdup(buf.c_str());
297 buf = pcrd->externalPotentialProvider;
298 serializer->doString(&buf);
303 pcrd->externalPotentialProvider.clear();
308 if (serializer->reading())
310 pcrd->externalPotentialProvider.clear();
313 /* Note that we try to support adding new geometries without
314 * changing the tpx version. This requires checks when printing the
315 * geometry string and a check and fatal_error in init_pull.
317 serializer->doInt(&pcrd->eGeom);
318 serializer->doInt(&pcrd->ngroup);
319 if (pcrd->ngroup <= c_pullCoordNgroupMax)
321 serializer->doIntArray(pcrd->group.data(), pcrd->ngroup);
325 /* More groups in file than supported, this must be a new geometry
326 * that is not supported by our current code. Since we will not
327 * use the groups for this coord (checks in the pull and WHAM code
328 * ensure this), we can ignore the groups and set ngroup=0.
331 snew(dum, pcrd->ngroup);
332 serializer->doIntArray(dum, pcrd->ngroup);
337 serializer->doIvec(&pcrd->dim.as_vec());
342 serializer->doInt(&pcrd->group[0]);
343 serializer->doInt(&pcrd->group[1]);
344 if (file_version >= tpxv_PullCoordTypeGeom)
346 pcrd->ngroup = (pcrd->eGeom == epullgDIRRELATIVE ? 4 : 2);
347 serializer->doInt(&pcrd->eType);
348 serializer->doInt(&pcrd->eGeom);
349 if (pcrd->ngroup == 4)
351 serializer->doInt(&pcrd->group[2]);
352 serializer->doInt(&pcrd->group[3]);
354 serializer->doIvec(&pcrd->dim.as_vec());
358 pcrd->eType = ePullOld;
359 pcrd->eGeom = eGeomOld;
360 copy_ivec(dimOld, pcrd->dim);
363 serializer->doRvec(&pcrd->origin.as_vec());
364 serializer->doRvec(&pcrd->vec.as_vec());
365 if (file_version >= tpxv_PullCoordTypeGeom)
367 serializer->doBool(&pcrd->bStart);
371 /* This parameter is only printed, but not actually used by mdrun */
372 pcrd->bStart = FALSE;
374 serializer->doReal(&pcrd->init);
375 serializer->doReal(&pcrd->rate);
376 serializer->doReal(&pcrd->k);
377 serializer->doReal(&pcrd->kB);
380 static void do_expandedvals(gmx::ISerializer* serializer, t_expanded* expand, t_lambda* fepvals, int file_version)
382 int n_lambda = fepvals->n_lambda;
384 /* reset the lambda calculation window */
385 fepvals->lambda_start_n = 0;
386 fepvals->lambda_stop_n = n_lambda;
387 if (file_version >= 79)
391 if (serializer->reading())
393 snew(expand->init_lambda_weights, n_lambda);
395 serializer->doRealArray(expand->init_lambda_weights, n_lambda);
396 serializer->doBool(&expand->bInit_weights);
399 serializer->doInt(&expand->nstexpanded);
400 serializer->doInt(&expand->elmcmove);
401 serializer->doInt(&expand->elamstats);
402 serializer->doInt(&expand->lmc_repeats);
403 serializer->doInt(&expand->gibbsdeltalam);
404 serializer->doInt(&expand->lmc_forced_nstart);
405 serializer->doInt(&expand->lmc_seed);
406 serializer->doReal(&expand->mc_temp);
407 serializer->doBool(&expand->bSymmetrizedTMatrix);
408 serializer->doInt(&expand->nstTij);
409 serializer->doInt(&expand->minvarmin);
410 serializer->doInt(&expand->c_range);
411 serializer->doReal(&expand->wl_scale);
412 serializer->doReal(&expand->wl_ratio);
413 serializer->doReal(&expand->init_wl_delta);
414 serializer->doBool(&expand->bWLoneovert);
415 serializer->doInt(&expand->elmceq);
416 serializer->doInt(&expand->equil_steps);
417 serializer->doInt(&expand->equil_samples);
418 serializer->doInt(&expand->equil_n_at_lam);
419 serializer->doReal(&expand->equil_wl_delta);
420 serializer->doReal(&expand->equil_ratio);
424 static void do_simtempvals(gmx::ISerializer* serializer, t_simtemp* simtemp, int n_lambda, int file_version)
426 if (file_version >= 79)
428 serializer->doInt(&simtemp->eSimTempScale);
429 serializer->doReal(&simtemp->simtemp_high);
430 serializer->doReal(&simtemp->simtemp_low);
433 if (serializer->reading())
435 snew(simtemp->temperatures, n_lambda);
437 serializer->doRealArray(simtemp->temperatures, n_lambda);
442 static void do_imd(gmx::ISerializer* serializer, t_IMD* imd)
444 serializer->doInt(&imd->nat);
445 if (serializer->reading())
447 snew(imd->ind, imd->nat);
449 serializer->doIntArray(imd->ind, imd->nat);
452 static void do_fepvals(gmx::ISerializer* serializer, t_lambda* fepvals, int file_version)
454 /* i is defined in the ndo_double macro; use g to iterate. */
458 /* free energy values */
460 if (file_version >= 79)
462 serializer->doInt(&fepvals->init_fep_state);
463 serializer->doDouble(&fepvals->init_lambda);
464 serializer->doDouble(&fepvals->delta_lambda);
466 else if (file_version >= 59)
468 serializer->doDouble(&fepvals->init_lambda);
469 serializer->doDouble(&fepvals->delta_lambda);
473 serializer->doReal(&rdum);
474 fepvals->init_lambda = rdum;
475 serializer->doReal(&rdum);
476 fepvals->delta_lambda = rdum;
478 if (file_version >= 79)
480 serializer->doInt(&fepvals->n_lambda);
481 if (serializer->reading())
483 snew(fepvals->all_lambda, efptNR);
485 for (g = 0; g < efptNR; g++)
487 if (fepvals->n_lambda > 0)
489 if (serializer->reading())
491 snew(fepvals->all_lambda[g], fepvals->n_lambda);
493 serializer->doDoubleArray(fepvals->all_lambda[g], fepvals->n_lambda);
494 serializer->doBoolArray(fepvals->separate_dvdl, efptNR);
496 else if (fepvals->init_lambda >= 0)
498 fepvals->separate_dvdl[efptFEP] = TRUE;
502 else if (file_version >= 64)
504 serializer->doInt(&fepvals->n_lambda);
505 if (serializer->reading())
509 snew(fepvals->all_lambda, efptNR);
510 /* still allocate the all_lambda array's contents. */
511 for (g = 0; g < efptNR; g++)
513 if (fepvals->n_lambda > 0)
515 snew(fepvals->all_lambda[g], fepvals->n_lambda);
519 serializer->doDoubleArray(fepvals->all_lambda[efptFEP], fepvals->n_lambda);
520 if (fepvals->init_lambda >= 0)
524 fepvals->separate_dvdl[efptFEP] = TRUE;
526 if (serializer->reading())
528 /* copy the contents of the efptFEP lambda component to all
529 the other components */
530 for (g = 0; g < efptNR; g++)
532 for (h = 0; h < fepvals->n_lambda; h++)
536 fepvals->all_lambda[g][h] = fepvals->all_lambda[efptFEP][h];
545 fepvals->n_lambda = 0;
546 fepvals->all_lambda = nullptr;
547 if (fepvals->init_lambda >= 0)
549 fepvals->separate_dvdl[efptFEP] = TRUE;
552 serializer->doReal(&fepvals->sc_alpha);
553 serializer->doInt(&fepvals->sc_power);
554 if (file_version >= 79)
556 serializer->doReal(&fepvals->sc_r_power);
560 fepvals->sc_r_power = 6.0;
562 if (fepvals->sc_r_power != 6.0)
564 gmx_fatal(FARGS, "sc-r-power=48 is no longer supported");
566 serializer->doReal(&fepvals->sc_sigma);
567 if (serializer->reading())
569 if (file_version >= 71)
571 fepvals->sc_sigma_min = fepvals->sc_sigma;
575 fepvals->sc_sigma_min = 0;
578 if (file_version >= 79)
580 serializer->doBool(&fepvals->bScCoul);
584 fepvals->bScCoul = TRUE;
586 if (file_version >= 64)
588 serializer->doInt(&fepvals->nstdhdl);
592 fepvals->nstdhdl = 1;
595 if (file_version >= 73)
597 serializer->doInt(&fepvals->separate_dhdl_file);
598 serializer->doInt(&fepvals->dhdl_derivatives);
602 fepvals->separate_dhdl_file = esepdhdlfileYES;
603 fepvals->dhdl_derivatives = edhdlderivativesYES;
605 if (file_version >= 71)
607 serializer->doInt(&fepvals->dh_hist_size);
608 serializer->doDouble(&fepvals->dh_hist_spacing);
612 fepvals->dh_hist_size = 0;
613 fepvals->dh_hist_spacing = 0.1;
615 if (file_version >= 79)
617 serializer->doInt(&fepvals->edHdLPrintEnergy);
621 fepvals->edHdLPrintEnergy = edHdLPrintEnergyNO;
624 /* handle lambda_neighbors */
625 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
627 serializer->doInt(&fepvals->lambda_neighbors);
628 if ((fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0)
629 && (fepvals->init_lambda < 0))
631 fepvals->lambda_start_n = (fepvals->init_fep_state - fepvals->lambda_neighbors);
632 fepvals->lambda_stop_n = (fepvals->init_fep_state + fepvals->lambda_neighbors + 1);
633 if (fepvals->lambda_start_n < 0)
635 fepvals->lambda_start_n = 0;
637 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
639 fepvals->lambda_stop_n = fepvals->n_lambda;
644 fepvals->lambda_start_n = 0;
645 fepvals->lambda_stop_n = fepvals->n_lambda;
650 fepvals->lambda_start_n = 0;
651 fepvals->lambda_stop_n = fepvals->n_lambda;
655 static void do_awhBias(gmx::ISerializer* serializer, gmx::AwhBiasParams* awhBiasParams, int gmx_unused file_version)
657 serializer->doInt(&awhBiasParams->eTarget);
658 serializer->doDouble(&awhBiasParams->targetBetaScaling);
659 serializer->doDouble(&awhBiasParams->targetCutoff);
660 serializer->doInt(&awhBiasParams->eGrowth);
661 serializer->doInt(&awhBiasParams->bUserData);
662 serializer->doDouble(&awhBiasParams->errorInitial);
663 serializer->doInt(&awhBiasParams->ndim);
664 serializer->doInt(&awhBiasParams->shareGroup);
665 serializer->doBool(&awhBiasParams->equilibrateHistogram);
667 if (serializer->reading())
669 snew(awhBiasParams->dimParams, awhBiasParams->ndim);
672 for (int d = 0; d < awhBiasParams->ndim; d++)
674 gmx::AwhDimParams* dimParams = &awhBiasParams->dimParams[d];
676 serializer->doInt(&dimParams->eCoordProvider);
677 serializer->doInt(&dimParams->coordIndex);
678 serializer->doDouble(&dimParams->origin);
679 serializer->doDouble(&dimParams->end);
680 serializer->doDouble(&dimParams->period);
681 serializer->doDouble(&dimParams->forceConstant);
682 serializer->doDouble(&dimParams->diffusion);
683 serializer->doDouble(&dimParams->coordValueInit);
684 serializer->doDouble(&dimParams->coverDiameter);
688 static void do_awh(gmx::ISerializer* serializer, gmx::AwhParams* awhParams, int gmx_unused file_version)
690 serializer->doInt(&awhParams->numBias);
691 serializer->doInt(&awhParams->nstOut);
692 serializer->doInt64(&awhParams->seed);
693 serializer->doInt(&awhParams->nstSampleCoord);
694 serializer->doInt(&awhParams->numSamplesUpdateFreeEnergy);
695 serializer->doInt(&awhParams->ePotential);
696 serializer->doBool(&awhParams->shareBiasMultisim);
698 if (awhParams->numBias > 0)
700 if (serializer->reading())
702 snew(awhParams->awhBiasParams, awhParams->numBias);
705 for (int k = 0; k < awhParams->numBias; k++)
707 do_awhBias(serializer, &awhParams->awhBiasParams[k], file_version);
712 static void do_pull(gmx::ISerializer* serializer, pull_params_t* pull, int file_version, int ePullOld)
718 if (file_version >= 95)
720 serializer->doInt(&pull->ngroup);
722 serializer->doInt(&pull->ncoord);
723 if (file_version < 95)
725 pull->ngroup = pull->ncoord + 1;
727 if (file_version < tpxv_PullCoordTypeGeom)
731 serializer->doInt(&eGeomOld);
732 serializer->doIvec(&dimOld);
733 /* The inner cylinder radius, now removed */
734 serializer->doReal(&dum);
736 serializer->doReal(&pull->cylinder_r);
737 serializer->doReal(&pull->constr_tol);
738 if (file_version >= 95)
740 serializer->doBool(&pull->bPrintCOM);
741 /* With file_version < 95 this value is set below */
743 if (file_version >= tpxv_ReplacePullPrintCOM12)
745 serializer->doBool(&pull->bPrintRefValue);
746 serializer->doBool(&pull->bPrintComp);
748 else if (file_version >= tpxv_PullCoordTypeGeom)
751 serializer->doInt(&idum); /* used to be bPrintCOM2 */
752 serializer->doBool(&pull->bPrintRefValue);
753 serializer->doBool(&pull->bPrintComp);
757 pull->bPrintRefValue = FALSE;
758 pull->bPrintComp = TRUE;
760 serializer->doInt(&pull->nstxout);
761 serializer->doInt(&pull->nstfout);
762 if (file_version >= tpxv_PullPrevStepCOMAsReference)
764 serializer->doBool(&pull->bSetPbcRefToPrevStepCOM);
768 pull->bSetPbcRefToPrevStepCOM = FALSE;
770 pull->group.resize(pull->ngroup);
771 pull->coord.resize(pull->ncoord);
772 if (file_version < 95)
774 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
775 if (eGeomOld == epullgDIRPBC)
777 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
779 if (eGeomOld > epullgDIRPBC)
784 for (g = 0; g < pull->ngroup; g++)
786 /* We read and ignore a pull coordinate for group 0 */
787 do_pullgrp_tpx_pre95(serializer, &pull->group[g], &pull->coord[std::max(g - 1, 0)]);
790 pull->coord[g - 1].group[0] = 0;
791 pull->coord[g - 1].group[1] = g;
795 pull->bPrintCOM = (!pull->group[0].ind.empty());
799 for (g = 0; g < pull->ngroup; g++)
801 do_pull_group(serializer, &pull->group[g]);
803 for (g = 0; g < pull->ncoord; g++)
805 do_pull_coord(serializer, &pull->coord[g], file_version, ePullOld, eGeomOld, dimOld);
808 if (file_version >= tpxv_PullAverage)
812 v = pull->bXOutAverage;
813 serializer->doBool(&v);
814 pull->bXOutAverage = v;
815 v = pull->bFOutAverage;
816 serializer->doBool(&v);
817 pull->bFOutAverage = v;
822 static void do_rotgrp(gmx::ISerializer* serializer, t_rotgrp* rotg)
824 serializer->doInt(&rotg->eType);
825 serializer->doInt(&rotg->bMassW);
826 serializer->doInt(&rotg->nat);
827 if (serializer->reading())
829 snew(rotg->ind, rotg->nat);
831 serializer->doIntArray(rotg->ind, rotg->nat);
832 if (serializer->reading())
834 snew(rotg->x_ref, rotg->nat);
836 serializer->doRvecArray(rotg->x_ref, rotg->nat);
837 serializer->doRvec(&rotg->inputVec);
838 serializer->doRvec(&rotg->pivot);
839 serializer->doReal(&rotg->rate);
840 serializer->doReal(&rotg->k);
841 serializer->doReal(&rotg->slab_dist);
842 serializer->doReal(&rotg->min_gaussian);
843 serializer->doReal(&rotg->eps);
844 serializer->doInt(&rotg->eFittype);
845 serializer->doInt(&rotg->PotAngle_nstep);
846 serializer->doReal(&rotg->PotAngle_step);
849 static void do_rot(gmx::ISerializer* serializer, t_rot* rot)
853 serializer->doInt(&rot->ngrp);
854 serializer->doInt(&rot->nstrout);
855 serializer->doInt(&rot->nstsout);
856 if (serializer->reading())
858 snew(rot->grp, rot->ngrp);
860 for (g = 0; g < rot->ngrp; g++)
862 do_rotgrp(serializer, &rot->grp[g]);
867 static void do_swapgroup(gmx::ISerializer* serializer, t_swapGroup* g)
870 /* Name of the group or molecule */
872 if (serializer->reading())
874 serializer->doString(&buf);
875 g->molname = gmx_strdup(buf.c_str());
880 serializer->doString(&buf);
883 /* Number of atoms in the group */
884 serializer->doInt(&g->nat);
886 /* The group's atom indices */
887 if (serializer->reading())
889 snew(g->ind, g->nat);
891 serializer->doIntArray(g->ind, g->nat);
893 /* Requested counts for compartments A and B */
894 serializer->doIntArray(g->nmolReq, eCompNR);
897 static void do_swapcoords_tpx(gmx::ISerializer* serializer, t_swapcoords* swap, int file_version)
899 /* Enums for better readability of the code */
912 if (file_version >= tpxv_CompElPolyatomicIonsAndMultipleIonTypes)
914 /* The total number of swap groups is the sum of the fixed groups
915 * (split0, split1, solvent), and the user-defined groups (2+ types of ions)
917 serializer->doInt(&swap->ngrp);
918 if (serializer->reading())
920 snew(swap->grp, swap->ngrp);
922 for (int ig = 0; ig < swap->ngrp; ig++)
924 do_swapgroup(serializer, &swap->grp[ig]);
926 serializer->doBool(&swap->massw_split[eChannel0]);
927 serializer->doBool(&swap->massw_split[eChannel1]);
928 serializer->doInt(&swap->nstswap);
929 serializer->doInt(&swap->nAverage);
930 serializer->doReal(&swap->threshold);
931 serializer->doReal(&swap->cyl0r);
932 serializer->doReal(&swap->cyl0u);
933 serializer->doReal(&swap->cyl0l);
934 serializer->doReal(&swap->cyl1r);
935 serializer->doReal(&swap->cyl1u);
936 serializer->doReal(&swap->cyl1l);
940 /*** Support reading older CompEl .tpr files ***/
942 /* In the original CompEl .tpr files, we always have 5 groups: */
944 snew(swap->grp, swap->ngrp);
946 swap->grp[eGrpSplit0].molname = gmx_strdup("split0"); // group 0: split0
947 swap->grp[eGrpSplit1].molname = gmx_strdup("split1"); // group 1: split1
948 swap->grp[eGrpSolvent].molname = gmx_strdup("solvent"); // group 2: solvent
949 swap->grp[3].molname = gmx_strdup("anions"); // group 3: anions
950 swap->grp[4].molname = gmx_strdup("cations"); // group 4: cations
952 serializer->doInt(&swap->grp[3].nat);
953 serializer->doInt(&swap->grp[eGrpSolvent].nat);
954 serializer->doInt(&swap->grp[eGrpSplit0].nat);
955 serializer->doBool(&swap->massw_split[eChannel0]);
956 serializer->doInt(&swap->grp[eGrpSplit1].nat);
957 serializer->doBool(&swap->massw_split[eChannel1]);
958 serializer->doInt(&swap->nstswap);
959 serializer->doInt(&swap->nAverage);
960 serializer->doReal(&swap->threshold);
961 serializer->doReal(&swap->cyl0r);
962 serializer->doReal(&swap->cyl0u);
963 serializer->doReal(&swap->cyl0l);
964 serializer->doReal(&swap->cyl1r);
965 serializer->doReal(&swap->cyl1u);
966 serializer->doReal(&swap->cyl1l);
968 // The order[] array keeps compatibility with older .tpr files
969 // by reading in the groups in the classic order
971 const int order[4] = { 3, eGrpSolvent, eGrpSplit0, eGrpSplit1 };
973 for (int ig = 0; ig < 4; ig++)
976 snew(swap->grp[g].ind, swap->grp[g].nat);
977 serializer->doIntArray(swap->grp[g].ind, swap->grp[g].nat);
981 for (int j = eCompA; j <= eCompB; j++)
983 serializer->doInt(&swap->grp[3].nmolReq[j]); // group 3 = anions
984 serializer->doInt(&swap->grp[4].nmolReq[j]); // group 4 = cations
986 } /* End support reading older CompEl .tpr files */
988 if (file_version >= tpxv_CompElWithSwapLayerOffset)
990 serializer->doReal(&swap->bulkOffset[eCompA]);
991 serializer->doReal(&swap->bulkOffset[eCompB]);
995 static void do_legacy_efield(gmx::ISerializer* serializer, gmx::KeyValueTreeObjectBuilder* root)
997 const char* const dimName[] = { "x", "y", "z" };
999 auto appliedForcesObj = root->addObject("applied-forces");
1000 auto efieldObj = appliedForcesObj.addObject("electric-field");
1001 // The content of the tpr file for this feature has
1002 // been the same since gromacs 4.0 that was used for
1004 for (int j = 0; j < DIM; ++j)
1007 serializer->doInt(&n);
1008 serializer->doInt(&nt);
1009 std::vector<real> aa(n + 1), phi(nt + 1), at(nt + 1), phit(nt + 1);
1010 serializer->doRealArray(aa.data(), n);
1011 serializer->doRealArray(phi.data(), n);
1012 serializer->doRealArray(at.data(), nt);
1013 serializer->doRealArray(phit.data(), nt);
1016 if (n > 1 || nt > 1)
1019 "Can not handle tpr files with more than one electric field term per "
1022 auto dimObj = efieldObj.addObject(dimName[j]);
1023 dimObj.addValue<real>("E0", aa[0]);
1024 dimObj.addValue<real>("omega", at[0]);
1025 dimObj.addValue<real>("t0", phi[0]);
1026 dimObj.addValue<real>("sigma", phit[0]);
1032 static void do_inputrec(gmx::ISerializer* serializer, t_inputrec* ir, int file_version)
1034 int i, j, k, idum = 0;
1036 gmx_bool bdum = false;
1038 if (file_version != tpx_version)
1040 /* Give a warning about features that are not accessible */
1041 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n", file_version, tpx_version);
1044 if (file_version == 0)
1049 gmx::KeyValueTreeBuilder paramsBuilder;
1050 gmx::KeyValueTreeObjectBuilder paramsObj = paramsBuilder.rootObject();
1052 /* Basic inputrec stuff */
1053 serializer->doInt(&ir->eI);
1054 if (file_version >= 62)
1056 serializer->doInt64(&ir->nsteps);
1060 serializer->doInt(&idum);
1064 if (file_version >= 62)
1066 serializer->doInt64(&ir->init_step);
1070 serializer->doInt(&idum);
1071 ir->init_step = idum;
1074 serializer->doInt(&ir->simulation_part);
1076 if (file_version >= tpxv_MTS)
1078 serializer->doBool(&ir->useMts);
1079 int numLevels = ir->mtsLevels.size();
1082 serializer->doInt(&numLevels);
1084 ir->mtsLevels.resize(numLevels);
1085 for (auto& mtsLevel : ir->mtsLevels)
1087 int forceGroups = mtsLevel.forceGroups.to_ulong();
1088 serializer->doInt(&forceGroups);
1089 mtsLevel.forceGroups = std::bitset<static_cast<int>(gmx::MtsForceGroups::Count)>(forceGroups);
1090 serializer->doInt(&mtsLevel.stepFactor);
1096 ir->mtsLevels.clear();
1099 if (file_version >= 67)
1101 serializer->doInt(&ir->nstcalcenergy);
1105 ir->nstcalcenergy = 1;
1107 if (file_version >= 81)
1109 serializer->doInt(&ir->cutoff_scheme);
1110 if (file_version < 94)
1112 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
1117 ir->cutoff_scheme = ecutsGROUP;
1119 serializer->doInt(&idum); /* used to be ns_type; not used anymore */
1120 serializer->doInt(&ir->nstlist);
1121 serializer->doInt(&idum); /* used to be ndelta; not used anymore */
1123 serializer->doReal(&ir->rtpi);
1125 serializer->doInt(&ir->nstcomm);
1126 serializer->doInt(&ir->comm_mode);
1128 /* ignore nstcheckpoint */
1129 if (file_version < tpxv_RemoveObsoleteParameters1)
1131 serializer->doInt(&idum);
1134 serializer->doInt(&ir->nstcgsteep);
1136 serializer->doInt(&ir->nbfgscorr);
1138 serializer->doInt(&ir->nstlog);
1139 serializer->doInt(&ir->nstxout);
1140 serializer->doInt(&ir->nstvout);
1141 serializer->doInt(&ir->nstfout);
1142 serializer->doInt(&ir->nstenergy);
1143 serializer->doInt(&ir->nstxout_compressed);
1144 if (file_version >= 59)
1146 serializer->doDouble(&ir->init_t);
1147 serializer->doDouble(&ir->delta_t);
1151 serializer->doReal(&rdum);
1153 serializer->doReal(&rdum);
1156 serializer->doReal(&ir->x_compression_precision);
1157 if (file_version >= 81)
1159 serializer->doReal(&ir->verletbuf_tol);
1163 ir->verletbuf_tol = 0;
1165 serializer->doReal(&ir->rlist);
1166 if (file_version >= 67 && file_version < tpxv_RemoveTwinRange)
1168 if (serializer->reading())
1170 // Reading such a file version could invoke the twin-range
1171 // scheme, about which mdrun should give a fatal error.
1172 real dummy_rlistlong = -1;
1173 serializer->doReal(&dummy_rlistlong);
1175 ir->useTwinRange = (ir->rlist > 0 && (dummy_rlistlong == 0 || dummy_rlistlong > ir->rlist));
1176 // When true, this forces mdrun to issue an error (regardless of
1177 // ir->cutoff_scheme).
1179 // Otherwise, grompp used to set rlistlong actively. Users
1180 // were probably also confused and set rlistlong == rlist.
1181 // However, in all remaining cases, it is safe to let
1182 // mdrun proceed normally.
1187 // No need to read or write anything
1188 ir->useTwinRange = false;
1190 if (file_version >= 82 && file_version != 90)
1192 // Multiple time-stepping is no longer enabled, but the old
1193 // support required the twin-range scheme, for which mdrun
1194 // already emits a fatal error.
1195 int dummy_nstcalclr = -1;
1196 serializer->doInt(&dummy_nstcalclr);
1198 serializer->doInt(&ir->coulombtype);
1199 if (file_version >= 81)
1201 serializer->doInt(&ir->coulomb_modifier);
1205 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1207 serializer->doReal(&ir->rcoulomb_switch);
1208 serializer->doReal(&ir->rcoulomb);
1209 serializer->doInt(&ir->vdwtype);
1210 if (file_version >= 81)
1212 serializer->doInt(&ir->vdw_modifier);
1216 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1218 serializer->doReal(&ir->rvdw_switch);
1219 serializer->doReal(&ir->rvdw);
1220 serializer->doInt(&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->doInt(&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->doInt(&ir->ljpme_combination_rule);
1292 serializer->doBool(&ir->bContinuation);
1293 serializer->doInt(&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->doInt(&ir->epc);
1311 serializer->doInt(&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->doInt(&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->doInt(&ir->efep);
1343 do_fepvals(serializer, ir->fepvals, file_version);
1345 if (file_version >= 79)
1347 serializer->doBool(&ir->bSimTemp);
1350 ir->bSimTemp = TRUE;
1355 ir->bSimTemp = FALSE;
1359 do_simtempvals(serializer, ir->simtempvals, ir->fepvals->n_lambda, file_version);
1362 if (file_version >= 79)
1364 serializer->doBool(&ir->bExpanded);
1368 do_expandedvals(serializer, ir->expandedvals, ir->fepvals, file_version);
1371 serializer->doInt(&ir->eDisre);
1372 serializer->doInt(&ir->eDisreWeighting);
1373 serializer->doBool(&ir->bDisreMixed);
1374 serializer->doReal(&ir->dr_fc);
1375 serializer->doReal(&ir->dr_tau);
1376 serializer->doInt(&ir->nstdisreout);
1377 serializer->doReal(&ir->orires_fc);
1378 serializer->doReal(&ir->orires_tau);
1379 serializer->doInt(&ir->nstorireout);
1381 /* ignore dihre_fc */
1382 if (file_version < 79)
1384 serializer->doReal(&rdum);
1387 serializer->doReal(&ir->em_stepsize);
1388 serializer->doReal(&ir->em_tol);
1389 serializer->doBool(&ir->bShakeSOR);
1390 serializer->doInt(&ir->niter);
1391 serializer->doReal(&ir->fc_stepsize);
1392 serializer->doInt(&ir->eConstrAlg);
1393 serializer->doInt(&ir->nProjOrder);
1394 serializer->doReal(&ir->LincsWarnAngle);
1395 serializer->doInt(&ir->nLincsIter);
1396 serializer->doReal(&ir->bd_fric);
1397 if (file_version >= tpxv_Use64BitRandomSeed)
1399 serializer->doInt64(&ir->ld_seed);
1403 serializer->doInt(&idum);
1407 for (i = 0; i < DIM; i++)
1409 serializer->doRvec(&ir->deform[i]);
1411 serializer->doReal(&ir->cos_accel);
1413 serializer->doInt(&ir->userint1);
1414 serializer->doInt(&ir->userint2);
1415 serializer->doInt(&ir->userint3);
1416 serializer->doInt(&ir->userint4);
1417 serializer->doReal(&ir->userreal1);
1418 serializer->doReal(&ir->userreal2);
1419 serializer->doReal(&ir->userreal3);
1420 serializer->doReal(&ir->userreal4);
1422 /* AdResS is removed, but we need to be able to read old files,
1423 and let mdrun refuse to run them */
1424 if (file_version >= 77 && file_version < tpxv_RemoveAdress)
1426 serializer->doBool(&ir->bAdress);
1429 int idum, numThermoForceGroups, numEnergyGroups;
1432 serializer->doInt(&idum);
1433 serializer->doReal(&rdum);
1434 serializer->doReal(&rdum);
1435 serializer->doReal(&rdum);
1436 serializer->doInt(&idum);
1437 serializer->doInt(&idum);
1438 serializer->doRvec(&rvecdum);
1439 serializer->doInt(&numThermoForceGroups);
1440 serializer->doReal(&rdum);
1441 serializer->doInt(&numEnergyGroups);
1442 serializer->doInt(&idum);
1444 if (numThermoForceGroups > 0)
1446 std::vector<int> idumn(numThermoForceGroups);
1447 serializer->doIntArray(idumn.data(), idumn.size());
1449 if (numEnergyGroups > 0)
1451 std::vector<int> idumn(numEnergyGroups);
1452 serializer->doIntArray(idumn.data(), idumn.size());
1458 ir->bAdress = FALSE;
1465 if (file_version >= tpxv_PullCoordTypeGeom)
1467 serializer->doBool(&ir->bPull);
1471 serializer->doInt(&ePullOld);
1472 ir->bPull = (ePullOld > 0);
1473 /* We removed the first ePull=ePullNo for the enum */
1478 if (serializer->reading())
1480 ir->pull = std::make_unique<pull_params_t>();
1482 do_pull(serializer, ir->pull.get(), file_version, ePullOld);
1486 if (file_version >= tpxv_AcceleratedWeightHistogram)
1488 serializer->doBool(&ir->bDoAwh);
1492 if (serializer->reading())
1494 snew(ir->awhParams, 1);
1496 do_awh(serializer, ir->awhParams, file_version);
1504 /* Enforced rotation */
1505 if (file_version >= 74)
1507 serializer->doBool(&ir->bRot);
1510 if (serializer->reading())
1514 do_rot(serializer, ir->rot);
1522 /* Interactive molecular dynamics */
1523 if (file_version >= tpxv_InteractiveMolecularDynamics)
1525 serializer->doBool(&ir->bIMD);
1528 if (serializer->reading())
1532 do_imd(serializer, ir->imd);
1537 /* We don't support IMD sessions for old .tpr files */
1542 serializer->doInt(&ir->opts.ngtc);
1543 if (file_version >= 69)
1545 serializer->doInt(&ir->opts.nhchainlength);
1549 ir->opts.nhchainlength = 1;
1551 serializer->doInt(&ir->opts.ngacc);
1552 serializer->doInt(&ir->opts.ngfrz);
1553 serializer->doInt(&ir->opts.ngener);
1555 if (serializer->reading())
1557 snew(ir->opts.nrdf, ir->opts.ngtc);
1558 snew(ir->opts.ref_t, ir->opts.ngtc);
1559 snew(ir->opts.annealing, ir->opts.ngtc);
1560 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1561 snew(ir->opts.anneal_time, ir->opts.ngtc);
1562 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1563 snew(ir->opts.tau_t, ir->opts.ngtc);
1564 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1565 snew(ir->opts.acc, ir->opts.ngacc);
1566 snew(ir->opts.egp_flags, ir->opts.ngener * ir->opts.ngener);
1568 if (ir->opts.ngtc > 0)
1570 serializer->doRealArray(ir->opts.nrdf, ir->opts.ngtc);
1571 serializer->doRealArray(ir->opts.ref_t, ir->opts.ngtc);
1572 serializer->doRealArray(ir->opts.tau_t, ir->opts.ngtc);
1574 if (ir->opts.ngfrz > 0)
1576 serializer->doIvecArray(ir->opts.nFreeze, ir->opts.ngfrz);
1578 if (ir->opts.ngacc > 0)
1580 serializer->doRvecArray(ir->opts.acc, ir->opts.ngacc);
1582 serializer->doIntArray(ir->opts.egp_flags, ir->opts.ngener * ir->opts.ngener);
1584 /* First read the lists with annealing and npoints for each group */
1585 serializer->doIntArray(ir->opts.annealing, ir->opts.ngtc);
1586 serializer->doIntArray(ir->opts.anneal_npoints, ir->opts.ngtc);
1587 for (j = 0; j < (ir->opts.ngtc); j++)
1589 k = ir->opts.anneal_npoints[j];
1590 if (serializer->reading())
1592 snew(ir->opts.anneal_time[j], k);
1593 snew(ir->opts.anneal_temp[j], k);
1595 serializer->doRealArray(ir->opts.anneal_time[j], k);
1596 serializer->doRealArray(ir->opts.anneal_temp[j], k);
1600 serializer->doInt(&ir->nwall);
1601 serializer->doInt(&ir->wall_type);
1602 serializer->doReal(&ir->wall_r_linpot);
1603 serializer->doInt(&ir->wall_atomtype[0]);
1604 serializer->doInt(&ir->wall_atomtype[1]);
1605 serializer->doReal(&ir->wall_density[0]);
1606 serializer->doReal(&ir->wall_density[1]);
1607 serializer->doReal(&ir->wall_ewald_zfac);
1610 /* Cosine stuff for electric fields */
1611 if (file_version < tpxv_GenericParamsForElectricField)
1613 do_legacy_efield(serializer, ¶msObj);
1617 if (file_version >= tpxv_ComputationalElectrophysiology)
1619 serializer->doInt(&ir->eSwapCoords);
1620 if (ir->eSwapCoords != eswapNO)
1622 if (serializer->reading())
1626 do_swapcoords_tpx(serializer, ir->swap, file_version);
1630 /* QMMM reading - despite defunct we require reading for backwards
1631 * compability and correct serialization
1634 serializer->doBool(&ir->bQMMM);
1636 serializer->doInt(&qmmmScheme);
1637 real unusedScalefactor;
1638 serializer->doReal(&unusedScalefactor);
1640 // this is still used in Mimic
1641 serializer->doInt(&ir->opts.ngQM);
1642 if (ir->opts.ngQM > 0 && ir->bQMMM)
1644 /* We leave in dummy i/o for removed parameters to avoid
1645 * changing the tpr format.
1647 std::vector<int> dummyIntVec(4 * ir->opts.ngQM, 0);
1648 serializer->doIntArray(dummyIntVec.data(), dummyIntVec.size());
1649 dummyIntVec.clear();
1651 // std::vector<bool> has no data()
1652 std::vector<char> dummyBoolVec(ir->opts.ngQM * sizeof(bool) / sizeof(char));
1653 serializer->doBoolArray(reinterpret_cast<bool*>(dummyBoolVec.data()), dummyBoolVec.size());
1654 dummyBoolVec.clear();
1656 dummyIntVec.resize(2 * ir->opts.ngQM, 0);
1657 serializer->doIntArray(dummyIntVec.data(), dummyIntVec.size());
1658 dummyIntVec.clear();
1660 std::vector<real> dummyRealVec(2 * ir->opts.ngQM, 0);
1661 serializer->doRealArray(dummyRealVec.data(), dummyRealVec.size());
1662 dummyRealVec.clear();
1664 dummyIntVec.resize(3 * ir->opts.ngQM, 0);
1665 serializer->doIntArray(dummyIntVec.data(), dummyIntVec.size());
1666 dummyIntVec.clear();
1668 /* end of QMMM stuff */
1671 if (file_version >= tpxv_GenericParamsForElectricField)
1673 if (serializer->reading())
1675 paramsObj.mergeObject(gmx::deserializeKeyValueTree(serializer));
1679 GMX_RELEASE_ASSERT(ir->params != nullptr,
1680 "Parameters should be present when writing inputrec");
1681 gmx::serializeKeyValueTree(*ir->params, serializer);
1684 if (serializer->reading())
1686 ir->params = new gmx::KeyValueTreeObject(paramsBuilder.build());
1687 // Initialize internal parameters to an empty kvt for all tpr versions
1688 ir->internalParameters = std::make_unique<gmx::KeyValueTreeObject>();
1691 if (file_version >= tpxv_GenericInternalParameters)
1693 if (serializer->reading())
1695 ir->internalParameters =
1696 std::make_unique<gmx::KeyValueTreeObject>(gmx::deserializeKeyValueTree(serializer));
1700 GMX_RELEASE_ASSERT(ir->internalParameters != nullptr,
1701 "Parameters should be present when writing inputrec");
1702 gmx::serializeKeyValueTree(*ir->internalParameters, serializer);
1708 static void do_harm(gmx::ISerializer* serializer, t_iparams* iparams)
1710 serializer->doReal(&iparams->harmonic.rA);
1711 serializer->doReal(&iparams->harmonic.krA);
1712 serializer->doReal(&iparams->harmonic.rB);
1713 serializer->doReal(&iparams->harmonic.krB);
1716 static void do_iparams(gmx::ISerializer* serializer, t_functype ftype, t_iparams* iparams, int file_version)
1729 do_harm(serializer, iparams);
1730 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && serializer->reading())
1732 /* Correct incorrect storage of parameters */
1733 iparams->pdihs.phiB = iparams->pdihs.phiA;
1734 iparams->pdihs.cpB = iparams->pdihs.cpA;
1738 serializer->doReal(&iparams->harmonic.rA);
1739 serializer->doReal(&iparams->harmonic.krA);
1741 case F_LINEAR_ANGLES:
1742 serializer->doReal(&iparams->linangle.klinA);
1743 serializer->doReal(&iparams->linangle.aA);
1744 serializer->doReal(&iparams->linangle.klinB);
1745 serializer->doReal(&iparams->linangle.aB);
1748 serializer->doReal(&iparams->fene.bm);
1749 serializer->doReal(&iparams->fene.kb);
1753 serializer->doReal(&iparams->restraint.lowA);
1754 serializer->doReal(&iparams->restraint.up1A);
1755 serializer->doReal(&iparams->restraint.up2A);
1756 serializer->doReal(&iparams->restraint.kA);
1757 serializer->doReal(&iparams->restraint.lowB);
1758 serializer->doReal(&iparams->restraint.up1B);
1759 serializer->doReal(&iparams->restraint.up2B);
1760 serializer->doReal(&iparams->restraint.kB);
1766 serializer->doReal(&iparams->tab.kA);
1767 serializer->doInt(&iparams->tab.table);
1768 serializer->doReal(&iparams->tab.kB);
1770 case F_CROSS_BOND_BONDS:
1771 serializer->doReal(&iparams->cross_bb.r1e);
1772 serializer->doReal(&iparams->cross_bb.r2e);
1773 serializer->doReal(&iparams->cross_bb.krr);
1775 case F_CROSS_BOND_ANGLES:
1776 serializer->doReal(&iparams->cross_ba.r1e);
1777 serializer->doReal(&iparams->cross_ba.r2e);
1778 serializer->doReal(&iparams->cross_ba.r3e);
1779 serializer->doReal(&iparams->cross_ba.krt);
1781 case F_UREY_BRADLEY:
1782 serializer->doReal(&iparams->u_b.thetaA);
1783 serializer->doReal(&iparams->u_b.kthetaA);
1784 serializer->doReal(&iparams->u_b.r13A);
1785 serializer->doReal(&iparams->u_b.kUBA);
1786 if (file_version >= 79)
1788 serializer->doReal(&iparams->u_b.thetaB);
1789 serializer->doReal(&iparams->u_b.kthetaB);
1790 serializer->doReal(&iparams->u_b.r13B);
1791 serializer->doReal(&iparams->u_b.kUBB);
1795 iparams->u_b.thetaB = iparams->u_b.thetaA;
1796 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1797 iparams->u_b.r13B = iparams->u_b.r13A;
1798 iparams->u_b.kUBB = iparams->u_b.kUBA;
1801 case F_QUARTIC_ANGLES:
1802 serializer->doReal(&iparams->qangle.theta);
1803 serializer->doRealArray(iparams->qangle.c, 5);
1806 serializer->doReal(&iparams->bham.a);
1807 serializer->doReal(&iparams->bham.b);
1808 serializer->doReal(&iparams->bham.c);
1811 serializer->doReal(&iparams->morse.b0A);
1812 serializer->doReal(&iparams->morse.cbA);
1813 serializer->doReal(&iparams->morse.betaA);
1814 if (file_version >= 79)
1816 serializer->doReal(&iparams->morse.b0B);
1817 serializer->doReal(&iparams->morse.cbB);
1818 serializer->doReal(&iparams->morse.betaB);
1822 iparams->morse.b0B = iparams->morse.b0A;
1823 iparams->morse.cbB = iparams->morse.cbA;
1824 iparams->morse.betaB = iparams->morse.betaA;
1828 serializer->doReal(&iparams->cubic.b0);
1829 serializer->doReal(&iparams->cubic.kb);
1830 serializer->doReal(&iparams->cubic.kcub);
1832 case F_CONNBONDS: break;
1833 case F_POLARIZATION: serializer->doReal(&iparams->polarize.alpha); break;
1835 serializer->doReal(&iparams->anharm_polarize.alpha);
1836 serializer->doReal(&iparams->anharm_polarize.drcut);
1837 serializer->doReal(&iparams->anharm_polarize.khyp);
1840 serializer->doReal(&iparams->wpol.al_x);
1841 serializer->doReal(&iparams->wpol.al_y);
1842 serializer->doReal(&iparams->wpol.al_z);
1843 serializer->doReal(&iparams->wpol.rOH);
1844 serializer->doReal(&iparams->wpol.rHH);
1845 serializer->doReal(&iparams->wpol.rOD);
1848 serializer->doReal(&iparams->thole.a);
1849 serializer->doReal(&iparams->thole.alpha1);
1850 serializer->doReal(&iparams->thole.alpha2);
1851 serializer->doReal(&iparams->thole.rfac);
1854 serializer->doReal(&iparams->lj.c6);
1855 serializer->doReal(&iparams->lj.c12);
1858 serializer->doReal(&iparams->lj14.c6A);
1859 serializer->doReal(&iparams->lj14.c12A);
1860 serializer->doReal(&iparams->lj14.c6B);
1861 serializer->doReal(&iparams->lj14.c12B);
1864 serializer->doReal(&iparams->ljc14.fqq);
1865 serializer->doReal(&iparams->ljc14.qi);
1866 serializer->doReal(&iparams->ljc14.qj);
1867 serializer->doReal(&iparams->ljc14.c6);
1868 serializer->doReal(&iparams->ljc14.c12);
1870 case F_LJC_PAIRS_NB:
1871 serializer->doReal(&iparams->ljcnb.qi);
1872 serializer->doReal(&iparams->ljcnb.qj);
1873 serializer->doReal(&iparams->ljcnb.c6);
1874 serializer->doReal(&iparams->ljcnb.c12);
1880 serializer->doReal(&iparams->pdihs.phiA);
1881 serializer->doReal(&iparams->pdihs.cpA);
1882 serializer->doReal(&iparams->pdihs.phiB);
1883 serializer->doReal(&iparams->pdihs.cpB);
1884 serializer->doInt(&iparams->pdihs.mult);
1887 serializer->doReal(&iparams->pdihs.phiA);
1888 serializer->doReal(&iparams->pdihs.cpA);
1891 serializer->doInt(&iparams->disres.label);
1892 serializer->doInt(&iparams->disres.type);
1893 serializer->doReal(&iparams->disres.low);
1894 serializer->doReal(&iparams->disres.up1);
1895 serializer->doReal(&iparams->disres.up2);
1896 serializer->doReal(&iparams->disres.kfac);
1899 serializer->doInt(&iparams->orires.ex);
1900 serializer->doInt(&iparams->orires.label);
1901 serializer->doInt(&iparams->orires.power);
1902 serializer->doReal(&iparams->orires.c);
1903 serializer->doReal(&iparams->orires.obs);
1904 serializer->doReal(&iparams->orires.kfac);
1907 if (file_version < 82)
1909 serializer->doInt(&idum);
1910 serializer->doInt(&idum);
1912 serializer->doReal(&iparams->dihres.phiA);
1913 serializer->doReal(&iparams->dihres.dphiA);
1914 serializer->doReal(&iparams->dihres.kfacA);
1915 if (file_version >= 82)
1917 serializer->doReal(&iparams->dihres.phiB);
1918 serializer->doReal(&iparams->dihres.dphiB);
1919 serializer->doReal(&iparams->dihres.kfacB);
1923 iparams->dihres.phiB = iparams->dihres.phiA;
1924 iparams->dihres.dphiB = iparams->dihres.dphiA;
1925 iparams->dihres.kfacB = iparams->dihres.kfacA;
1929 serializer->doRvec(&iparams->posres.pos0A);
1930 serializer->doRvec(&iparams->posres.fcA);
1931 serializer->doRvec(&iparams->posres.pos0B);
1932 serializer->doRvec(&iparams->posres.fcB);
1935 serializer->doInt(&iparams->fbposres.geom);
1936 serializer->doRvec(&iparams->fbposres.pos0);
1937 serializer->doReal(&iparams->fbposres.r);
1938 serializer->doReal(&iparams->fbposres.k);
1940 case F_CBTDIHS: serializer->doRealArray(iparams->cbtdihs.cbtcA, NR_CBTDIHS); break;
1942 // Fall-through intended
1944 /* Fourier dihedrals are internally represented
1945 * as Ryckaert-Bellemans since those are faster to compute.
1947 serializer->doRealArray(iparams->rbdihs.rbcA, NR_RBDIHS);
1948 serializer->doRealArray(iparams->rbdihs.rbcB, NR_RBDIHS);
1952 serializer->doReal(&iparams->constr.dA);
1953 serializer->doReal(&iparams->constr.dB);
1956 serializer->doReal(&iparams->settle.doh);
1957 serializer->doReal(&iparams->settle.dhh);
1959 case F_VSITE1: break; // VSite1 has 0 parameters
1961 case F_VSITE2FD: serializer->doReal(&iparams->vsite.a); break;
1965 serializer->doReal(&iparams->vsite.a);
1966 serializer->doReal(&iparams->vsite.b);
1971 serializer->doReal(&iparams->vsite.a);
1972 serializer->doReal(&iparams->vsite.b);
1973 serializer->doReal(&iparams->vsite.c);
1976 serializer->doInt(&iparams->vsiten.n);
1977 serializer->doReal(&iparams->vsiten.a);
1979 case F_GB12_NOLONGERUSED:
1980 case F_GB13_NOLONGERUSED:
1981 case F_GB14_NOLONGERUSED:
1982 // Implicit solvent parameters can still be read, but never used
1983 if (serializer->reading())
1985 if (file_version < 68)
1987 serializer->doReal(&rdum);
1988 serializer->doReal(&rdum);
1989 serializer->doReal(&rdum);
1990 serializer->doReal(&rdum);
1992 if (file_version < tpxv_RemoveImplicitSolvation)
1994 serializer->doReal(&rdum);
1995 serializer->doReal(&rdum);
1996 serializer->doReal(&rdum);
1997 serializer->doReal(&rdum);
1998 serializer->doReal(&rdum);
2003 serializer->doInt(&iparams->cmap.cmapA);
2004 serializer->doInt(&iparams->cmap.cmapB);
2008 "unknown function type %d (%s) in %s line %d",
2010 interaction_function[ftype].name,
2016 static void do_ilist(gmx::ISerializer* serializer, InteractionList* ilist)
2018 int nr = ilist->size();
2019 serializer->doInt(&nr);
2020 if (serializer->reading())
2022 ilist->iatoms.resize(nr);
2024 serializer->doIntArray(ilist->iatoms.data(), ilist->size());
2027 static void do_ffparams(gmx::ISerializer* serializer, gmx_ffparams_t* ffparams, int file_version)
2029 serializer->doInt(&ffparams->atnr);
2030 int numTypes = ffparams->numTypes();
2031 serializer->doInt(&numTypes);
2032 if (serializer->reading())
2034 ffparams->functype.resize(numTypes);
2035 ffparams->iparams.resize(numTypes);
2037 /* Read/write all the function types */
2038 serializer->doIntArray(ffparams->functype.data(), ffparams->functype.size());
2040 if (file_version >= 66)
2042 serializer->doDouble(&ffparams->reppow);
2046 ffparams->reppow = 12.0;
2049 serializer->doReal(&ffparams->fudgeQQ);
2051 /* Check whether all these function types are supported by the code.
2052 * In practice the code is backwards compatible, which means that the
2053 * numbering may have to be altered from old numbering to new numbering
2055 for (int i = 0; i < ffparams->numTypes(); i++)
2057 if (serializer->reading())
2059 /* Loop over file versions */
2060 for (int k = 0; k < NFTUPD; k++)
2062 /* Compare the read file_version to the update table */
2063 if ((file_version < ftupd[k].fvnr) && (ffparams->functype[i] >= ftupd[k].ftype))
2065 ffparams->functype[i] += 1;
2070 do_iparams(serializer, ffparams->functype[i], &ffparams->iparams[i], file_version);
2074 static void add_settle_atoms(InteractionList* ilist)
2078 /* Settle used to only store the first atom: add the other two */
2079 ilist->iatoms.resize(2 * ilist->size());
2080 for (i = ilist->size() / 4 - 1; i >= 0; i--)
2082 ilist->iatoms[4 * i + 0] = ilist->iatoms[2 * i + 0];
2083 ilist->iatoms[4 * i + 1] = ilist->iatoms[2 * i + 1];
2084 ilist->iatoms[4 * i + 2] = ilist->iatoms[2 * i + 1] + 1;
2085 ilist->iatoms[4 * i + 3] = ilist->iatoms[2 * i + 1] + 2;
2089 static void do_ilists(gmx::ISerializer* serializer, InteractionLists* ilists, int file_version)
2091 GMX_RELEASE_ASSERT(ilists, "Need a valid ilists object");
2092 GMX_RELEASE_ASSERT(ilists->size() == F_NRE,
2093 "The code needs to be in sync with InteractionLists");
2095 for (int j = 0; j < F_NRE; j++)
2097 InteractionList& ilist = (*ilists)[j];
2098 gmx_bool bClear = FALSE;
2099 if (serializer->reading())
2101 for (int k = 0; k < NFTUPD; k++)
2103 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2111 ilist.iatoms.clear();
2115 do_ilist(serializer, &ilist);
2116 if (file_version < 78 && j == F_SETTLE && !ilist.empty())
2118 add_settle_atoms(&ilist);
2124 static void do_block(gmx::ISerializer* serializer, t_block* block)
2126 serializer->doInt(&block->nr);
2127 if (serializer->reading())
2129 if ((block->nalloc_index > 0) && (nullptr != block->index))
2131 sfree(block->index);
2133 block->nalloc_index = block->nr + 1;
2134 snew(block->index, block->nalloc_index);
2136 serializer->doIntArray(block->index, block->nr + 1);
2139 static void doListOfLists(gmx::ISerializer* serializer, gmx::ListOfLists<int>* listOfLists)
2141 int numLists = listOfLists->ssize();
2142 serializer->doInt(&numLists);
2143 int numElements = listOfLists->elementsView().ssize();
2144 serializer->doInt(&numElements);
2145 if (serializer->reading())
2147 std::vector<int> listRanges(numLists + 1);
2148 serializer->doIntArray(listRanges.data(), numLists + 1);
2149 std::vector<int> elements(numElements);
2150 serializer->doIntArray(elements.data(), numElements);
2151 *listOfLists = gmx::ListOfLists<int>(std::move(listRanges), std::move(elements));
2155 serializer->doIntArray(const_cast<int*>(listOfLists->listRangesView().data()), numLists + 1);
2156 serializer->doIntArray(const_cast<int*>(listOfLists->elementsView().data()), numElements);
2160 /* This is a primitive routine to make it possible to translate atomic numbers
2161 * to element names when reading TPR files, without making the Gromacs library
2162 * directory a dependency on mdrun (which is the case if we need elements.dat).
2164 static const char* atomicnumber_to_element(int atomicnumber)
2168 /* This does not have to be complete, so we only include elements likely
2169 * to occur in PDB files.
2171 switch (atomicnumber)
2173 case 1: p = "H"; break;
2174 case 5: p = "B"; break;
2175 case 6: p = "C"; break;
2176 case 7: p = "N"; break;
2177 case 8: p = "O"; break;
2178 case 9: p = "F"; break;
2179 case 11: p = "Na"; break;
2180 case 12: p = "Mg"; break;
2181 case 15: p = "P"; break;
2182 case 16: p = "S"; break;
2183 case 17: p = "Cl"; break;
2184 case 18: p = "Ar"; break;
2185 case 19: p = "K"; break;
2186 case 20: p = "Ca"; break;
2187 case 25: p = "Mn"; break;
2188 case 26: p = "Fe"; break;
2189 case 28: p = "Ni"; break;
2190 case 29: p = "Cu"; break;
2191 case 30: p = "Zn"; break;
2192 case 35: p = "Br"; break;
2193 case 47: p = "Ag"; break;
2194 default: p = ""; break;
2200 static void do_atom(gmx::ISerializer* serializer, t_atom* atom)
2202 serializer->doReal(&atom->m);
2203 serializer->doReal(&atom->q);
2204 serializer->doReal(&atom->mB);
2205 serializer->doReal(&atom->qB);
2206 serializer->doUShort(&atom->type);
2207 serializer->doUShort(&atom->typeB);
2208 serializer->doInt(&atom->ptype);
2209 serializer->doInt(&atom->resind);
2210 serializer->doInt(&atom->atomnumber);
2211 if (serializer->reading())
2213 /* Set element string from atomic number if present.
2214 * This routine returns an empty string if the name is not found.
2216 std::strncpy(atom->elem, atomicnumber_to_element(atom->atomnumber), 4);
2217 /* avoid warnings about potentially unterminated string */
2218 atom->elem[3] = '\0';
2222 static void do_grps(gmx::ISerializer* serializer, gmx::ArrayRef<AtomGroupIndices> grps)
2224 for (auto& group : grps)
2226 int size = group.size();
2227 serializer->doInt(&size);
2228 if (serializer->reading())
2232 serializer->doIntArray(group.data(), size);
2236 static void do_symstr(gmx::ISerializer* serializer, char*** nm, t_symtab* symtab)
2240 if (serializer->reading())
2242 serializer->doInt(&ls);
2243 *nm = get_symtab_handle(symtab, ls);
2247 ls = lookup_symtab(symtab, *nm);
2248 serializer->doInt(&ls);
2252 static void do_strstr(gmx::ISerializer* serializer, int nstr, char*** nm, t_symtab* symtab)
2256 for (j = 0; (j < nstr); j++)
2258 do_symstr(serializer, &(nm[j]), symtab);
2262 static void do_resinfo(gmx::ISerializer* serializer, int n, t_resinfo* ri, t_symtab* symtab, int file_version)
2266 for (j = 0; (j < n); j++)
2268 do_symstr(serializer, &(ri[j].name), symtab);
2269 if (file_version >= 63)
2271 serializer->doInt(&ri[j].nr);
2272 serializer->doUChar(&ri[j].ic);
2282 static void do_atoms(gmx::ISerializer* serializer, t_atoms* atoms, t_symtab* symtab, int file_version)
2286 serializer->doInt(&atoms->nr);
2287 serializer->doInt(&atoms->nres);
2288 if (serializer->reading())
2290 /* Since we have always written all t_atom properties in the tpr file
2291 * (at least for all backward compatible versions), we don't store
2292 * but simple set the booleans here.
2294 atoms->haveMass = TRUE;
2295 atoms->haveCharge = TRUE;
2296 atoms->haveType = TRUE;
2297 atoms->haveBState = TRUE;
2298 atoms->havePdbInfo = FALSE;
2300 snew(atoms->atom, atoms->nr);
2301 snew(atoms->atomname, atoms->nr);
2302 snew(atoms->atomtype, atoms->nr);
2303 snew(atoms->atomtypeB, atoms->nr);
2304 snew(atoms->resinfo, atoms->nres);
2305 atoms->pdbinfo = nullptr;
2309 GMX_RELEASE_ASSERT(atoms->haveMass && atoms->haveCharge && atoms->haveType && atoms->haveBState,
2310 "Mass, charge, atomtype and B-state parameters should be present in "
2311 "t_atoms when writing a tpr file");
2313 for (i = 0; (i < atoms->nr); i++)
2315 do_atom(serializer, &atoms->atom[i]);
2317 do_strstr(serializer, atoms->nr, atoms->atomname, symtab);
2318 do_strstr(serializer, atoms->nr, atoms->atomtype, symtab);
2319 do_strstr(serializer, atoms->nr, atoms->atomtypeB, symtab);
2321 do_resinfo(serializer, atoms->nres, atoms->resinfo, symtab, file_version);
2324 static void do_groups(gmx::ISerializer* serializer, SimulationGroups* groups, t_symtab* symtab)
2326 do_grps(serializer, groups->groups);
2327 int numberOfGroupNames = groups->groupNames.size();
2328 serializer->doInt(&numberOfGroupNames);
2329 if (serializer->reading())
2331 groups->groupNames.resize(numberOfGroupNames);
2333 do_strstr(serializer, numberOfGroupNames, groups->groupNames.data(), symtab);
2334 for (auto group : gmx::keysOf(groups->groupNumbers))
2336 int numberOfGroupNumbers = groups->numberOfGroupNumbers(group);
2337 serializer->doInt(&numberOfGroupNumbers);
2338 if (numberOfGroupNumbers != 0)
2340 if (serializer->reading())
2342 groups->groupNumbers[group].resize(numberOfGroupNumbers);
2344 serializer->doUCharArray(groups->groupNumbers[group].data(), numberOfGroupNumbers);
2349 static void do_atomtypes(gmx::ISerializer* serializer, t_atomtypes* atomtypes, int file_version)
2353 serializer->doInt(&atomtypes->nr);
2355 if (serializer->reading())
2357 snew(atomtypes->atomnumber, j);
2359 if (serializer->reading() && file_version < tpxv_RemoveImplicitSolvation)
2361 std::vector<real> dummy(atomtypes->nr, 0);
2362 serializer->doRealArray(dummy.data(), dummy.size());
2363 serializer->doRealArray(dummy.data(), dummy.size());
2364 serializer->doRealArray(dummy.data(), dummy.size());
2366 serializer->doIntArray(atomtypes->atomnumber, j);
2368 if (serializer->reading() && file_version >= 60 && file_version < tpxv_RemoveImplicitSolvation)
2370 std::vector<real> dummy(atomtypes->nr, 0);
2371 serializer->doRealArray(dummy.data(), dummy.size());
2372 serializer->doRealArray(dummy.data(), dummy.size());
2376 static void do_symtab(gmx::ISerializer* serializer, t_symtab* symtab)
2381 serializer->doInt(&symtab->nr);
2383 if (serializer->reading())
2385 snew(symtab->symbuf, 1);
2386 symbuf = symtab->symbuf;
2387 symbuf->bufsize = nr;
2388 snew(symbuf->buf, nr);
2389 for (i = 0; (i < nr); i++)
2392 serializer->doString(&buf);
2393 symbuf->buf[i] = gmx_strdup(buf.c_str());
2398 symbuf = symtab->symbuf;
2399 while (symbuf != nullptr)
2401 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2403 std::string buf = symbuf->buf[i];
2404 serializer->doString(&buf);
2407 symbuf = symbuf->next;
2411 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2416 static void do_cmap(gmx::ISerializer* serializer, gmx_cmap_t* cmap_grid)
2419 int ngrid = cmap_grid->cmapdata.size();
2420 serializer->doInt(&ngrid);
2421 serializer->doInt(&cmap_grid->grid_spacing);
2423 int gs = cmap_grid->grid_spacing;
2424 int nelem = gs * gs;
2426 if (serializer->reading())
2428 cmap_grid->cmapdata.resize(ngrid);
2430 for (int i = 0; i < ngrid; i++)
2432 cmap_grid->cmapdata[i].cmap.resize(4 * nelem);
2436 for (int i = 0; i < ngrid; i++)
2438 for (int j = 0; j < nelem; j++)
2440 serializer->doReal(&cmap_grid->cmapdata[i].cmap[j * 4]);
2441 serializer->doReal(&cmap_grid->cmapdata[i].cmap[j * 4 + 1]);
2442 serializer->doReal(&cmap_grid->cmapdata[i].cmap[j * 4 + 2]);
2443 serializer->doReal(&cmap_grid->cmapdata[i].cmap[j * 4 + 3]);
2449 static void do_moltype(gmx::ISerializer* serializer, gmx_moltype_t* molt, t_symtab* symtab, int file_version)
2451 do_symstr(serializer, &(molt->name), symtab);
2453 do_atoms(serializer, &molt->atoms, symtab, file_version);
2455 do_ilists(serializer, &molt->ilist, file_version);
2457 /* TODO: Remove the obsolete charge group index from the file */
2459 cgs.nr = molt->atoms.nr;
2460 cgs.nalloc_index = cgs.nr + 1;
2461 snew(cgs.index, cgs.nalloc_index);
2462 for (int i = 0; i < cgs.nr + 1; i++)
2466 do_block(serializer, &cgs);
2469 /* This used to be in the atoms struct */
2470 doListOfLists(serializer, &molt->excls);
2473 static void do_molblock(gmx::ISerializer* serializer, gmx_molblock_t* molb, int numAtomsPerMolecule)
2475 serializer->doInt(&molb->type);
2476 serializer->doInt(&molb->nmol);
2477 /* To maintain forward topology reading compatibility, we store #atoms.
2478 * TODO: Change this to conditional reading of a dummy int when we
2479 * increase tpx_generation.
2481 serializer->doInt(&numAtomsPerMolecule);
2482 /* Position restraint coordinates */
2483 int numPosres_xA = molb->posres_xA.size();
2484 serializer->doInt(&numPosres_xA);
2485 if (numPosres_xA > 0)
2487 if (serializer->reading())
2489 molb->posres_xA.resize(numPosres_xA);
2491 serializer->doRvecArray(as_rvec_array(molb->posres_xA.data()), numPosres_xA);
2493 int numPosres_xB = molb->posres_xB.size();
2494 serializer->doInt(&numPosres_xB);
2495 if (numPosres_xB > 0)
2497 if (serializer->reading())
2499 molb->posres_xB.resize(numPosres_xB);
2501 serializer->doRvecArray(as_rvec_array(molb->posres_xB.data()), numPosres_xB);
2505 static void set_disres_npair(gmx_mtop_t* mtop)
2507 gmx_mtop_ilistloop_t iloop;
2510 gmx::ArrayRef<t_iparams> ip = mtop->ffparams.iparams;
2512 iloop = gmx_mtop_ilistloop_init(mtop);
2513 while (const InteractionLists* ilist = gmx_mtop_ilistloop_next(iloop, &nmol))
2515 const InteractionList& il = (*ilist)[F_DISRES];
2519 gmx::ArrayRef<const int> a = il.iatoms;
2521 for (int i = 0; i < il.size(); i += 3)
2524 if (i + 3 == il.size() || ip[a[i]].disres.label != ip[a[i + 3]].disres.label)
2526 ip[a[i]].disres.npair = npair;
2534 static void do_mtop(gmx::ISerializer* serializer, gmx_mtop_t* mtop, int file_version)
2536 do_symtab(serializer, &(mtop->symtab));
2538 do_symstr(serializer, &(mtop->name), &(mtop->symtab));
2540 do_ffparams(serializer, &mtop->ffparams, file_version);
2542 int nmoltype = mtop->moltype.size();
2543 serializer->doInt(&nmoltype);
2544 if (serializer->reading())
2546 mtop->moltype.resize(nmoltype);
2548 for (gmx_moltype_t& moltype : mtop->moltype)
2550 do_moltype(serializer, &moltype, &mtop->symtab, file_version);
2553 int nmolblock = mtop->molblock.size();
2554 serializer->doInt(&nmolblock);
2555 if (serializer->reading())
2557 mtop->molblock.resize(nmolblock);
2559 for (gmx_molblock_t& molblock : mtop->molblock)
2561 int numAtomsPerMolecule = (serializer->reading() ? 0 : mtop->moltype[molblock.type].atoms.nr);
2562 do_molblock(serializer, &molblock, numAtomsPerMolecule);
2564 serializer->doInt(&mtop->natoms);
2566 if (file_version >= tpxv_IntermolecularBondeds)
2568 serializer->doBool(&mtop->bIntermolecularInteractions);
2569 if (mtop->bIntermolecularInteractions)
2571 if (serializer->reading())
2573 mtop->intermolecular_ilist = std::make_unique<InteractionLists>();
2575 do_ilists(serializer, mtop->intermolecular_ilist.get(), file_version);
2580 mtop->bIntermolecularInteractions = FALSE;
2583 do_atomtypes(serializer, &(mtop->atomtypes), file_version);
2585 if (file_version >= 65)
2587 do_cmap(serializer, &mtop->ffparams.cmap_grid);
2591 mtop->ffparams.cmap_grid.grid_spacing = 0;
2592 mtop->ffparams.cmap_grid.cmapdata.clear();
2595 do_groups(serializer, &mtop->groups, &(mtop->symtab));
2597 mtop->haveMoleculeIndices = true;
2599 if (file_version >= tpxv_StoreNonBondedInteractionExclusionGroup)
2601 std::int64_t intermolecularExclusionGroupSize = gmx::ssize(mtop->intermolecularExclusionGroup);
2602 serializer->doInt64(&intermolecularExclusionGroupSize);
2603 GMX_RELEASE_ASSERT(intermolecularExclusionGroupSize >= 0,
2604 "Number of atoms with excluded intermolecular non-bonded interactions "
2606 mtop->intermolecularExclusionGroup.resize(intermolecularExclusionGroupSize); // no effect when writing
2607 serializer->doIntArray(mtop->intermolecularExclusionGroup.data(),
2608 mtop->intermolecularExclusionGroup.size());
2611 if (serializer->reading())
2613 close_symtab(&(mtop->symtab));
2618 * Read the first part of the TPR file to find general system information.
2620 * If \p TopOnlyOK is true then we can read even future versions
2621 * of tpx files, provided the \p fileGeneration hasn't changed.
2622 * If it is false, we need the \p ir too, and bail out
2623 * if the file is newer than the program.
2625 * The version and generation of the topology (see top of this file)
2626 * are returned in the two last arguments, if those arguments are non-nullptr.
2628 * If possible, we will read the \p ir even when \p TopOnlyOK is true.
2630 * \param[in,out] serializer The serializer used to handle header processing.
2631 * \param[in,out] tpx File header datastructure.
2632 * \param[in] filename The name of the file being read/written
2633 * \param[in,out] fio File handle.
2634 * \param[in] TopOnlyOK If not reading \p ir is fine or not.
2636 static void do_tpxheader(gmx::FileIOXdrSerializer* serializer,
2638 const char* filename,
2646 /* XDR binary topology file */
2647 precision = sizeof(real);
2649 std::string fileTag;
2650 if (serializer->reading())
2652 serializer->doString(&buf);
2653 if (std::strncmp(buf.c_str(), "VERSION", 7) != 0)
2657 "Can not read file %s,\n"
2658 " this file is from a GROMACS version which is older than 2.0\n"
2659 " Make a new one with grompp or use a gro or pdb file, if possible",
2662 // We need to know the precision used to write the TPR file, to match it
2663 // to the precision of the currently running binary. If the precisions match
2664 // there is no problem, but mismatching precision needs to be accounted for
2665 // by reading into temporary variables of the correct precision instead
2666 // of the desired target datastructures.
2667 serializer->doInt(&precision);
2668 tpx->isDouble = (precision == sizeof(double));
2669 if ((precision != sizeof(float)) && !tpx->isDouble)
2672 "Unknown precision in file %s: real is %d bytes "
2673 "instead of %zu or %zu",
2679 gmx_fio_setprecision(fio, tpx->isDouble);
2681 "Reading file %s, %s (%s precision)\n",
2684 tpx->isDouble ? "double" : "single");
2688 buf = gmx::formatString("VERSION %s", gmx_version());
2689 serializer->doString(&buf);
2690 gmx_fio_setprecision(fio, tpx->isDouble);
2691 serializer->doInt(&precision);
2692 fileTag = gmx::formatString("%s", tpx_tag);
2695 /* Check versions! */
2696 serializer->doInt(&tpx->fileVersion);
2698 /* This is for backward compatibility with development versions 77-79
2699 * where the tag was, mistakenly, placed before the generation,
2700 * which would cause a segv instead of a proper error message
2701 * when reading the topology only from tpx with <77 code.
2703 if (tpx->fileVersion >= 77 && tpx->fileVersion <= 79)
2705 serializer->doString(&fileTag);
2708 serializer->doInt(&tpx->fileGeneration);
2710 if (tpx->fileVersion >= 81)
2712 serializer->doString(&fileTag);
2714 if (serializer->reading())
2716 if (tpx->fileVersion < 77)
2718 /* Versions before 77 don't have the tag, set it to release */
2719 fileTag = gmx::formatString("%s", TPX_TAG_RELEASE);
2722 if (fileTag != tpx_tag)
2724 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n", fileTag.c_str(), tpx_tag);
2726 /* We only support reading tpx files with the same tag as the code
2727 * or tpx files with the release tag and with lower version number.
2729 if (fileTag != TPX_TAG_RELEASE && tpx->fileVersion < tpx_version)
2732 "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' "
2733 "with program for tpx version %d, tag '%s'",
2743 if ((tpx->fileVersion <= tpx_incompatible_version)
2744 || ((tpx->fileVersion > tpx_version) && !TopOnlyOK) || (tpx->fileGeneration > tpx_generation)
2745 || tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
2748 "reading tpx file (%s) version %d with version %d program",
2754 serializer->doInt(&tpx->natoms);
2755 serializer->doInt(&tpx->ngtc);
2757 if (tpx->fileVersion < 62)
2759 serializer->doInt(&idum);
2760 serializer->doReal(&rdum);
2762 if (tpx->fileVersion >= 79)
2764 serializer->doInt(&tpx->fep_state);
2766 serializer->doReal(&tpx->lambda);
2767 serializer->doBool(&tpx->bIr);
2768 serializer->doBool(&tpx->bTop);
2769 serializer->doBool(&tpx->bX);
2770 serializer->doBool(&tpx->bV);
2771 serializer->doBool(&tpx->bF);
2772 serializer->doBool(&tpx->bBox);
2774 if (tpx->fileVersion >= tpxv_AddSizeField && tpx->fileGeneration >= 27)
2776 if (!serializer->reading())
2778 GMX_RELEASE_ASSERT(tpx->sizeOfTprBody != 0,
2779 "Not possible to write new file with zero TPR body size");
2781 serializer->doInt64(&tpx->sizeOfTprBody);
2784 if ((tpx->fileGeneration > tpx_generation))
2786 /* This can only happen if TopOnlyOK=TRUE */
2791 #define do_test(serializer, b, p) \
2792 if ((serializer)->reading() && ((p) != nullptr) && !(b)) \
2793 gmx_fatal(FARGS, "No %s in input file", #p)
2796 * Process the first part of the TPR into the state datastructure.
2798 * Due to the structure of the legacy code, it is necessary
2799 * to split up the state reading into two parts, with the
2800 * box and legacy temperature coupling processed before the
2801 * topology datastructures.
2803 * See the documentation for do_tpx_body for the correct order of
2804 * the operations for reading a tpr file.
2806 * \param[in] serializer Abstract serializer used to read/write data.
2807 * \param[in] tpx The file header data.
2808 * \param[in, out] state Global state data.
2810 static void do_tpx_state_first(gmx::ISerializer* serializer, TpxFileHeader* tpx, t_state* state)
2812 if (serializer->reading())
2815 init_gtc_state(state, tpx->ngtc, 0, 0);
2817 do_test(serializer, tpx->bBox, state->box);
2820 serializer->doRvecArray(state->box, DIM);
2821 if (tpx->fileVersion >= 51)
2823 serializer->doRvecArray(state->box_rel, DIM);
2827 /* We initialize box_rel after reading the inputrec */
2828 clear_mat(state->box_rel);
2830 serializer->doRvecArray(state->boxv, DIM);
2831 if (tpx->fileVersion < 56)
2834 serializer->doRvecArray(mdum, DIM);
2838 if (state->ngtc > 0)
2841 snew(dumv, state->ngtc);
2842 if (tpx->fileVersion < 69)
2844 serializer->doRealArray(dumv, state->ngtc);
2846 /* These used to be the Berendsen tcoupl_lambda's */
2847 serializer->doRealArray(dumv, state->ngtc);
2853 * Process global topology data.
2855 * See the documentation for do_tpx_body for the correct order of
2856 * the operations for reading a tpr file.
2858 * \param[in] serializer Abstract serializer used to read/write data.
2859 * \param[in] tpx The file header data.
2860 * \param[in,out] mtop Global topology.
2862 static void do_tpx_mtop(gmx::ISerializer* serializer, TpxFileHeader* tpx, gmx_mtop_t* mtop)
2864 do_test(serializer, tpx->bTop, mtop);
2869 do_mtop(serializer, mtop, tpx->fileVersion);
2870 set_disres_npair(mtop);
2876 do_mtop(serializer, &dum_top, tpx->fileVersion);
2881 * Process coordinate vectors for state data.
2883 * Main part of state gets processed here.
2885 * See the documentation for do_tpx_body for the correct order of
2886 * the operations for reading a tpr file.
2888 * \param[in] serializer Abstract serializer used to read/write data.
2889 * \param[in] tpx The file header data.
2890 * \param[in,out] state Global state data.
2891 * \param[in,out] x Individual coordinates for processing, deprecated.
2892 * \param[in,out] v Individual velocities for processing, deprecated.
2894 static void do_tpx_state_second(gmx::ISerializer* serializer, TpxFileHeader* tpx, t_state* state, rvec* x, rvec* v)
2896 if (!serializer->reading())
2899 x == nullptr && v == nullptr,
2900 "Passing separate x and v pointers to do_tpx() is not supported when writing");
2904 GMX_RELEASE_ASSERT(!(x == nullptr && v != nullptr),
2905 "Passing x==NULL and v!=NULL is not supported");
2908 if (serializer->reading())
2912 // v is also nullptr by the above assertion, so we may
2913 // need to make memory in state for storing the contents
2917 state->flags |= (1 << estX);
2921 state->flags |= (1 << estV);
2923 state_change_natoms(state, tpx->natoms);
2929 x = state->x.rvec_array();
2930 v = state->v.rvec_array();
2932 do_test(serializer, tpx->bX, x);
2935 if (serializer->reading())
2937 state->flags |= (1 << estX);
2939 serializer->doRvecArray(x, tpx->natoms);
2942 do_test(serializer, tpx->bV, v);
2945 if (serializer->reading())
2947 state->flags |= (1 << estV);
2951 std::vector<gmx::RVec> dummyVelocities(tpx->natoms);
2952 serializer->doRvecArray(as_rvec_array(dummyVelocities.data()), tpx->natoms);
2956 serializer->doRvecArray(v, tpx->natoms);
2960 // No need to run do_test when the last argument is NULL
2963 std::vector<gmx::RVec> dummyForces(state->natoms);
2964 serializer->doRvecArray(as_rvec_array(dummyForces.data()), tpx->natoms);
2968 * Process simulation parameters.
2970 * See the documentation for do_tpx_body for the correct order of
2971 * the operations for reading a tpr file.
2973 * \param[in] serializer Abstract serializer used to read/write data.
2974 * \param[in] tpx The file header data.
2975 * \param[in,out] ir Datastructure with simulation parameters.
2977 static PbcType do_tpx_ir(gmx::ISerializer* serializer, TpxFileHeader* tpx, t_inputrec* ir)
2980 gmx_bool bPeriodicMols;
2982 /* Starting with tpx version 26, we have the inputrec
2983 * at the end of the file, so we can ignore it
2984 * if the file is never than the software (but still the
2985 * same generation - see comments at the top of this file.
2989 pbcType = PbcType::Unset;
2990 bPeriodicMols = FALSE;
2992 do_test(serializer, tpx->bIr, ir);
2995 if (tpx->fileVersion >= 53)
2997 /* Removed the pbc info from do_inputrec, since we always want it */
2998 if (!serializer->reading())
3000 pbcType = ir->pbcType;
3001 bPeriodicMols = ir->bPeriodicMols;
3003 serializer->doInt(reinterpret_cast<int*>(&pbcType));
3004 serializer->doBool(&bPeriodicMols);
3006 if (tpx->fileGeneration <= tpx_generation && ir)
3008 do_inputrec(serializer, ir, tpx->fileVersion);
3009 if (tpx->fileVersion < 53)
3011 pbcType = ir->pbcType;
3012 bPeriodicMols = ir->bPeriodicMols;
3015 if (serializer->reading() && ir && tpx->fileVersion >= 53)
3017 /* We need to do this after do_inputrec, since that initializes ir */
3018 ir->pbcType = pbcType;
3019 ir->bPeriodicMols = bPeriodicMols;
3026 * Correct and finalize read information.
3028 * If \p state is nullptr, skip the parts dependent on it.
3030 * See the documentation for do_tpx_body for the correct order of
3031 * the operations for reading a tpr file.
3033 * \param[in] tpx The file header used to check version numbers.
3034 * \param[out] ir Input rec that needs correction.
3035 * \param[out] state State needing correction.
3036 * \param[out] mtop Topology to finalize.
3038 static void do_tpx_finalize(TpxFileHeader* tpx, t_inputrec* ir, t_state* state, gmx_mtop_t* mtop)
3040 if (tpx->fileVersion < 51 && state)
3042 set_box_rel(ir, state);
3046 if (state && state->ngtc == 0)
3048 /* Reading old version without tcoupl state data: set it */
3049 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3051 if (tpx->bTop && mtop)
3053 if (tpx->fileVersion < 57)
3055 if (!mtop->moltype[0].ilist[F_DISRES].empty())
3057 ir->eDisre = edrSimple;
3061 ir->eDisre = edrNone;
3069 * Process TPR data for file reading/writing.
3071 * The TPR file gets processed in in four stages due to the organization
3072 * of the data within it.
3074 * First, state data for the box is processed in do_tpx_state_first.
3075 * This is followed by processing the topology in do_tpx_mtop.
3076 * Coordinate and velocity vectors are handled next in do_tpx_state_second.
3077 * The last file information processed is the collection of simulation parameters in do_tpx_ir.
3078 * When reading, a final processing step is undertaken at the end.
3080 * \param[in] serializer Abstract serializer used to read/write data.
3081 * \param[in] tpx The file header data.
3082 * \param[in,out] ir Datastructures with simulation parameters.
3083 * \param[in,out] state Global state data.
3084 * \param[in,out] x Individual coordinates for processing, deprecated.
3085 * \param[in,out] v Individual velocities for processing, deprecated.
3086 * \param[in,out] mtop Global topology.
3088 static PbcType do_tpx_body(gmx::ISerializer* serializer,
3098 do_tpx_state_first(serializer, tpx, state);
3100 do_tpx_mtop(serializer, tpx, mtop);
3103 do_tpx_state_second(serializer, tpx, state, x, v);
3105 PbcType pbcType = do_tpx_ir(serializer, tpx, ir);
3106 if (serializer->reading())
3108 do_tpx_finalize(tpx, ir, state, mtop);
3114 * Overload for do_tpx_body that defaults to state vectors being nullptr.
3116 * \param[in] serializer Abstract serializer used to read/write data.
3117 * \param[in] tpx The file header data.
3118 * \param[in,out] ir Datastructures with simulation parameters.
3119 * \param[in,out] mtop Global topology.
3121 static PbcType do_tpx_body(gmx::ISerializer* serializer, TpxFileHeader* tpx, t_inputrec* ir, gmx_mtop_t* mtop)
3123 return do_tpx_body(serializer, tpx, ir, nullptr, nullptr, nullptr, mtop);
3126 static t_fileio* open_tpx(const char* fn, const char* mode)
3128 return gmx_fio_open(fn, mode);
3131 static void close_tpx(t_fileio* fio)
3137 * Fill information into the header only from state before writing.
3139 * Populating the header needs to be independent from writing the information
3140 * to file to allow things like writing the raw byte stream.
3142 * \param[in] state The current simulation state. Can't write without it.
3143 * \param[in] ir Parameter and system information.
3144 * \param[in] mtop Global topology.
3145 * \returns Fully populated header.
3147 static TpxFileHeader populateTpxHeader(const t_state& state, const t_inputrec* ir, const gmx_mtop_t* mtop)
3149 TpxFileHeader header;
3150 header.natoms = state.natoms;
3151 header.ngtc = state.ngtc;
3152 header.fep_state = state.fep_state;
3153 header.lambda = state.lambda[efptFEP];
3154 header.bIr = ir != nullptr;
3155 header.bTop = mtop != nullptr;
3156 header.bX = (state.flags & (1 << estX)) != 0;
3157 header.bV = (state.flags & (1 << estV)) != 0;
3160 header.fileVersion = tpx_version;
3161 header.fileGeneration = tpx_generation;
3162 header.isDouble = (sizeof(real) == sizeof(double));
3168 * Process the body of a TPR file as an opaque data buffer.
3170 * Reads/writes the information in \p buffer from/to the \p serializer
3171 * provided to the function. Does not interact with the actual
3172 * TPR datastructures but with an in memory representation of the
3173 * data, so that this data can be efficiently read or written from/to
3174 * an original source.
3176 * \param[in] serializer The abstract serializer used for reading or writing
3177 * the information in \p buffer.
3178 * \param[in,out] buffer Information from TPR file as char buffer.
3180 static void doTpxBodyBuffer(gmx::ISerializer* serializer, gmx::ArrayRef<char> buffer)
3182 serializer->doOpaque(buffer.data(), buffer.size());
3186 * Populates simulation datastructures.
3188 * Here the information from the serialization interface \p serializer
3189 * is used to first populate the datastructures containing the simulation
3190 * information. Depending on the version found in the header \p tpx,
3191 * this is done using the new reading of the data as one block from disk,
3192 * followed by complete deserialization of the information read from there.
3193 * Otherwise, the datastructures are populated as before one by one from disk.
3194 * The second version is the default for the legacy tools that read the
3195 * coordinates and velocities separate from the state.
3197 * After reading in the data, a separate buffer is populated from them
3198 * containing only \p ir and \p mtop that can be communicated directly
3199 * to nodes needing the information to set up a simulation.
3201 * \param[in] tpx The file header.
3202 * \param[in] serializer The Serialization interface used to read the TPR.
3203 * \param[out] ir Input rec to populate.
3204 * \param[out] state State vectors to populate.
3205 * \param[out] x Coordinates to populate if needed.
3206 * \param[out] v Velocities to populate if needed.
3207 * \param[out] mtop Global topology to populate.
3209 * \returns Partial de-serialized TPR used for communication to nodes.
3211 static PartialDeserializedTprFile readTpxBody(TpxFileHeader* tpx,
3212 gmx::ISerializer* serializer,
3219 PartialDeserializedTprFile partialDeserializedTpr;
3220 if (tpx->fileVersion >= tpxv_AddSizeField && tpx->fileGeneration >= 27)
3222 partialDeserializedTpr.body.resize(tpx->sizeOfTprBody);
3223 partialDeserializedTpr.header = *tpx;
3224 doTpxBodyBuffer(serializer, partialDeserializedTpr.body);
3226 partialDeserializedTpr.pbcType =
3227 completeTprDeserialization(&partialDeserializedTpr, ir, state, x, v, mtop);
3231 partialDeserializedTpr.pbcType = do_tpx_body(serializer, tpx, ir, state, x, v, mtop);
3233 // Update header to system info for communication to nodes.
3234 // As we only need to communicate the inputrec and mtop to other nodes,
3235 // we prepare a new char buffer with the information we have already read
3237 partialDeserializedTpr.header = populateTpxHeader(*state, ir, mtop);
3238 // Long-term we should move to use little endian in files to avoid extra byte swapping,
3239 // but since we just used the default XDR format (which is big endian) for the TPR
3240 // header it would cause third-party libraries reading our raw data to tear their hair
3241 // if we swap the endian in the middle of the file, so we stick to big endian in the
3242 // TPR file for now - and thus we ask the serializer to swap if this host is little endian.
3243 gmx::InMemorySerializer tprBodySerializer(gmx::EndianSwapBehavior::SwapIfHostIsLittleEndian);
3244 do_tpx_body(&tprBodySerializer, &partialDeserializedTpr.header, ir, mtop);
3245 partialDeserializedTpr.body = tprBodySerializer.finishAndGetBuffer();
3247 return partialDeserializedTpr;
3250 /************************************************************
3252 * The following routines are the exported ones
3254 ************************************************************/
3256 TpxFileHeader readTpxHeader(const char* fileName, bool canReadTopologyOnly)
3260 fio = open_tpx(fileName, "r");
3261 gmx::FileIOXdrSerializer serializer(fio);
3264 do_tpxheader(&serializer, &tpx, fileName, fio, canReadTopologyOnly);
3269 void write_tpx_state(const char* fn, const t_inputrec* ir, const t_state* state, const gmx_mtop_t* mtop)
3271 /* To write a state, we first need to write the state information to a buffer before
3272 * we append the raw bytes to the file. For this, the header information needs to be
3273 * populated before we write the main body because it has some information that is
3274 * otherwise not available.
3279 TpxFileHeader tpx = populateTpxHeader(*state, ir, mtop);
3280 // Long-term we should move to use little endian in files to avoid extra byte swapping,
3281 // but since we just used the default XDR format (which is big endian) for the TPR
3282 // header it would cause third-party libraries reading our raw data to tear their hair
3283 // if we swap the endian in the middle of the file, so we stick to big endian in the
3284 // TPR file for now - and thus we ask the serializer to swap if this host is little endian.
3285 gmx::InMemorySerializer tprBodySerializer(gmx::EndianSwapBehavior::SwapIfHostIsLittleEndian);
3287 do_tpx_body(&tprBodySerializer,
3289 const_cast<t_inputrec*>(ir),
3290 const_cast<t_state*>(state),
3293 const_cast<gmx_mtop_t*>(mtop));
3295 std::vector<char> tprBody = tprBodySerializer.finishAndGetBuffer();
3296 tpx.sizeOfTprBody = tprBody.size();
3298 fio = open_tpx(fn, "w");
3299 gmx::FileIOXdrSerializer serializer(fio);
3300 do_tpxheader(&serializer, &tpx, fn, fio, ir == nullptr);
3301 doTpxBodyBuffer(&serializer, tprBody);
3306 PbcType completeTprDeserialization(PartialDeserializedTprFile* partialDeserializedTpr,
3313 // Long-term we should move to use little endian in files to avoid extra byte swapping,
3314 // but since we just used the default XDR format (which is big endian) for the TPR
3315 // header it would cause third-party libraries reading our raw data to tear their hair
3316 // if we swap the endian in the middle of the file, so we stick to big endian in the
3317 // TPR file for now - and thus we ask the serializer to swap if this host is little endian.
3318 gmx::InMemoryDeserializer tprBodyDeserializer(partialDeserializedTpr->body,
3319 partialDeserializedTpr->header.isDouble,
3320 gmx::EndianSwapBehavior::SwapIfHostIsLittleEndian);
3321 return do_tpx_body(&tprBodyDeserializer, &partialDeserializedTpr->header, ir, state, x, v, mtop);
3324 PbcType completeTprDeserialization(PartialDeserializedTprFile* partialDeserializedTpr,
3328 return completeTprDeserialization(partialDeserializedTpr, ir, nullptr, nullptr, nullptr, mtop);
3331 PartialDeserializedTprFile read_tpx_state(const char* fn, t_inputrec* ir, t_state* state, gmx_mtop_t* mtop)
3334 fio = open_tpx(fn, "r");
3335 gmx::FileIOXdrSerializer serializer(fio);
3336 PartialDeserializedTprFile partialDeserializedTpr;
3337 do_tpxheader(&serializer, &partialDeserializedTpr.header, fn, fio, ir == nullptr);
3338 partialDeserializedTpr =
3339 readTpxBody(&partialDeserializedTpr.header, &serializer, ir, state, nullptr, nullptr, mtop);
3341 return partialDeserializedTpr;
3344 PbcType read_tpx(const char* fn, t_inputrec* ir, matrix box, int* natoms, rvec* x, rvec* v, gmx_mtop_t* mtop)
3350 fio = open_tpx(fn, "r");
3351 gmx::FileIOXdrSerializer serializer(fio);
3352 do_tpxheader(&serializer, &tpx, fn, fio, ir == nullptr);
3353 PartialDeserializedTprFile partialDeserializedTpr =
3354 readTpxBody(&tpx, &serializer, ir, &state, x, v, mtop);
3356 if (mtop != nullptr && natoms != nullptr)
3358 *natoms = mtop->natoms;
3362 copy_mat(state.box, box);
3364 return partialDeserializedTpr.pbcType;
3367 PbcType read_tpx_top(const char* fn, t_inputrec* ir, matrix box, int* natoms, rvec* x, rvec* v, t_topology* top)
3372 pbcType = read_tpx(fn, ir, box, natoms, x, v, &mtop);
3374 *top = gmx_mtop_t_to_t_topology(&mtop, true);
3379 gmx_bool fn2bTPX(const char* file)
3381 return (efTPR == fn2ftp(file));
3384 void pr_tpxheader(FILE* fp, int indent, const char* title, const TpxFileHeader* sh)
3386 if (available(fp, sh, indent, title))
3388 indent = pr_title(fp, indent, title);
3389 pr_indent(fp, indent);
3390 fprintf(fp, "bIr = %spresent\n", sh->bIr ? "" : "not ");
3391 pr_indent(fp, indent);
3392 fprintf(fp, "bBox = %spresent\n", sh->bBox ? "" : "not ");
3393 pr_indent(fp, indent);
3394 fprintf(fp, "bTop = %spresent\n", sh->bTop ? "" : "not ");
3395 pr_indent(fp, indent);
3396 fprintf(fp, "bX = %spresent\n", sh->bX ? "" : "not ");
3397 pr_indent(fp, indent);
3398 fprintf(fp, "bV = %spresent\n", sh->bV ? "" : "not ");
3399 pr_indent(fp, indent);
3400 fprintf(fp, "bF = %spresent\n", sh->bF ? "" : "not ");
3402 pr_indent(fp, indent);
3403 fprintf(fp, "natoms = %d\n", sh->natoms);
3404 pr_indent(fp, indent);
3405 fprintf(fp, "lambda = %e\n", sh->lambda);
3406 pr_indent(fp, indent);
3407 fprintf(fp, "buffer size = %" PRId64 "\n", sh->sizeOfTprBody);