2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 /* This file is completely threadsafe - keep it that way! */
51 #include "gromacs/fileio/filetypes.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/gmxfio_xdr.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdtypes/awh_history.h"
57 #include "gromacs/mdtypes/awh_params.h"
58 #include "gromacs/mdtypes/inputrec.h"
59 #include "gromacs/mdtypes/md_enums.h"
60 #include "gromacs/mdtypes/pull_params.h"
61 #include "gromacs/mdtypes/state.h"
62 #include "gromacs/pbcutil/boxutilities.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/topology/block.h"
65 #include "gromacs/topology/ifunc.h"
66 #include "gromacs/topology/mtop_util.h"
67 #include "gromacs/topology/symtab.h"
68 #include "gromacs/topology/topology.h"
69 #include "gromacs/utility/arraysize.h"
70 #include "gromacs/utility/baseversion.h"
71 #include "gromacs/utility/cstringutil.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/futil.h"
74 #include "gromacs/utility/gmxassert.h"
75 #include "gromacs/utility/keyvaluetreebuilder.h"
76 #include "gromacs/utility/keyvaluetreeserializer.h"
77 #include "gromacs/utility/smalloc.h"
78 #include "gromacs/utility/snprintf.h"
79 #include "gromacs/utility/txtdump.h"
81 #define TPX_TAG_RELEASE "release"
83 /*! \brief Tag string for the file format written to run input files
84 * written by this version of the code.
86 * Change this if you want to change the run input format in a feature
87 * branch. This ensures that there will not be different run input
88 * formats around which cannot be distinguished, while not causing
89 * problems rebasing the feature branch onto upstream changes. When
90 * merging with mainstream GROMACS, set this tag string back to
91 * TPX_TAG_RELEASE, and instead add an element to tpxv.
93 static const char *tpx_tag = TPX_TAG_RELEASE;
95 /*! \brief Enum of values that describe the contents of a tpr file
96 * whose format matches a version number
98 * The enum helps the code be more self-documenting and ensure merges
99 * do not silently resolve when two patches make the same bump. When
100 * adding new functionality, add a new element just above tpxv_Count
101 * in this enumeration, and write code below that does the right thing
102 * according to the value of file_version.
106 tpxv_ComputationalElectrophysiology = 96, /**< support for ion/water position swaps (computational electrophysiology) */
107 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to int64_t */
108 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
109 tpxv_InteractiveMolecularDynamics, /**< interactive molecular dynamics (IMD) */
110 tpxv_RemoveObsoleteParameters1, /**< remove optimize_fft, dihre_fc, nstcheckpoint */
111 tpxv_PullCoordTypeGeom, /**< add pull type and geometry per group and flat-bottom */
112 tpxv_PullGeomDirRel, /**< add pull geometry direction-relative */
113 tpxv_IntermolecularBondeds, /**< permit inter-molecular bonded interactions in the topology */
114 tpxv_CompElWithSwapLayerOffset, /**< added parameters for improved CompEl setups */
115 tpxv_CompElPolyatomicIonsAndMultipleIonTypes, /**< CompEl now can handle polyatomic ions and more than two types of ions */
116 tpxv_RemoveAdress, /**< removed support for AdResS */
117 tpxv_PullCoordNGroup, /**< add ngroup to pull coord */
118 tpxv_RemoveTwinRange, /**< removed support for twin-range interactions */
119 tpxv_ReplacePullPrintCOM12, /**< Replaced print-com-1, 2 with pull-print-com */
120 tpxv_PullExternalPotential, /**< Added pull type external potential */
121 tpxv_GenericParamsForElectricField, /**< Introduced KeyValueTree and moved electric field parameters */
122 tpxv_AcceleratedWeightHistogram, /**< sampling with accelerated weight histogram method (AWH) */
123 tpxv_RemoveImplicitSolvation, /**< removed support for implicit solvation */
124 tpxv_PullPrevStepCOMAsReference, /**< Enabled using the COM of the pull group of the last frame as reference for PBC */
125 tpxv_MimicQMMM, /**< Introduced support for MiMiC QM/MM interface */
126 tpxv_PullAverage, /**< Added possibility to output average pull force and position */
127 tpxv_GenericInternalParameters, /**< Added internal parameters for mdrun modules*/
128 tpxv_Count /**< the total number of tpxv versions */
131 /*! \brief Version number of the file format written to run input
132 * files by this version of the code.
134 * The tpx_version increases whenever the file format in the main
135 * development branch changes, due to an extension of the tpxv enum above.
136 * Backward compatibility for reading old run input files is maintained
137 * by checking this version number against that of the file and then using
138 * the correct code path.
140 * When developing a feature branch that needs to change the run input
141 * file format, change tpx_tag instead. */
142 static const int tpx_version = tpxv_Count - 1;
145 /* This number should only be increased when you edit the TOPOLOGY section
146 * or the HEADER of the tpx format.
147 * This way we can maintain forward compatibility too for all analysis tools
148 * and/or external programs that only need to know the atom/residue names,
149 * charges, and bond connectivity.
151 * It first appeared in tpx version 26, when I also moved the inputrecord
152 * to the end of the tpx file, so we can just skip it if we only
155 * In particular, it must be increased when adding new elements to
156 * ftupd, so that old code can read new .tpr files.
158 static const int tpx_generation = 26;
160 /* This number should be the most recent backwards incompatible version
161 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
163 static const int tpx_incompatible_version = 57; // GMX4.0 has version 58
167 /* Struct used to maintain tpx compatibility when function types are added */
169 int fvnr; /* file version number in which the function type first appeared */
170 int ftype; /* function type */
174 * TODO The following three lines make little sense, please clarify if
175 * you've had to work out how ftupd works.
177 * The entries should be ordered in:
178 * 1. ascending function type number
179 * 2. ascending file version number
181 * Because we support reading of old .tpr file versions (even when
182 * mdrun can no longer run the simulation), we need to be able to read
183 * obsolete t_interaction_function types. Any data read from such
184 * fields is discarded. Their names have _NOLONGERUSED appended to
185 * them to make things clear.
187 static const t_ftupd ftupd[] = {
188 { 70, F_RESTRBONDS },
189 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
190 { 76, F_LINEAR_ANGLES },
191 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
192 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
194 { 60, F_GB12_NOLONGERUSED },
195 { 61, F_GB13_NOLONGERUSED },
196 { 61, F_GB14_NOLONGERUSED },
197 { 72, F_GBPOL_NOLONGERUSED },
198 { 72, F_NPSOLVATION_NOLONGERUSED },
201 { 69, F_VTEMP_NOLONGERUSED},
203 { 76, F_ANHARM_POL },
206 { 79, F_DVDL_BONDED, },
207 { 79, F_DVDL_RESTRAINT },
208 { 79, F_DVDL_TEMPERATURE },
209 { tpxv_GenericInternalParameters, F_DENSITYFITTING },
211 #define NFTUPD asize(ftupd)
213 /* Needed for backward compatibility */
216 /**************************************************************
218 * Now the higer level routines that do io of the structures and arrays
220 **************************************************************/
221 static void do_pullgrp_tpx_pre95(gmx::ISerializer *serializer,
227 serializer->doInt(&pgrp->nat);
228 if (serializer->reading())
230 snew(pgrp->ind, pgrp->nat);
232 serializer->doIntArray(pgrp->ind, pgrp->nat);
233 serializer->doInt(&pgrp->nweight);
234 if (serializer->reading())
236 snew(pgrp->weight, pgrp->nweight);
238 serializer->doRealArray(pgrp->weight, pgrp->nweight);
239 serializer->doInt(&pgrp->pbcatom);
240 serializer->doRvec(&pcrd->vec);
241 clear_rvec(pcrd->origin);
242 serializer->doRvec(&tmp);
244 serializer->doReal(&pcrd->rate);
245 serializer->doReal(&pcrd->k);
246 serializer->doReal(&pcrd->kB);
249 static void do_pull_group(gmx::ISerializer *serializer, t_pull_group *pgrp)
251 serializer->doInt(&pgrp->nat);
252 if (serializer->reading())
254 snew(pgrp->ind, pgrp->nat);
256 serializer->doIntArray(pgrp->ind, pgrp->nat);
257 serializer->doInt(&pgrp->nweight);
258 if (serializer->reading())
260 snew(pgrp->weight, pgrp->nweight);
262 serializer->doRealArray(pgrp->weight, pgrp->nweight);
263 serializer->doInt(&pgrp->pbcatom);
266 static void do_pull_coord(gmx::ISerializer *serializer,
269 int ePullOld, int eGeomOld, ivec dimOld)
271 if (file_version >= tpxv_PullCoordNGroup)
273 serializer->doInt(&pcrd->eType);
274 if (file_version >= tpxv_PullExternalPotential)
276 if (pcrd->eType == epullEXTERNAL)
279 if (serializer->reading())
281 serializer->doString(&buf);
282 pcrd->externalPotentialProvider = gmx_strdup(buf.c_str());
286 buf = pcrd->externalPotentialProvider;
287 serializer->doString(&buf);
292 pcrd->externalPotentialProvider = nullptr;
297 if (serializer->reading())
299 pcrd->externalPotentialProvider = nullptr;
302 /* Note that we try to support adding new geometries without
303 * changing the tpx version. This requires checks when printing the
304 * geometry string and a check and fatal_error in init_pull.
306 serializer->doInt(&pcrd->eGeom);
307 serializer->doInt(&pcrd->ngroup);
308 if (pcrd->ngroup <= c_pullCoordNgroupMax)
310 serializer->doIntArray(pcrd->group, pcrd->ngroup);
314 /* More groups in file than supported, this must be a new geometry
315 * that is not supported by our current code. Since we will not
316 * use the groups for this coord (checks in the pull and WHAM code
317 * ensure this), we can ignore the groups and set ngroup=0.
320 snew(dum, pcrd->ngroup);
321 serializer->doIntArray(dum, pcrd->ngroup);
326 serializer->doIvec(&pcrd->dim);
331 serializer->doInt(&pcrd->group[0]);
332 serializer->doInt(&pcrd->group[1]);
333 if (file_version >= tpxv_PullCoordTypeGeom)
335 pcrd->ngroup = (pcrd->eGeom == epullgDIRRELATIVE ? 4 : 2);
336 serializer->doInt(&pcrd->eType);
337 serializer->doInt(&pcrd->eGeom);
338 if (pcrd->ngroup == 4)
340 serializer->doInt(&pcrd->group[2]);
341 serializer->doInt(&pcrd->group[3]);
343 serializer->doIvec(&pcrd->dim);
347 pcrd->eType = ePullOld;
348 pcrd->eGeom = eGeomOld;
349 copy_ivec(dimOld, pcrd->dim);
352 serializer->doRvec(&pcrd->origin);
353 serializer->doRvec(&pcrd->vec);
354 if (file_version >= tpxv_PullCoordTypeGeom)
356 serializer->doBool(&pcrd->bStart);
360 /* This parameter is only printed, but not actually used by mdrun */
361 pcrd->bStart = FALSE;
363 serializer->doReal(&pcrd->init);
364 serializer->doReal(&pcrd->rate);
365 serializer->doReal(&pcrd->k);
366 serializer->doReal(&pcrd->kB);
369 static void do_expandedvals(gmx::ISerializer *serializer,
374 int n_lambda = fepvals->n_lambda;
376 /* reset the lambda calculation window */
377 fepvals->lambda_start_n = 0;
378 fepvals->lambda_stop_n = n_lambda;
379 if (file_version >= 79)
383 if (serializer->reading())
385 snew(expand->init_lambda_weights, n_lambda);
387 serializer->doRealArray(expand->init_lambda_weights, n_lambda);
388 serializer->doBool(&expand->bInit_weights);
391 serializer->doInt(&expand->nstexpanded);
392 serializer->doInt(&expand->elmcmove);
393 serializer->doInt(&expand->elamstats);
394 serializer->doInt(&expand->lmc_repeats);
395 serializer->doInt(&expand->gibbsdeltalam);
396 serializer->doInt(&expand->lmc_forced_nstart);
397 serializer->doInt(&expand->lmc_seed);
398 serializer->doReal(&expand->mc_temp);
399 serializer->doBool(&expand->bSymmetrizedTMatrix);
400 serializer->doInt(&expand->nstTij);
401 serializer->doInt(&expand->minvarmin);
402 serializer->doInt(&expand->c_range);
403 serializer->doReal(&expand->wl_scale);
404 serializer->doReal(&expand->wl_ratio);
405 serializer->doReal(&expand->init_wl_delta);
406 serializer->doBool(&expand->bWLoneovert);
407 serializer->doInt(&expand->elmceq);
408 serializer->doInt(&expand->equil_steps);
409 serializer->doInt(&expand->equil_samples);
410 serializer->doInt(&expand->equil_n_at_lam);
411 serializer->doReal(&expand->equil_wl_delta);
412 serializer->doReal(&expand->equil_ratio);
416 static void do_simtempvals(gmx::ISerializer *serializer,
421 if (file_version >= 79)
423 serializer->doInt(&simtemp->eSimTempScale);
424 serializer->doReal(&simtemp->simtemp_high);
425 serializer->doReal(&simtemp->simtemp_low);
428 if (serializer->reading())
430 snew(simtemp->temperatures, n_lambda);
432 serializer->doRealArray(simtemp->temperatures, n_lambda);
437 static void do_imd(gmx::ISerializer *serializer, t_IMD *imd)
439 serializer->doInt(&imd->nat);
440 if (serializer->reading())
442 snew(imd->ind, imd->nat);
444 serializer->doIntArray(imd->ind, imd->nat);
447 static void do_fepvals(gmx::ISerializer *serializer,
451 /* i is defined in the ndo_double macro; use g to iterate. */
455 /* free energy values */
457 if (file_version >= 79)
459 serializer->doInt(&fepvals->init_fep_state);
460 serializer->doDouble(&fepvals->init_lambda);
461 serializer->doDouble(&fepvals->delta_lambda);
463 else if (file_version >= 59)
465 serializer->doDouble(&fepvals->init_lambda);
466 serializer->doDouble(&fepvals->delta_lambda);
470 serializer->doReal(&rdum);
471 fepvals->init_lambda = rdum;
472 serializer->doReal(&rdum);
473 fepvals->delta_lambda = rdum;
475 if (file_version >= 79)
477 serializer->doInt(&fepvals->n_lambda);
478 if (serializer->reading())
480 snew(fepvals->all_lambda, efptNR);
482 for (g = 0; g < efptNR; g++)
484 if (fepvals->n_lambda > 0)
486 if (serializer->reading())
488 snew(fepvals->all_lambda[g], fepvals->n_lambda);
490 serializer->doDoubleArray(fepvals->all_lambda[g], fepvals->n_lambda);
491 serializer->doBoolArray(fepvals->separate_dvdl, efptNR);
493 else if (fepvals->init_lambda >= 0)
495 fepvals->separate_dvdl[efptFEP] = TRUE;
499 else if (file_version >= 64)
501 serializer->doInt(&fepvals->n_lambda);
502 if (serializer->reading())
506 snew(fepvals->all_lambda, efptNR);
507 /* still allocate the all_lambda array's contents. */
508 for (g = 0; g < efptNR; g++)
510 if (fepvals->n_lambda > 0)
512 snew(fepvals->all_lambda[g], fepvals->n_lambda);
516 serializer->doDoubleArray(fepvals->all_lambda[efptFEP], fepvals->n_lambda);
517 if (fepvals->init_lambda >= 0)
521 fepvals->separate_dvdl[efptFEP] = TRUE;
523 if (serializer->reading())
525 /* copy the contents of the efptFEP lambda component to all
526 the other components */
527 for (g = 0; g < efptNR; g++)
529 for (h = 0; h < fepvals->n_lambda; h++)
533 fepvals->all_lambda[g][h] =
534 fepvals->all_lambda[efptFEP][h];
543 fepvals->n_lambda = 0;
544 fepvals->all_lambda = nullptr;
545 if (fepvals->init_lambda >= 0)
547 fepvals->separate_dvdl[efptFEP] = TRUE;
550 serializer->doReal(&fepvals->sc_alpha);
551 serializer->doInt(&fepvals->sc_power);
552 if (file_version >= 79)
554 serializer->doReal(&fepvals->sc_r_power);
558 fepvals->sc_r_power = 6.0;
560 serializer->doReal(&fepvals->sc_sigma);
561 if (serializer->reading())
563 if (file_version >= 71)
565 fepvals->sc_sigma_min = fepvals->sc_sigma;
569 fepvals->sc_sigma_min = 0;
572 if (file_version >= 79)
574 serializer->doBool(&fepvals->bScCoul);
578 fepvals->bScCoul = TRUE;
580 if (file_version >= 64)
582 serializer->doInt(&fepvals->nstdhdl);
586 fepvals->nstdhdl = 1;
589 if (file_version >= 73)
591 serializer->doInt(&fepvals->separate_dhdl_file);
592 serializer->doInt(&fepvals->dhdl_derivatives);
596 fepvals->separate_dhdl_file = esepdhdlfileYES;
597 fepvals->dhdl_derivatives = edhdlderivativesYES;
599 if (file_version >= 71)
601 serializer->doInt(&fepvals->dh_hist_size);
602 serializer->doDouble(&fepvals->dh_hist_spacing);
606 fepvals->dh_hist_size = 0;
607 fepvals->dh_hist_spacing = 0.1;
609 if (file_version >= 79)
611 serializer->doInt(&fepvals->edHdLPrintEnergy);
615 fepvals->edHdLPrintEnergy = edHdLPrintEnergyNO;
618 /* handle lambda_neighbors */
619 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
621 serializer->doInt(&fepvals->lambda_neighbors);
622 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
623 (fepvals->init_lambda < 0) )
625 fepvals->lambda_start_n = (fepvals->init_fep_state -
626 fepvals->lambda_neighbors);
627 fepvals->lambda_stop_n = (fepvals->init_fep_state +
628 fepvals->lambda_neighbors + 1);
629 if (fepvals->lambda_start_n < 0)
631 fepvals->lambda_start_n = 0;;
633 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
635 fepvals->lambda_stop_n = fepvals->n_lambda;
640 fepvals->lambda_start_n = 0;
641 fepvals->lambda_stop_n = fepvals->n_lambda;
646 fepvals->lambda_start_n = 0;
647 fepvals->lambda_stop_n = fepvals->n_lambda;
651 static void do_awhBias(gmx::ISerializer *serializer, gmx::AwhBiasParams *awhBiasParams,
652 int gmx_unused file_version)
654 serializer->doInt(&awhBiasParams->eTarget);
655 serializer->doDouble(&awhBiasParams->targetBetaScaling);
656 serializer->doDouble(&awhBiasParams->targetCutoff);
657 serializer->doInt(&awhBiasParams->eGrowth);
658 serializer->doInt(&awhBiasParams->bUserData);
659 serializer->doDouble(&awhBiasParams->errorInitial);
660 serializer->doInt(&awhBiasParams->ndim);
661 serializer->doInt(&awhBiasParams->shareGroup);
662 serializer->doBool(&awhBiasParams->equilibrateHistogram);
664 if (serializer->reading())
666 snew(awhBiasParams->dimParams, awhBiasParams->ndim);
669 for (int d = 0; d < awhBiasParams->ndim; d++)
671 gmx::AwhDimParams *dimParams = &awhBiasParams->dimParams[d];
673 serializer->doInt(&dimParams->eCoordProvider);
674 serializer->doInt(&dimParams->coordIndex);
675 serializer->doDouble(&dimParams->origin);
676 serializer->doDouble(&dimParams->end);
677 serializer->doDouble(&dimParams->period);
678 serializer->doDouble(&dimParams->forceConstant);
679 serializer->doDouble(&dimParams->diffusion);
680 serializer->doDouble(&dimParams->coordValueInit);
681 serializer->doDouble(&dimParams->coverDiameter);
685 static void do_awh(gmx::ISerializer *serializer,
686 gmx::AwhParams *awhParams,
687 int gmx_unused file_version)
689 serializer->doInt(&awhParams->numBias);
690 serializer->doInt(&awhParams->nstOut);
691 serializer->doInt64(&awhParams->seed);
692 serializer->doInt(&awhParams->nstSampleCoord);
693 serializer->doInt(&awhParams->numSamplesUpdateFreeEnergy);
694 serializer->doInt(&awhParams->ePotential);
695 serializer->doBool(&awhParams->shareBiasMultisim);
697 if (awhParams->numBias > 0)
699 if (serializer->reading())
701 snew(awhParams->awhBiasParams, awhParams->numBias);
704 for (int k = 0; k < awhParams->numBias; k++)
706 do_awhBias(serializer, &awhParams->awhBiasParams[k], file_version);
711 static void do_pull(gmx::ISerializer *serializer,
720 if (file_version >= 95)
722 serializer->doInt(&pull->ngroup);
724 serializer->doInt(&pull->ncoord);
725 if (file_version < 95)
727 pull->ngroup = pull->ncoord + 1;
729 if (file_version < tpxv_PullCoordTypeGeom)
733 serializer->doInt(&eGeomOld);
734 serializer->doIvec(&dimOld);
735 /* The inner cylinder radius, now removed */
736 serializer->doReal(&dum);
738 serializer->doReal(&pull->cylinder_r);
739 serializer->doReal(&pull->constr_tol);
740 if (file_version >= 95)
742 serializer->doBool(&pull->bPrintCOM);
743 /* With file_version < 95 this value is set below */
745 if (file_version >= tpxv_ReplacePullPrintCOM12)
747 serializer->doBool(&pull->bPrintRefValue);
748 serializer->doBool(&pull->bPrintComp);
750 else if (file_version >= tpxv_PullCoordTypeGeom)
753 serializer->doInt(&idum); /* used to be bPrintCOM2 */
754 serializer->doBool(&pull->bPrintRefValue);
755 serializer->doBool(&pull->bPrintComp);
759 pull->bPrintRefValue = FALSE;
760 pull->bPrintComp = TRUE;
762 serializer->doInt(&pull->nstxout);
763 serializer->doInt(&pull->nstfout);
764 if (file_version >= tpxv_PullPrevStepCOMAsReference)
766 serializer->doBool(&pull->bSetPbcRefToPrevStepCOM);
770 pull->bSetPbcRefToPrevStepCOM = FALSE;
772 if (serializer->reading())
774 snew(pull->group, pull->ngroup);
775 snew(pull->coord, pull->ncoord);
777 if (file_version < 95)
779 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
780 if (eGeomOld == epullgDIRPBC)
782 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
784 if (eGeomOld > epullgDIRPBC)
789 for (g = 0; g < pull->ngroup; g++)
791 /* We read and ignore a pull coordinate for group 0 */
792 do_pullgrp_tpx_pre95(serializer, &pull->group[g], &pull->coord[std::max(g-1, 0)]);
795 pull->coord[g-1].group[0] = 0;
796 pull->coord[g-1].group[1] = g;
800 pull->bPrintCOM = (pull->group[0].nat > 0);
804 for (g = 0; g < pull->ngroup; g++)
806 do_pull_group(serializer, &pull->group[g]);
808 for (g = 0; g < pull->ncoord; g++)
810 do_pull_coord(serializer, &pull->coord[g],
811 file_version, ePullOld, eGeomOld, dimOld);
814 if (file_version >= tpxv_PullAverage)
818 v = pull->bXOutAverage;
819 serializer->doBool(&v);
820 pull->bXOutAverage = v;
821 v = pull->bFOutAverage;
822 serializer->doBool(&v);
823 pull->bFOutAverage = v;
828 static void do_rotgrp(gmx::ISerializer *serializer,
831 serializer->doInt(&rotg->eType);
832 serializer->doInt(&rotg->bMassW);
833 serializer->doInt(&rotg->nat);
834 if (serializer->reading())
836 snew(rotg->ind, rotg->nat);
838 serializer->doIntArray(rotg->ind, rotg->nat);
839 if (serializer->reading())
841 snew(rotg->x_ref, rotg->nat);
843 serializer->doRvecArray(rotg->x_ref, rotg->nat);
844 serializer->doRvec(&rotg->inputVec);
845 serializer->doRvec(&rotg->pivot);
846 serializer->doReal(&rotg->rate);
847 serializer->doReal(&rotg->k);
848 serializer->doReal(&rotg->slab_dist);
849 serializer->doReal(&rotg->min_gaussian);
850 serializer->doReal(&rotg->eps);
851 serializer->doInt(&rotg->eFittype);
852 serializer->doInt(&rotg->PotAngle_nstep);
853 serializer->doReal(&rotg->PotAngle_step);
856 static void do_rot(gmx::ISerializer *serializer,
861 serializer->doInt(&rot->ngrp);
862 serializer->doInt(&rot->nstrout);
863 serializer->doInt(&rot->nstsout);
864 if (serializer->reading())
866 snew(rot->grp, rot->ngrp);
868 for (g = 0; g < rot->ngrp; g++)
870 do_rotgrp(serializer, &rot->grp[g]);
875 static void do_swapgroup(gmx::ISerializer *serializer,
879 /* Name of the group or molecule */
881 if (serializer->reading())
883 serializer->doString(&buf);
884 g->molname = gmx_strdup(buf.c_str());
889 serializer->doString(&buf);
892 /* Number of atoms in the group */
893 serializer->doInt(&g->nat);
895 /* The group's atom indices */
896 if (serializer->reading())
898 snew(g->ind, g->nat);
900 serializer->doIntArray(g->ind, g->nat);
902 /* Requested counts for compartments A and B */
903 serializer->doIntArray(g->nmolReq, eCompNR);
906 static void do_swapcoords_tpx(gmx::ISerializer *serializer,
910 /* Enums for better readability of the code */
915 eChannel0 = 0, eChannel1
919 if (file_version >= tpxv_CompElPolyatomicIonsAndMultipleIonTypes)
921 /* The total number of swap groups is the sum of the fixed groups
922 * (split0, split1, solvent), and the user-defined groups (2+ types of ions)
924 serializer->doInt(&swap->ngrp);
925 if (serializer->reading())
927 snew(swap->grp, swap->ngrp);
929 for (int ig = 0; ig < swap->ngrp; ig++)
931 do_swapgroup(serializer, &swap->grp[ig]);
933 serializer->doBool(&swap->massw_split[eChannel0]);
934 serializer->doBool(&swap->massw_split[eChannel1]);
935 serializer->doInt(&swap->nstswap);
936 serializer->doInt(&swap->nAverage);
937 serializer->doReal(&swap->threshold);
938 serializer->doReal(&swap->cyl0r);
939 serializer->doReal(&swap->cyl0u);
940 serializer->doReal(&swap->cyl0l);
941 serializer->doReal(&swap->cyl1r);
942 serializer->doReal(&swap->cyl1u);
943 serializer->doReal(&swap->cyl1l);
947 /*** Support reading older CompEl .tpr files ***/
949 /* In the original CompEl .tpr files, we always have 5 groups: */
951 snew(swap->grp, swap->ngrp);
953 swap->grp[eGrpSplit0 ].molname = gmx_strdup("split0" ); // group 0: split0
954 swap->grp[eGrpSplit1 ].molname = gmx_strdup("split1" ); // group 1: split1
955 swap->grp[eGrpSolvent].molname = gmx_strdup("solvent"); // group 2: solvent
956 swap->grp[3 ].molname = gmx_strdup("anions" ); // group 3: anions
957 swap->grp[4 ].molname = gmx_strdup("cations"); // group 4: cations
959 serializer->doInt(&swap->grp[3].nat);
960 serializer->doInt(&swap->grp[eGrpSolvent].nat);
961 serializer->doInt(&swap->grp[eGrpSplit0].nat);
962 serializer->doBool(&swap->massw_split[eChannel0]);
963 serializer->doInt(&swap->grp[eGrpSplit1].nat);
964 serializer->doBool(&swap->massw_split[eChannel1]);
965 serializer->doInt(&swap->nstswap);
966 serializer->doInt(&swap->nAverage);
967 serializer->doReal(&swap->threshold);
968 serializer->doReal(&swap->cyl0r);
969 serializer->doReal(&swap->cyl0u);
970 serializer->doReal(&swap->cyl0l);
971 serializer->doReal(&swap->cyl1r);
972 serializer->doReal(&swap->cyl1u);
973 serializer->doReal(&swap->cyl1l);
975 // The order[] array keeps compatibility with older .tpr files
976 // by reading in the groups in the classic order
978 const int order[4] = {3, eGrpSolvent, eGrpSplit0, eGrpSplit1 };
980 for (int ig = 0; ig < 4; ig++)
983 snew(swap->grp[g].ind, swap->grp[g].nat);
984 serializer->doIntArray(swap->grp[g].ind, swap->grp[g].nat);
988 for (int j = eCompA; j <= eCompB; j++)
990 serializer->doInt(&swap->grp[3].nmolReq[j]); // group 3 = anions
991 serializer->doInt(&swap->grp[4].nmolReq[j]); // group 4 = cations
993 } /* End support reading older CompEl .tpr files */
995 if (file_version >= tpxv_CompElWithSwapLayerOffset)
997 serializer->doReal(&swap->bulkOffset[eCompA]);
998 serializer->doReal(&swap->bulkOffset[eCompB]);
1003 static void do_legacy_efield(gmx::ISerializer *serializer, gmx::KeyValueTreeObjectBuilder *root)
1005 const char *const dimName[] = { "x", "y", "z" };
1007 auto appliedForcesObj = root->addObject("applied-forces");
1008 auto efieldObj = appliedForcesObj.addObject("electric-field");
1009 // The content of the tpr file for this feature has
1010 // been the same since gromacs 4.0 that was used for
1012 for (int j = 0; j < DIM; ++j)
1015 serializer->doInt(&n);
1016 serializer->doInt(&nt);
1017 std::vector<real> aa(n+1), phi(nt+1), at(nt+1), phit(nt+1);
1018 serializer->doRealArray(aa.data(), n);
1019 serializer->doRealArray(phi.data(), n);
1020 serializer->doRealArray(at.data(), nt);
1021 serializer->doRealArray(phit.data(), nt);
1024 if (n > 1 || nt > 1)
1026 gmx_fatal(FARGS, "Can not handle tpr files with more than one electric field term per direction.");
1028 auto dimObj = efieldObj.addObject(dimName[j]);
1029 dimObj.addValue<real>("E0", aa[0]);
1030 dimObj.addValue<real>("omega", at[0]);
1031 dimObj.addValue<real>("t0", phi[0]);
1032 dimObj.addValue<real>("sigma", phit[0]);
1038 static void do_inputrec(gmx::ISerializer *serializer,
1042 int i, j, k, idum = 0;
1044 gmx_bool bdum = false;
1046 if (file_version != tpx_version)
1048 /* Give a warning about features that are not accessible */
1049 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
1050 file_version, tpx_version);
1053 if (file_version == 0)
1058 gmx::KeyValueTreeBuilder paramsBuilder;
1059 gmx::KeyValueTreeObjectBuilder paramsObj = paramsBuilder.rootObject();
1061 /* Basic inputrec stuff */
1062 serializer->doInt(&ir->eI);
1063 if (file_version >= 62)
1065 serializer->doInt64(&ir->nsteps);
1069 serializer->doInt(&idum);
1073 if (file_version >= 62)
1075 serializer->doInt64(&ir->init_step);
1079 serializer->doInt(&idum);
1080 ir->init_step = idum;
1083 serializer->doInt(&ir->simulation_part);
1085 if (file_version >= 67)
1087 serializer->doInt(&ir->nstcalcenergy);
1091 ir->nstcalcenergy = 1;
1093 if (file_version >= 81)
1095 serializer->doInt(&ir->cutoff_scheme);
1096 if (file_version < 94)
1098 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
1103 ir->cutoff_scheme = ecutsGROUP;
1105 serializer->doInt(&ir->ns_type);
1106 serializer->doInt(&ir->nstlist);
1107 serializer->doInt(&idum); /* used to be ndelta; not used anymore */
1109 serializer->doReal(&ir->rtpi);
1111 serializer->doInt(&ir->nstcomm);
1112 serializer->doInt(&ir->comm_mode);
1114 /* ignore nstcheckpoint */
1115 if (file_version < tpxv_RemoveObsoleteParameters1)
1117 serializer->doInt(&idum);
1120 serializer->doInt(&ir->nstcgsteep);
1122 serializer->doInt(&ir->nbfgscorr);
1124 serializer->doInt(&ir->nstlog);
1125 serializer->doInt(&ir->nstxout);
1126 serializer->doInt(&ir->nstvout);
1127 serializer->doInt(&ir->nstfout);
1128 serializer->doInt(&ir->nstenergy);
1129 serializer->doInt(&ir->nstxout_compressed);
1130 if (file_version >= 59)
1132 serializer->doDouble(&ir->init_t);
1133 serializer->doDouble(&ir->delta_t);
1137 serializer->doReal(&rdum);
1139 serializer->doReal(&rdum);
1142 serializer->doReal(&ir->x_compression_precision);
1143 if (file_version >= 81)
1145 serializer->doReal(&ir->verletbuf_tol);
1149 ir->verletbuf_tol = 0;
1151 serializer->doReal(&ir->rlist);
1152 if (file_version >= 67 && file_version < tpxv_RemoveTwinRange)
1154 if (serializer->reading())
1156 // Reading such a file version could invoke the twin-range
1157 // scheme, about which mdrun should give a fatal error.
1158 real dummy_rlistlong = -1;
1159 serializer->doReal(&dummy_rlistlong);
1161 if (ir->rlist > 0 && (dummy_rlistlong == 0 || dummy_rlistlong > ir->rlist))
1163 // Get mdrun to issue an error (regardless of
1164 // ir->cutoff_scheme).
1165 ir->useTwinRange = true;
1169 // grompp used to set rlistlong actively. Users were
1170 // probably also confused and set rlistlong == rlist.
1171 // However, in all remaining cases, it is safe to let
1172 // mdrun proceed normally.
1173 ir->useTwinRange = false;
1179 // No need to read or write anything
1180 ir->useTwinRange = false;
1182 if (file_version >= 82 && file_version != 90)
1184 // Multiple time-stepping is no longer enabled, but the old
1185 // support required the twin-range scheme, for which mdrun
1186 // already emits a fatal error.
1187 int dummy_nstcalclr = -1;
1188 serializer->doInt(&dummy_nstcalclr);
1190 serializer->doInt(&ir->coulombtype);
1191 if (file_version >= 81)
1193 serializer->doInt(&ir->coulomb_modifier);
1197 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1199 serializer->doReal(&ir->rcoulomb_switch);
1200 serializer->doReal(&ir->rcoulomb);
1201 serializer->doInt(&ir->vdwtype);
1202 if (file_version >= 81)
1204 serializer->doInt(&ir->vdw_modifier);
1208 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1210 serializer->doReal(&ir->rvdw_switch);
1211 serializer->doReal(&ir->rvdw);
1212 serializer->doInt(&ir->eDispCorr);
1213 serializer->doReal(&ir->epsilon_r);
1214 serializer->doReal(&ir->epsilon_rf);
1215 serializer->doReal(&ir->tabext);
1217 // This permits reading a .tpr file that used implicit solvent,
1218 // and later permitting mdrun to refuse to run it.
1219 if (serializer->reading())
1221 if (file_version < tpxv_RemoveImplicitSolvation)
1223 serializer->doInt(&idum);
1224 serializer->doInt(&idum);
1225 serializer->doReal(&rdum);
1226 serializer->doReal(&rdum);
1227 serializer->doInt(&idum);
1228 ir->implicit_solvent = (idum > 0);
1232 ir->implicit_solvent = false;
1234 if (file_version < tpxv_RemoveImplicitSolvation)
1236 serializer->doReal(&rdum);
1237 serializer->doReal(&rdum);
1238 serializer->doReal(&rdum);
1239 serializer->doReal(&rdum);
1240 if (file_version >= 60)
1242 serializer->doReal(&rdum);
1243 serializer->doInt(&idum);
1245 serializer->doReal(&rdum);
1249 if (file_version >= 81)
1251 serializer->doReal(&ir->fourier_spacing);
1255 ir->fourier_spacing = 0.0;
1257 serializer->doInt(&ir->nkx);
1258 serializer->doInt(&ir->nky);
1259 serializer->doInt(&ir->nkz);
1260 serializer->doInt(&ir->pme_order);
1261 serializer->doReal(&ir->ewald_rtol);
1263 if (file_version >= 93)
1265 serializer->doReal(&ir->ewald_rtol_lj);
1269 ir->ewald_rtol_lj = ir->ewald_rtol;
1271 serializer->doInt(&ir->ewald_geometry);
1272 serializer->doReal(&ir->epsilon_surface);
1274 /* ignore bOptFFT */
1275 if (file_version < tpxv_RemoveObsoleteParameters1)
1277 serializer->doBool(&bdum);
1280 if (file_version >= 93)
1282 serializer->doInt(&ir->ljpme_combination_rule);
1284 serializer->doBool(&ir->bContinuation);
1285 serializer->doInt(&ir->etc);
1286 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1287 * but the values 0 and 1 still mean no and
1288 * berendsen temperature coupling, respectively.
1290 if (file_version >= 79)
1292 serializer->doBool(&ir->bPrintNHChains);
1294 if (file_version >= 71)
1296 serializer->doInt(&ir->nsttcouple);
1300 ir->nsttcouple = ir->nstcalcenergy;
1302 serializer->doInt(&ir->epc);
1303 serializer->doInt(&ir->epct);
1304 if (file_version >= 71)
1306 serializer->doInt(&ir->nstpcouple);
1310 ir->nstpcouple = ir->nstcalcenergy;
1312 serializer->doReal(&ir->tau_p);
1313 serializer->doRvec(&ir->ref_p[XX]);
1314 serializer->doRvec(&ir->ref_p[YY]);
1315 serializer->doRvec(&ir->ref_p[ZZ]);
1316 serializer->doRvec(&ir->compress[XX]);
1317 serializer->doRvec(&ir->compress[YY]);
1318 serializer->doRvec(&ir->compress[ZZ]);
1319 serializer->doInt(&ir->refcoord_scaling);
1320 serializer->doRvec(&ir->posres_com);
1321 serializer->doRvec(&ir->posres_comB);
1323 if (file_version < 79)
1325 serializer->doInt(&ir->andersen_seed);
1329 ir->andersen_seed = 0;
1332 serializer->doReal(&ir->shake_tol);
1334 serializer->doInt(&ir->efep);
1335 do_fepvals(serializer, ir->fepvals, file_version);
1337 if (file_version >= 79)
1339 serializer->doBool(&ir->bSimTemp);
1342 ir->bSimTemp = TRUE;
1347 ir->bSimTemp = FALSE;
1351 do_simtempvals(serializer, ir->simtempvals, ir->fepvals->n_lambda, file_version);
1354 if (file_version >= 79)
1356 serializer->doBool(&ir->bExpanded);
1359 ir->bExpanded = TRUE;
1363 ir->bExpanded = FALSE;
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())
1482 do_pull(serializer, ir->pull, 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,
1583 ir->opts.ngener*ir->opts.ngener);
1585 /* First read the lists with annealing and npoints for each group */
1586 serializer->doIntArray(ir->opts.annealing, ir->opts.ngtc);
1587 serializer->doIntArray(ir->opts.anneal_npoints, ir->opts.ngtc);
1588 for (j = 0; j < (ir->opts.ngtc); j++)
1590 k = ir->opts.anneal_npoints[j];
1591 if (serializer->reading())
1593 snew(ir->opts.anneal_time[j], k);
1594 snew(ir->opts.anneal_temp[j], k);
1596 serializer->doRealArray(ir->opts.anneal_time[j], k);
1597 serializer->doRealArray(ir->opts.anneal_temp[j], k);
1601 serializer->doInt(&ir->nwall);
1602 serializer->doInt(&ir->wall_type);
1603 serializer->doReal(&ir->wall_r_linpot);
1604 serializer->doInt(&ir->wall_atomtype[0]);
1605 serializer->doInt(&ir->wall_atomtype[1]);
1606 serializer->doReal(&ir->wall_density[0]);
1607 serializer->doReal(&ir->wall_density[1]);
1608 serializer->doReal(&ir->wall_ewald_zfac);
1611 /* Cosine stuff for electric fields */
1612 if (file_version < tpxv_GenericParamsForElectricField)
1614 do_legacy_efield(serializer, ¶msObj);
1618 if (file_version >= tpxv_ComputationalElectrophysiology)
1620 serializer->doInt(&ir->eSwapCoords);
1621 if (ir->eSwapCoords != eswapNO)
1623 if (serializer->reading())
1627 do_swapcoords_tpx(serializer, ir->swap, file_version);
1633 serializer->doBool(&ir->bQMMM);
1634 serializer->doInt(&ir->QMMMscheme);
1635 serializer->doReal(&ir->scalefactor);
1636 serializer->doInt(&ir->opts.ngQM);
1637 if (serializer->reading())
1639 snew(ir->opts.QMmethod, ir->opts.ngQM);
1640 snew(ir->opts.QMbasis, ir->opts.ngQM);
1641 snew(ir->opts.QMcharge, ir->opts.ngQM);
1642 snew(ir->opts.QMmult, ir->opts.ngQM);
1643 snew(ir->opts.bSH, ir->opts.ngQM);
1644 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1645 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1646 snew(ir->opts.SAon, ir->opts.ngQM);
1647 snew(ir->opts.SAoff, ir->opts.ngQM);
1648 snew(ir->opts.SAsteps, ir->opts.ngQM);
1650 if (ir->opts.ngQM > 0 && ir->bQMMM)
1652 serializer->doIntArray(ir->opts.QMmethod, ir->opts.ngQM);
1653 serializer->doIntArray(ir->opts.QMbasis, ir->opts.ngQM);
1654 serializer->doIntArray(ir->opts.QMcharge, ir->opts.ngQM);
1655 serializer->doIntArray(ir->opts.QMmult, ir->opts.ngQM);
1656 serializer->doBoolArray(ir->opts.bSH, ir->opts.ngQM);
1657 serializer->doIntArray(ir->opts.CASorbitals, ir->opts.ngQM);
1658 serializer->doIntArray(ir->opts.CASelectrons, ir->opts.ngQM);
1659 serializer->doRealArray(ir->opts.SAon, ir->opts.ngQM);
1660 serializer->doRealArray(ir->opts.SAoff, ir->opts.ngQM);
1661 serializer->doIntArray(ir->opts.SAsteps, ir->opts.ngQM);
1662 /* We leave in dummy i/o for removed parameters to avoid
1663 * changing the tpr format for every QMMM change.
1665 std::vector<int> dummy(ir->opts.ngQM, 0);
1666 serializer->doIntArray(dummy.data(), ir->opts.ngQM);
1667 serializer->doIntArray(dummy.data(), ir->opts.ngQM);
1669 /* end of QMMM stuff */
1672 if (file_version >= tpxv_GenericParamsForElectricField)
1674 if (serializer->reading())
1676 paramsObj.mergeObject(
1677 gmx::deserializeKeyValueTree(serializer));
1681 GMX_RELEASE_ASSERT(ir->params != nullptr,
1682 "Parameters should be present when writing inputrec");
1683 gmx::serializeKeyValueTree(*ir->params, serializer);
1686 if (serializer->reading())
1688 ir->params = new gmx::KeyValueTreeObject(paramsBuilder.build());
1689 // Initialize internal parameters to an empty kvt for all tpr versions
1690 ir->internalParameters = std::make_unique<gmx::KeyValueTreeObject>();
1693 if (file_version >= tpxv_GenericInternalParameters)
1695 if (serializer->reading())
1697 ir->internalParameters = std::make_unique<gmx::KeyValueTreeObject>(gmx::deserializeKeyValueTree(serializer));
1701 GMX_RELEASE_ASSERT(ir->internalParameters != nullptr,
1702 "Parameters should be present when writing inputrec");
1703 gmx::serializeKeyValueTree(*ir->internalParameters, serializer);
1709 static void do_harm(gmx::ISerializer *serializer, t_iparams *iparams)
1711 serializer->doReal(&iparams->harmonic.rA);
1712 serializer->doReal(&iparams->harmonic.krA);
1713 serializer->doReal(&iparams->harmonic.rB);
1714 serializer->doReal(&iparams->harmonic.krB);
1717 static void do_iparams(gmx::ISerializer *serializer,
1733 do_harm(serializer, iparams);
1734 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && serializer->reading())
1736 /* Correct incorrect storage of parameters */
1737 iparams->pdihs.phiB = iparams->pdihs.phiA;
1738 iparams->pdihs.cpB = iparams->pdihs.cpA;
1742 serializer->doReal(&iparams->harmonic.rA);
1743 serializer->doReal(&iparams->harmonic.krA);
1745 case F_LINEAR_ANGLES:
1746 serializer->doReal(&iparams->linangle.klinA);
1747 serializer->doReal(&iparams->linangle.aA);
1748 serializer->doReal(&iparams->linangle.klinB);
1749 serializer->doReal(&iparams->linangle.aB);
1752 serializer->doReal(&iparams->fene.bm);
1753 serializer->doReal(&iparams->fene.kb);
1757 serializer->doReal(&iparams->restraint.lowA);
1758 serializer->doReal(&iparams->restraint.up1A);
1759 serializer->doReal(&iparams->restraint.up2A);
1760 serializer->doReal(&iparams->restraint.kA);
1761 serializer->doReal(&iparams->restraint.lowB);
1762 serializer->doReal(&iparams->restraint.up1B);
1763 serializer->doReal(&iparams->restraint.up2B);
1764 serializer->doReal(&iparams->restraint.kB);
1770 serializer->doReal(&iparams->tab.kA);
1771 serializer->doInt(&iparams->tab.table);
1772 serializer->doReal(&iparams->tab.kB);
1774 case F_CROSS_BOND_BONDS:
1775 serializer->doReal(&iparams->cross_bb.r1e);
1776 serializer->doReal(&iparams->cross_bb.r2e);
1777 serializer->doReal(&iparams->cross_bb.krr);
1779 case F_CROSS_BOND_ANGLES:
1780 serializer->doReal(&iparams->cross_ba.r1e);
1781 serializer->doReal(&iparams->cross_ba.r2e);
1782 serializer->doReal(&iparams->cross_ba.r3e);
1783 serializer->doReal(&iparams->cross_ba.krt);
1785 case F_UREY_BRADLEY:
1786 serializer->doReal(&iparams->u_b.thetaA);
1787 serializer->doReal(&iparams->u_b.kthetaA);
1788 serializer->doReal(&iparams->u_b.r13A);
1789 serializer->doReal(&iparams->u_b.kUBA);
1790 if (file_version >= 79)
1792 serializer->doReal(&iparams->u_b.thetaB);
1793 serializer->doReal(&iparams->u_b.kthetaB);
1794 serializer->doReal(&iparams->u_b.r13B);
1795 serializer->doReal(&iparams->u_b.kUBB);
1799 iparams->u_b.thetaB = iparams->u_b.thetaA;
1800 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1801 iparams->u_b.r13B = iparams->u_b.r13A;
1802 iparams->u_b.kUBB = iparams->u_b.kUBA;
1805 case F_QUARTIC_ANGLES:
1806 serializer->doReal(&iparams->qangle.theta);
1807 serializer->doRealArray(iparams->qangle.c, 5);
1810 serializer->doReal(&iparams->bham.a);
1811 serializer->doReal(&iparams->bham.b);
1812 serializer->doReal(&iparams->bham.c);
1815 serializer->doReal(&iparams->morse.b0A);
1816 serializer->doReal(&iparams->morse.cbA);
1817 serializer->doReal(&iparams->morse.betaA);
1818 if (file_version >= 79)
1820 serializer->doReal(&iparams->morse.b0B);
1821 serializer->doReal(&iparams->morse.cbB);
1822 serializer->doReal(&iparams->morse.betaB);
1826 iparams->morse.b0B = iparams->morse.b0A;
1827 iparams->morse.cbB = iparams->morse.cbA;
1828 iparams->morse.betaB = iparams->morse.betaA;
1832 serializer->doReal(&iparams->cubic.b0);
1833 serializer->doReal(&iparams->cubic.kb);
1834 serializer->doReal(&iparams->cubic.kcub);
1838 case F_POLARIZATION:
1839 serializer->doReal(&iparams->polarize.alpha);
1842 serializer->doReal(&iparams->anharm_polarize.alpha);
1843 serializer->doReal(&iparams->anharm_polarize.drcut);
1844 serializer->doReal(&iparams->anharm_polarize.khyp);
1847 serializer->doReal(&iparams->wpol.al_x);
1848 serializer->doReal(&iparams->wpol.al_y);
1849 serializer->doReal(&iparams->wpol.al_z);
1850 serializer->doReal(&iparams->wpol.rOH);
1851 serializer->doReal(&iparams->wpol.rHH);
1852 serializer->doReal(&iparams->wpol.rOD);
1855 serializer->doReal(&iparams->thole.a);
1856 serializer->doReal(&iparams->thole.alpha1);
1857 serializer->doReal(&iparams->thole.alpha2);
1858 serializer->doReal(&iparams->thole.rfac);
1861 serializer->doReal(&iparams->lj.c6);
1862 serializer->doReal(&iparams->lj.c12);
1865 serializer->doReal(&iparams->lj14.c6A);
1866 serializer->doReal(&iparams->lj14.c12A);
1867 serializer->doReal(&iparams->lj14.c6B);
1868 serializer->doReal(&iparams->lj14.c12B);
1871 serializer->doReal(&iparams->ljc14.fqq);
1872 serializer->doReal(&iparams->ljc14.qi);
1873 serializer->doReal(&iparams->ljc14.qj);
1874 serializer->doReal(&iparams->ljc14.c6);
1875 serializer->doReal(&iparams->ljc14.c12);
1877 case F_LJC_PAIRS_NB:
1878 serializer->doReal(&iparams->ljcnb.qi);
1879 serializer->doReal(&iparams->ljcnb.qj);
1880 serializer->doReal(&iparams->ljcnb.c6);
1881 serializer->doReal(&iparams->ljcnb.c12);
1887 serializer->doReal(&iparams->pdihs.phiA);
1888 serializer->doReal(&iparams->pdihs.cpA);
1889 serializer->doReal(&iparams->pdihs.phiB);
1890 serializer->doReal(&iparams->pdihs.cpB);
1891 serializer->doInt(&iparams->pdihs.mult);
1894 serializer->doReal(&iparams->pdihs.phiA);
1895 serializer->doReal(&iparams->pdihs.cpA);
1898 serializer->doInt(&iparams->disres.label);
1899 serializer->doInt(&iparams->disres.type);
1900 serializer->doReal(&iparams->disres.low);
1901 serializer->doReal(&iparams->disres.up1);
1902 serializer->doReal(&iparams->disres.up2);
1903 serializer->doReal(&iparams->disres.kfac);
1906 serializer->doInt(&iparams->orires.ex);
1907 serializer->doInt(&iparams->orires.label);
1908 serializer->doInt(&iparams->orires.power);
1909 serializer->doReal(&iparams->orires.c);
1910 serializer->doReal(&iparams->orires.obs);
1911 serializer->doReal(&iparams->orires.kfac);
1914 if (file_version < 82)
1916 serializer->doInt(&idum);
1917 serializer->doInt(&idum);
1919 serializer->doReal(&iparams->dihres.phiA);
1920 serializer->doReal(&iparams->dihres.dphiA);
1921 serializer->doReal(&iparams->dihres.kfacA);
1922 if (file_version >= 82)
1924 serializer->doReal(&iparams->dihres.phiB);
1925 serializer->doReal(&iparams->dihres.dphiB);
1926 serializer->doReal(&iparams->dihres.kfacB);
1930 iparams->dihres.phiB = iparams->dihres.phiA;
1931 iparams->dihres.dphiB = iparams->dihres.dphiA;
1932 iparams->dihres.kfacB = iparams->dihres.kfacA;
1936 serializer->doRvec(&iparams->posres.pos0A);
1937 serializer->doRvec(&iparams->posres.fcA);
1938 serializer->doRvec(&iparams->posres.pos0B);
1939 serializer->doRvec(&iparams->posres.fcB);
1942 serializer->doInt(&iparams->fbposres.geom);
1943 serializer->doRvec(&iparams->fbposres.pos0);
1944 serializer->doReal(&iparams->fbposres.r);
1945 serializer->doReal(&iparams->fbposres.k);
1948 serializer->doRealArray(iparams->cbtdihs.cbtcA, NR_CBTDIHS);
1951 serializer->doRealArray(iparams->rbdihs.rbcA, NR_RBDIHS);
1952 serializer->doRealArray(iparams->rbdihs.rbcB, NR_RBDIHS);
1955 /* Fourier dihedrals are internally represented
1956 * as Ryckaert-Bellemans since those are faster to compute.
1958 serializer->doRealArray(iparams->rbdihs.rbcA, NR_RBDIHS);
1959 serializer->doRealArray(iparams->rbdihs.rbcB, NR_RBDIHS);
1963 serializer->doReal(&iparams->constr.dA);
1964 serializer->doReal(&iparams->constr.dB);
1967 serializer->doReal(&iparams->settle.doh);
1968 serializer->doReal(&iparams->settle.dhh);
1971 serializer->doReal(&iparams->vsite.a);
1976 serializer->doReal(&iparams->vsite.a);
1977 serializer->doReal(&iparams->vsite.b);
1982 serializer->doReal(&iparams->vsite.a);
1983 serializer->doReal(&iparams->vsite.b);
1984 serializer->doReal(&iparams->vsite.c);
1987 serializer->doInt(&iparams->vsiten.n);
1988 serializer->doReal(&iparams->vsiten.a);
1990 case F_GB12_NOLONGERUSED:
1991 case F_GB13_NOLONGERUSED:
1992 case F_GB14_NOLONGERUSED:
1993 // Implicit solvent parameters can still be read, but never used
1994 if (serializer->reading())
1996 if (file_version < 68)
1998 serializer->doReal(&rdum);
1999 serializer->doReal(&rdum);
2000 serializer->doReal(&rdum);
2001 serializer->doReal(&rdum);
2003 if (file_version < tpxv_RemoveImplicitSolvation)
2005 serializer->doReal(&rdum);
2006 serializer->doReal(&rdum);
2007 serializer->doReal(&rdum);
2008 serializer->doReal(&rdum);
2009 serializer->doReal(&rdum);
2014 serializer->doInt(&iparams->cmap.cmapA);
2015 serializer->doInt(&iparams->cmap.cmapB);
2018 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
2019 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
2023 static void do_ilist(gmx::ISerializer *serializer, InteractionList *ilist)
2025 int nr = ilist->size();
2026 serializer->doInt(&nr);
2027 if (serializer->reading())
2029 ilist->iatoms.resize(nr);
2031 serializer->doIntArray(ilist->iatoms.data(), ilist->size());
2034 static void do_ffparams(gmx::ISerializer *serializer,
2035 gmx_ffparams_t *ffparams,
2038 serializer->doInt(&ffparams->atnr);
2039 int numTypes = ffparams->numTypes();
2040 serializer->doInt(&numTypes);
2041 if (serializer->reading())
2043 ffparams->functype.resize(numTypes);
2044 ffparams->iparams.resize(numTypes);
2046 /* Read/write all the function types */
2047 serializer->doIntArray(ffparams->functype.data(), ffparams->functype.size());
2049 if (file_version >= 66)
2051 serializer->doDouble(&ffparams->reppow);
2055 ffparams->reppow = 12.0;
2058 serializer->doReal(&ffparams->fudgeQQ);
2060 /* Check whether all these function types are supported by the code.
2061 * In practice the code is backwards compatible, which means that the
2062 * numbering may have to be altered from old numbering to new numbering
2064 for (int i = 0; i < ffparams->numTypes(); i++)
2066 if (serializer->reading())
2068 /* Loop over file versions */
2069 for (int k = 0; k < NFTUPD; k++)
2071 /* Compare the read file_version to the update table */
2072 if ((file_version < ftupd[k].fvnr) &&
2073 (ffparams->functype[i] >= ftupd[k].ftype))
2075 ffparams->functype[i] += 1;
2080 do_iparams(serializer, ffparams->functype[i], &ffparams->iparams[i],
2085 static void add_settle_atoms(InteractionList *ilist)
2089 /* Settle used to only store the first atom: add the other two */
2090 ilist->iatoms.resize(2*ilist->size());
2091 for (i = ilist->size()/4 - 1; i >= 0; i--)
2093 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2094 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2095 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2096 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2100 static void do_ilists(gmx::ISerializer *serializer,
2101 InteractionLists *ilists,
2104 GMX_RELEASE_ASSERT(ilists, "Need a valid ilists object");
2105 GMX_RELEASE_ASSERT(ilists->size() == F_NRE, "The code needs to be in sync with InteractionLists");
2107 for (int j = 0; j < F_NRE; j++)
2109 InteractionList &ilist = (*ilists)[j];
2110 gmx_bool bClear = FALSE;
2111 if (serializer->reading())
2113 for (int k = 0; k < NFTUPD; k++)
2115 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2123 ilist.iatoms.clear();
2127 do_ilist(serializer, &ilist);
2128 if (file_version < 78 && j == F_SETTLE && ilist.size() > 0)
2130 add_settle_atoms(&ilist);
2136 static void do_block(gmx::ISerializer *serializer, t_block *block)
2138 serializer->doInt(&block->nr);
2139 if (serializer->reading())
2141 if ((block->nalloc_index > 0) && (nullptr != block->index))
2143 sfree(block->index);
2145 block->nalloc_index = block->nr+1;
2146 snew(block->index, block->nalloc_index);
2148 serializer->doIntArray(block->index, block->nr+1);
2151 static void do_blocka(gmx::ISerializer *serializer, t_blocka *block)
2153 serializer->doInt(&block->nr);
2154 serializer->doInt(&block->nra);
2155 if (serializer->reading())
2157 block->nalloc_index = block->nr+1;
2158 snew(block->index, block->nalloc_index);
2159 block->nalloc_a = block->nra;
2160 snew(block->a, block->nalloc_a);
2162 serializer->doIntArray(block->index, block->nr+1);
2163 serializer->doIntArray(block->a, block->nra);
2166 /* This is a primitive routine to make it possible to translate atomic numbers
2167 * to element names when reading TPR files, without making the Gromacs library
2168 * directory a dependency on mdrun (which is the case if we need elements.dat).
2171 atomicnumber_to_element(int atomicnumber)
2175 /* This does not have to be complete, so we only include elements likely
2176 * to occur in PDB files.
2178 switch (atomicnumber)
2180 case 1: p = "H"; break;
2181 case 5: p = "B"; break;
2182 case 6: p = "C"; break;
2183 case 7: p = "N"; break;
2184 case 8: p = "O"; break;
2185 case 9: p = "F"; break;
2186 case 11: p = "Na"; break;
2187 case 12: p = "Mg"; break;
2188 case 15: p = "P"; break;
2189 case 16: p = "S"; break;
2190 case 17: p = "Cl"; break;
2191 case 18: p = "Ar"; break;
2192 case 19: p = "K"; break;
2193 case 20: p = "Ca"; break;
2194 case 25: p = "Mn"; break;
2195 case 26: p = "Fe"; break;
2196 case 28: p = "Ni"; break;
2197 case 29: p = "Cu"; break;
2198 case 30: p = "Zn"; break;
2199 case 35: p = "Br"; break;
2200 case 47: p = "Ag"; break;
2201 default: p = ""; break;
2207 static void do_atom(gmx::ISerializer *serializer, t_atom *atom)
2209 serializer->doReal(&atom->m);
2210 serializer->doReal(&atom->q);
2211 serializer->doReal(&atom->mB);
2212 serializer->doReal(&atom->qB);
2213 serializer->doUShort(&atom->type);
2214 serializer->doUShort(&atom->typeB);
2215 serializer->doInt(&atom->ptype);
2216 serializer->doInt(&atom->resind);
2217 serializer->doInt(&atom->atomnumber);
2218 if (serializer->reading())
2220 /* Set element string from atomic number if present.
2221 * This routine returns an empty string if the name is not found.
2223 std::strncpy(atom->elem, atomicnumber_to_element(atom->atomnumber), 4);
2224 /* avoid warnings about potentially unterminated string */
2225 atom->elem[3] = '\0';
2229 static void do_grps(gmx::ISerializer *serializer,
2230 gmx::ArrayRef<AtomGroupIndices> grps)
2232 for (auto &group : grps)
2234 int size = group.size();
2235 serializer->doInt(&size);
2236 if (serializer->reading())
2240 serializer->doIntArray(group.data(), size);
2244 static void do_symstr(gmx::ISerializer *serializer, char ***nm, t_symtab *symtab)
2248 if (serializer->reading())
2250 serializer->doInt(&ls);
2251 *nm = get_symtab_handle(symtab, ls);
2255 ls = lookup_symtab(symtab, *nm);
2256 serializer->doInt(&ls);
2260 static void do_strstr(gmx::ISerializer *serializer,
2267 for (j = 0; (j < nstr); j++)
2269 do_symstr(serializer, &(nm[j]), symtab);
2273 static void do_resinfo(gmx::ISerializer *serializer,
2281 for (j = 0; (j < n); j++)
2283 do_symstr(serializer, &(ri[j].name), symtab);
2284 if (file_version >= 63)
2286 serializer->doInt(&ri[j].nr);
2287 serializer->doUChar(&ri[j].ic);
2297 static void do_atoms(gmx::ISerializer *serializer,
2304 serializer->doInt(&atoms->nr);
2305 serializer->doInt(&atoms->nres);
2306 if (serializer->reading())
2308 /* Since we have always written all t_atom properties in the tpr file
2309 * (at least for all backward compatible versions), we don't store
2310 * but simple set the booleans here.
2312 atoms->haveMass = TRUE;
2313 atoms->haveCharge = TRUE;
2314 atoms->haveType = TRUE;
2315 atoms->haveBState = TRUE;
2316 atoms->havePdbInfo = FALSE;
2318 snew(atoms->atom, atoms->nr);
2319 snew(atoms->atomname, atoms->nr);
2320 snew(atoms->atomtype, atoms->nr);
2321 snew(atoms->atomtypeB, atoms->nr);
2322 snew(atoms->resinfo, atoms->nres);
2323 atoms->pdbinfo = nullptr;
2327 GMX_RELEASE_ASSERT(atoms->haveMass && atoms->haveCharge && atoms->haveType && atoms->haveBState, "Mass, charge, atomtype and B-state parameters should be present in t_atoms when writing a tpr file");
2329 for (i = 0; (i < atoms->nr); i++)
2331 do_atom(serializer, &atoms->atom[i]);
2333 do_strstr(serializer, atoms->nr, atoms->atomname, symtab);
2334 do_strstr(serializer, atoms->nr, atoms->atomtype, symtab);
2335 do_strstr(serializer, atoms->nr, atoms->atomtypeB, symtab);
2337 do_resinfo(serializer, atoms->nres, atoms->resinfo, symtab, file_version);
2340 static void do_groups(gmx::ISerializer *serializer,
2341 SimulationGroups *groups,
2344 do_grps(serializer, groups->groups);
2345 int numberOfGroupNames = groups->groupNames.size();
2346 serializer->doInt(&numberOfGroupNames);
2347 if (serializer->reading())
2349 groups->groupNames.resize(numberOfGroupNames);
2351 do_strstr(serializer, numberOfGroupNames, groups->groupNames.data(), symtab);
2352 for (auto group : gmx::keysOf(groups->groupNumbers))
2354 int numberOfGroupNumbers = groups->numberOfGroupNumbers(group);
2355 serializer->doInt(&numberOfGroupNumbers);
2356 if (numberOfGroupNumbers != 0)
2358 if (serializer->reading())
2360 groups->groupNumbers[group].resize(numberOfGroupNumbers);
2362 serializer->doUCharArray(groups->groupNumbers[group].data(), numberOfGroupNumbers);
2367 static void do_atomtypes(gmx::ISerializer *serializer, t_atomtypes *atomtypes,
2372 serializer->doInt(&atomtypes->nr);
2374 if (serializer->reading())
2376 snew(atomtypes->atomnumber, j);
2378 if (serializer->reading() && file_version < tpxv_RemoveImplicitSolvation)
2380 std::vector<real> dummy(atomtypes->nr, 0);
2381 serializer->doRealArray(dummy.data(), dummy.size());
2382 serializer->doRealArray(dummy.data(), dummy.size());
2383 serializer->doRealArray(dummy.data(), dummy.size());
2385 serializer->doIntArray(atomtypes->atomnumber, j);
2387 if (serializer->reading() && file_version >= 60 && file_version < tpxv_RemoveImplicitSolvation)
2389 std::vector<real> dummy(atomtypes->nr, 0);
2390 serializer->doRealArray(dummy.data(), dummy.size());
2391 serializer->doRealArray(dummy.data(), dummy.size());
2395 static void do_symtab(gmx::ISerializer *serializer, t_symtab *symtab)
2400 serializer->doInt(&symtab->nr);
2402 if (serializer->reading())
2404 snew(symtab->symbuf, 1);
2405 symbuf = symtab->symbuf;
2406 symbuf->bufsize = nr;
2407 snew(symbuf->buf, nr);
2408 for (i = 0; (i < nr); i++)
2411 serializer->doString(&buf);
2412 symbuf->buf[i] = gmx_strdup(buf.c_str());
2417 symbuf = symtab->symbuf;
2418 while (symbuf != nullptr)
2420 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2422 std::string buf = symbuf->buf[i];
2423 serializer->doString(&buf);
2426 symbuf = symbuf->next;
2430 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2435 static void do_cmap(gmx::ISerializer *serializer, gmx_cmap_t *cmap_grid)
2438 int ngrid = cmap_grid->cmapdata.size();
2439 serializer->doInt(&ngrid);
2440 serializer->doInt(&cmap_grid->grid_spacing);
2442 int gs = cmap_grid->grid_spacing;
2443 int nelem = gs * gs;
2445 if (serializer->reading())
2447 cmap_grid->cmapdata.resize(ngrid);
2449 for (int i = 0; i < ngrid; i++)
2451 cmap_grid->cmapdata[i].cmap.resize(4*nelem);
2455 for (int i = 0; i < ngrid; i++)
2457 for (int j = 0; j < nelem; j++)
2459 serializer->doReal(&cmap_grid->cmapdata[i].cmap[j*4]);
2460 serializer->doReal(&cmap_grid->cmapdata[i].cmap[j*4+1]);
2461 serializer->doReal(&cmap_grid->cmapdata[i].cmap[j*4+2]);
2462 serializer->doReal(&cmap_grid->cmapdata[i].cmap[j*4+3]);
2468 static void do_moltype(gmx::ISerializer *serializer,
2469 gmx_moltype_t *molt,
2473 do_symstr(serializer, &(molt->name), symtab);
2475 do_atoms(serializer, &molt->atoms, symtab, file_version);
2477 do_ilists(serializer, &molt->ilist, file_version);
2479 do_block(serializer, &molt->cgs);
2481 /* This used to be in the atoms struct */
2482 do_blocka(serializer, &molt->excls);
2485 static void do_molblock(gmx::ISerializer *serializer,
2486 gmx_molblock_t *molb,
2487 int numAtomsPerMolecule)
2489 serializer->doInt(&molb->type);
2490 serializer->doInt(&molb->nmol);
2491 /* To maintain forward topology reading compatibility, we store #atoms.
2492 * TODO: Change this to conditional reading of a dummy int when we
2493 * increase tpx_generation.
2495 serializer->doInt(&numAtomsPerMolecule);
2496 /* Position restraint coordinates */
2497 int numPosres_xA = molb->posres_xA.size();
2498 serializer->doInt(&numPosres_xA);
2499 if (numPosres_xA > 0)
2501 if (serializer->reading())
2503 molb->posres_xA.resize(numPosres_xA);
2505 serializer->doRvecArray(as_rvec_array(molb->posres_xA.data()), numPosres_xA);
2507 int numPosres_xB = molb->posres_xB.size();
2508 serializer->doInt(&numPosres_xB);
2509 if (numPosres_xB > 0)
2511 if (serializer->reading())
2513 molb->posres_xB.resize(numPosres_xB);
2515 serializer->doRvecArray(as_rvec_array(molb->posres_xB.data()), numPosres_xB);
2520 static void set_disres_npair(gmx_mtop_t *mtop)
2522 gmx_mtop_ilistloop_t iloop;
2525 gmx::ArrayRef<t_iparams> ip = mtop->ffparams.iparams;
2527 iloop = gmx_mtop_ilistloop_init(mtop);
2528 while (const InteractionLists *ilist = gmx_mtop_ilistloop_next(iloop, &nmol))
2530 const InteractionList &il = (*ilist)[F_DISRES];
2534 gmx::ArrayRef<const int> a = il.iatoms;
2536 for (int i = 0; i < il.size(); i += 3)
2539 if (i+3 == il.size() || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2541 ip[a[i]].disres.npair = npair;
2549 static void do_mtop(gmx::ISerializer *serializer,
2553 do_symtab(serializer, &(mtop->symtab));
2555 do_symstr(serializer, &(mtop->name), &(mtop->symtab));
2557 do_ffparams(serializer, &mtop->ffparams, file_version);
2559 int nmoltype = mtop->moltype.size();
2560 serializer->doInt(&nmoltype);
2561 if (serializer->reading())
2563 mtop->moltype.resize(nmoltype);
2565 for (gmx_moltype_t &moltype : mtop->moltype)
2567 do_moltype(serializer, &moltype, &mtop->symtab, file_version);
2570 int nmolblock = mtop->molblock.size();
2571 serializer->doInt(&nmolblock);
2572 if (serializer->reading())
2574 mtop->molblock.resize(nmolblock);
2576 for (gmx_molblock_t &molblock : mtop->molblock)
2578 int numAtomsPerMolecule = (serializer->reading() ? 0 : mtop->moltype[molblock.type].atoms.nr);
2579 do_molblock(serializer, &molblock, numAtomsPerMolecule);
2581 serializer->doInt(&mtop->natoms);
2583 if (file_version >= tpxv_IntermolecularBondeds)
2585 serializer->doBool(&mtop->bIntermolecularInteractions);
2586 if (mtop->bIntermolecularInteractions)
2588 if (serializer->reading())
2590 mtop->intermolecular_ilist = std::make_unique<InteractionLists>();
2592 do_ilists(serializer, mtop->intermolecular_ilist.get(), file_version);
2597 mtop->bIntermolecularInteractions = FALSE;
2600 do_atomtypes (serializer, &(mtop->atomtypes), file_version);
2602 if (file_version >= 65)
2604 do_cmap(serializer, &mtop->ffparams.cmap_grid);
2608 mtop->ffparams.cmap_grid.grid_spacing = 0;
2609 mtop->ffparams.cmap_grid.cmapdata.clear();
2612 do_groups(serializer, &mtop->groups, &(mtop->symtab));
2614 mtop->haveMoleculeIndices = true;
2616 if (serializer->reading())
2618 close_symtab(&(mtop->symtab));
2623 * Read the first part of the TPR file to find general system information.
2625 * If \p TopOnlyOK is true then we can read even future versions
2626 * of tpx files, provided the \p fileGeneration hasn't changed.
2627 * If it is false, we need the \p ir too, and bail out
2628 * if the file is newer than the program.
2630 * The version and generation of the topology (see top of this file)
2631 * are returned in the two last arguments, if those arguments are non-nullptr.
2633 * If possible, we will read the \p ir even when \p TopOnlyOK is true.
2635 * \param[in,out] serializer The serializer used to handle header processing.
2636 * \param[in,out] tpx File header datastructure.
2637 * \param[in] filename The name of the file being read/written
2638 * \param[in,out] fio File handle.
2639 * \param[in] TopOnlyOK If not reading \p ir is fine or not.
2641 static void do_tpxheader(gmx::FileIOXdrSerializer *serializer,
2643 const char *filename,
2652 /* XDR binary topology file */
2653 precision = sizeof(real);
2655 std::string fileTag;
2656 if (serializer->reading())
2658 serializer->doString(&buf);
2659 if (std::strncmp(buf.c_str(), "VERSION", 7) != 0)
2661 gmx_fatal(FARGS, "Can not read file %s,\n"
2662 " this file is from a GROMACS version which is older than 2.0\n"
2663 " Make a new one with grompp or use a gro or pdb file, if possible",
2666 serializer->doInt(&precision);
2667 bDouble = (precision == sizeof(double));
2668 if ((precision != sizeof(float)) && !bDouble)
2670 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
2671 "instead of %zu or %zu",
2672 filename, precision, sizeof(float), sizeof(double));
2674 gmx_fio_setprecision(fio, bDouble);
2675 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
2676 filename, buf.c_str(), bDouble ? "double" : "single");
2680 buf = gmx::formatString("VERSION %s", gmx_version());
2681 serializer->doString(&buf);
2682 bDouble = (precision == sizeof(double));
2683 gmx_fio_setprecision(fio, bDouble);
2684 serializer->doInt(&precision);
2685 fileTag = gmx::formatString("%s", tpx_tag);
2688 /* Check versions! */
2689 serializer->doInt(&tpx->fileVersion);
2691 /* This is for backward compatibility with development versions 77-79
2692 * where the tag was, mistakenly, placed before the generation,
2693 * which would cause a segv instead of a proper error message
2694 * when reading the topology only from tpx with <77 code.
2696 if (tpx->fileVersion >= 77 && tpx->fileVersion <= 79)
2698 serializer->doString(&fileTag);
2701 serializer->doInt(&tpx->fileGeneration);
2703 if (tpx->fileVersion >= 81)
2705 serializer->doString(&fileTag);
2707 if (serializer->reading())
2709 if (tpx->fileVersion < 77)
2711 /* Versions before 77 don't have the tag, set it to release */
2712 fileTag = gmx::formatString("%s", TPX_TAG_RELEASE);
2715 if (fileTag != tpx_tag)
2717 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
2718 fileTag.c_str(), tpx_tag);
2720 /* We only support reading tpx files with the same tag as the code
2721 * or tpx files with the release tag and with lower version number.
2723 if (fileTag != TPX_TAG_RELEASE && tpx->fileVersion < tpx_version)
2725 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2726 filename, tpx->fileVersion, fileTag.c_str(),
2727 tpx_version, tpx_tag);
2732 if ((tpx->fileVersion <= tpx_incompatible_version) ||
2733 ((tpx->fileVersion > tpx_version) && !TopOnlyOK) ||
2734 (tpx->fileGeneration > tpx_generation) ||
2735 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
2737 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
2738 filename, tpx->fileVersion, tpx_version);
2741 serializer->doInt(&tpx->natoms);
2742 serializer->doInt(&tpx->ngtc);
2744 if (tpx->fileVersion < 62)
2746 serializer->doInt(&idum);
2747 serializer->doReal(&rdum);
2749 if (tpx->fileVersion >= 79)
2751 serializer->doInt(&tpx->fep_state);
2753 serializer->doReal(&tpx->lambda);
2754 serializer->doBool(&tpx->bIr);
2755 serializer->doBool(&tpx->bTop);
2756 serializer->doBool(&tpx->bX);
2757 serializer->doBool(&tpx->bV);
2758 serializer->doBool(&tpx->bF);
2759 serializer->doBool(&tpx->bBox);
2761 if ((tpx->fileGeneration > tpx_generation))
2763 /* This can only happen if TopOnlyOK=TRUE */
2768 static int do_tpx_body(gmx::ISerializer *serializer,
2777 gmx_bool bPeriodicMols;
2779 if (!serializer->reading())
2781 GMX_RELEASE_ASSERT(x == nullptr && v == nullptr, "Passing separate x and v pointers to do_tpx() is not supported when writing");
2785 GMX_RELEASE_ASSERT(!(x == nullptr && v != nullptr), "Passing x==NULL and v!=NULL is not supported");
2788 if (serializer->reading())
2791 init_gtc_state(state, tpx->ngtc, 0, 0);
2794 // v is also nullptr by the above assertion, so we may
2795 // need to make memory in state for storing the contents
2799 state->flags |= (1 << estX);
2803 state->flags |= (1 << estV);
2805 state_change_natoms(state, tpx->natoms);
2811 x = state->x.rvec_array();
2812 v = state->v.rvec_array();
2815 #define do_test(serializer, b, p) if ((serializer)->reading() && ((p) != nullptr) && !(b)) gmx_fatal(FARGS, "No %s in input file",#p)
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);
2851 do_test(serializer, tpx->bTop, mtop);
2856 do_mtop(serializer, mtop, tpx->fileVersion);
2861 do_mtop(serializer, &dum_top, tpx->fileVersion);
2864 do_test(serializer, tpx->bX, x);
2867 if (serializer->reading())
2869 state->flags |= (1<<estX);
2871 serializer->doRvecArray(x, tpx->natoms);
2874 do_test(serializer, tpx->bV, v);
2877 if (serializer->reading())
2879 state->flags |= (1<<estV);
2883 std::vector<gmx::RVec> dummyVelocities(tpx->natoms);
2884 serializer->doRvecArray(as_rvec_array(dummyVelocities.data()), tpx->natoms);
2888 serializer->doRvecArray(v, tpx->natoms);
2892 // No need to run do_test when the last argument is NULL
2895 std::vector<gmx::RVec> dummyForces(state->natoms);
2896 serializer->doRvecArray(as_rvec_array(dummyForces.data()), tpx->natoms);
2899 /* Starting with tpx version 26, we have the inputrec
2900 * at the end of the file, so we can ignore it
2901 * if the file is never than the software (but still the
2902 * same generation - see comments at the top of this file.
2907 bPeriodicMols = FALSE;
2909 do_test(serializer, tpx->bIr, ir);
2912 if (tpx->fileVersion >= 53)
2914 /* Removed the pbc info from do_inputrec, since we always want it */
2915 if (!serializer->reading())
2918 bPeriodicMols = ir->bPeriodicMols;
2920 serializer->doInt(&ePBC);
2921 serializer->doBool(&bPeriodicMols);
2923 if (tpx->fileGeneration <= tpx_generation && ir)
2925 do_inputrec(serializer, ir, tpx->fileVersion);
2926 if (tpx->fileVersion < 51)
2928 set_box_rel(ir, state);
2930 if (tpx->fileVersion < 53)
2933 bPeriodicMols = ir->bPeriodicMols;
2936 if (serializer->reading() && ir && tpx->fileVersion >= 53)
2938 /* We need to do this after do_inputrec, since that initializes ir */
2940 ir->bPeriodicMols = bPeriodicMols;
2944 if (serializer->reading())
2948 if (state->ngtc == 0)
2950 /* Reading old version without tcoupl state data: set it */
2951 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
2953 if (tpx->bTop && mtop)
2955 if (tpx->fileVersion < 57)
2957 if (mtop->moltype[0].ilist[F_DISRES].size() > 0)
2959 ir->eDisre = edrSimple;
2963 ir->eDisre = edrNone;
2966 set_disres_npair(mtop);
2970 if (tpx->bTop && mtop)
2972 gmx_mtop_finalize(mtop);
2979 static t_fileio *open_tpx(const char *fn, const char *mode)
2981 return gmx_fio_open(fn, mode);
2984 static void close_tpx(t_fileio *fio)
2990 * Fill information into the header only from state before writing.
2992 * Populating the header needs to be independent from writing the information
2993 * to file to allow things like writing the raw byte stream.
2995 * \param[in] state The current simulation state. Can't write without it.
2996 * \param[in] ir Parameter and system information.
2997 * \param[in] mtop Global topology.
2998 * \returns Fully populated header.
3000 static TpxFileHeader populateTpxHeader(const t_state &state,
3001 const t_inputrec *ir,
3002 const gmx_mtop_t *mtop)
3004 TpxFileHeader header;
3005 header.natoms = state.natoms;
3006 header.ngtc = state.ngtc;
3007 header.fep_state = state.fep_state;
3008 header.lambda = state.lambda[efptFEP];
3009 header.bIr = ir != nullptr;
3010 header.bTop = mtop != nullptr;
3011 header.bX = (state.flags & (1 << estX)) != 0;
3012 header.bV = (state.flags & (1 << estV)) != 0;
3015 header.fileVersion = tpx_version;
3016 header.fileGeneration = tpx_generation;
3021 /************************************************************
3023 * The following routines are the exported ones
3025 ************************************************************/
3027 TpxFileHeader readTpxHeader(const char *fileName, bool canReadTopologyOnly)
3031 fio = open_tpx(fileName, "r");
3032 gmx::FileIOXdrSerializer serializer(fio);
3035 do_tpxheader(&serializer, &tpx, fileName, fio, canReadTopologyOnly);
3040 void write_tpx_state(const char *fn,
3041 const t_inputrec *ir, const t_state *state, const gmx_mtop_t *mtop)
3045 fio = open_tpx(fn, "w");
3047 TpxFileHeader tpx = populateTpxHeader(*state,
3051 gmx::FileIOXdrSerializer serializer(fio);
3052 do_tpxheader(&serializer,
3057 do_tpx_body(&serializer,
3059 const_cast<t_inputrec *>(ir),
3060 const_cast<t_state *>(state), nullptr, nullptr,
3061 const_cast<gmx_mtop_t *>(mtop));
3065 void read_tpx_state(const char *fn,
3066 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3071 fio = open_tpx(fn, "r");
3072 gmx::FileIOXdrSerializer serializer(fio);
3073 do_tpxheader(&serializer,
3078 do_tpx_body(&serializer,
3088 int read_tpx(const char *fn,
3089 t_inputrec *ir, matrix box, int *natoms,
3090 rvec *x, rvec *v, gmx_mtop_t *mtop)
3097 fio = open_tpx(fn, "r");
3098 gmx::FileIOXdrSerializer serializer(fio);
3099 do_tpxheader(&serializer,
3104 ePBC = do_tpx_body(&serializer,
3106 ir, &state, x, v, mtop);
3108 if (mtop != nullptr && natoms != nullptr)
3110 *natoms = mtop->natoms;
3114 copy_mat(state.box, box);
3120 int read_tpx_top(const char *fn,
3121 t_inputrec *ir, matrix box, int *natoms,
3122 rvec *x, rvec *v, t_topology *top)
3127 ePBC = read_tpx(fn, ir, box, natoms, x, v, &mtop);
3129 *top = gmx_mtop_t_to_t_topology(&mtop, true);
3134 gmx_bool fn2bTPX(const char *file)
3136 return (efTPR == fn2ftp(file));
3139 void pr_tpxheader(FILE *fp, int indent, const char *title, const TpxFileHeader *sh)
3141 if (available(fp, sh, indent, title))
3143 indent = pr_title(fp, indent, title);
3144 pr_indent(fp, indent);
3145 fprintf(fp, "bIr = %spresent\n", sh->bIr ? "" : "not ");
3146 pr_indent(fp, indent);
3147 fprintf(fp, "bBox = %spresent\n", sh->bBox ? "" : "not ");
3148 pr_indent(fp, indent);
3149 fprintf(fp, "bTop = %spresent\n", sh->bTop ? "" : "not ");
3150 pr_indent(fp, indent);
3151 fprintf(fp, "bX = %spresent\n", sh->bX ? "" : "not ");
3152 pr_indent(fp, indent);
3153 fprintf(fp, "bV = %spresent\n", sh->bV ? "" : "not ");
3154 pr_indent(fp, indent);
3155 fprintf(fp, "bF = %spresent\n", sh->bF ? "" : "not ");
3157 pr_indent(fp, indent);
3158 fprintf(fp, "natoms = %d\n", sh->natoms);
3159 pr_indent(fp, indent);
3160 fprintf(fp, "lambda = %e\n", sh->lambda);