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
},
210 #define NFTUPD asize(ftupd)
212 /* Needed for backward compatibility */
215 /**************************************************************
217 * Now the higer level routines that do io of the structures and arrays
219 **************************************************************/
220 static void do_pullgrp_tpx_pre95(gmx::ISerializer
*serializer
,
226 serializer
->doInt(&pgrp
->nat
);
227 if (serializer
->reading())
229 snew(pgrp
->ind
, pgrp
->nat
);
231 serializer
->doIntArray(pgrp
->ind
, pgrp
->nat
);
232 serializer
->doInt(&pgrp
->nweight
);
233 if (serializer
->reading())
235 snew(pgrp
->weight
, pgrp
->nweight
);
237 serializer
->doRealArray(pgrp
->weight
, pgrp
->nweight
);
238 serializer
->doInt(&pgrp
->pbcatom
);
239 serializer
->doRvec(&pcrd
->vec
);
240 clear_rvec(pcrd
->origin
);
241 serializer
->doRvec(&tmp
);
243 serializer
->doReal(&pcrd
->rate
);
244 serializer
->doReal(&pcrd
->k
);
245 serializer
->doReal(&pcrd
->kB
);
248 static void do_pull_group(gmx::ISerializer
*serializer
, t_pull_group
*pgrp
)
250 serializer
->doInt(&pgrp
->nat
);
251 if (serializer
->reading())
253 snew(pgrp
->ind
, pgrp
->nat
);
255 serializer
->doIntArray(pgrp
->ind
, pgrp
->nat
);
256 serializer
->doInt(&pgrp
->nweight
);
257 if (serializer
->reading())
259 snew(pgrp
->weight
, pgrp
->nweight
);
261 serializer
->doRealArray(pgrp
->weight
, pgrp
->nweight
);
262 serializer
->doInt(&pgrp
->pbcatom
);
265 static void do_pull_coord(gmx::ISerializer
*serializer
,
268 int ePullOld
, int eGeomOld
, ivec dimOld
)
270 if (file_version
>= tpxv_PullCoordNGroup
)
272 serializer
->doInt(&pcrd
->eType
);
273 if (file_version
>= tpxv_PullExternalPotential
)
275 if (pcrd
->eType
== epullEXTERNAL
)
278 if (serializer
->reading())
280 serializer
->doString(&buf
);
281 pcrd
->externalPotentialProvider
= gmx_strdup(buf
.c_str());
285 buf
= pcrd
->externalPotentialProvider
;
286 serializer
->doString(&buf
);
291 pcrd
->externalPotentialProvider
= nullptr;
296 if (serializer
->reading())
298 pcrd
->externalPotentialProvider
= nullptr;
301 /* Note that we try to support adding new geometries without
302 * changing the tpx version. This requires checks when printing the
303 * geometry string and a check and fatal_error in init_pull.
305 serializer
->doInt(&pcrd
->eGeom
);
306 serializer
->doInt(&pcrd
->ngroup
);
307 if (pcrd
->ngroup
<= c_pullCoordNgroupMax
)
309 serializer
->doIntArray(pcrd
->group
, pcrd
->ngroup
);
313 /* More groups in file than supported, this must be a new geometry
314 * that is not supported by our current code. Since we will not
315 * use the groups for this coord (checks in the pull and WHAM code
316 * ensure this), we can ignore the groups and set ngroup=0.
319 snew(dum
, pcrd
->ngroup
);
320 serializer
->doIntArray(dum
, pcrd
->ngroup
);
325 serializer
->doIvec(&pcrd
->dim
);
330 serializer
->doInt(&pcrd
->group
[0]);
331 serializer
->doInt(&pcrd
->group
[1]);
332 if (file_version
>= tpxv_PullCoordTypeGeom
)
334 pcrd
->ngroup
= (pcrd
->eGeom
== epullgDIRRELATIVE
? 4 : 2);
335 serializer
->doInt(&pcrd
->eType
);
336 serializer
->doInt(&pcrd
->eGeom
);
337 if (pcrd
->ngroup
== 4)
339 serializer
->doInt(&pcrd
->group
[2]);
340 serializer
->doInt(&pcrd
->group
[3]);
342 serializer
->doIvec(&pcrd
->dim
);
346 pcrd
->eType
= ePullOld
;
347 pcrd
->eGeom
= eGeomOld
;
348 copy_ivec(dimOld
, pcrd
->dim
);
351 serializer
->doRvec(&pcrd
->origin
);
352 serializer
->doRvec(&pcrd
->vec
);
353 if (file_version
>= tpxv_PullCoordTypeGeom
)
355 serializer
->doBool(&pcrd
->bStart
);
359 /* This parameter is only printed, but not actually used by mdrun */
360 pcrd
->bStart
= FALSE
;
362 serializer
->doReal(&pcrd
->init
);
363 serializer
->doReal(&pcrd
->rate
);
364 serializer
->doReal(&pcrd
->k
);
365 serializer
->doReal(&pcrd
->kB
);
368 static void do_expandedvals(gmx::ISerializer
*serializer
,
373 int n_lambda
= fepvals
->n_lambda
;
375 /* reset the lambda calculation window */
376 fepvals
->lambda_start_n
= 0;
377 fepvals
->lambda_stop_n
= n_lambda
;
378 if (file_version
>= 79)
382 if (serializer
->reading())
384 snew(expand
->init_lambda_weights
, n_lambda
);
386 serializer
->doRealArray(expand
->init_lambda_weights
, n_lambda
);
387 serializer
->doBool(&expand
->bInit_weights
);
390 serializer
->doInt(&expand
->nstexpanded
);
391 serializer
->doInt(&expand
->elmcmove
);
392 serializer
->doInt(&expand
->elamstats
);
393 serializer
->doInt(&expand
->lmc_repeats
);
394 serializer
->doInt(&expand
->gibbsdeltalam
);
395 serializer
->doInt(&expand
->lmc_forced_nstart
);
396 serializer
->doInt(&expand
->lmc_seed
);
397 serializer
->doReal(&expand
->mc_temp
);
398 serializer
->doBool(&expand
->bSymmetrizedTMatrix
);
399 serializer
->doInt(&expand
->nstTij
);
400 serializer
->doInt(&expand
->minvarmin
);
401 serializer
->doInt(&expand
->c_range
);
402 serializer
->doReal(&expand
->wl_scale
);
403 serializer
->doReal(&expand
->wl_ratio
);
404 serializer
->doReal(&expand
->init_wl_delta
);
405 serializer
->doBool(&expand
->bWLoneovert
);
406 serializer
->doInt(&expand
->elmceq
);
407 serializer
->doInt(&expand
->equil_steps
);
408 serializer
->doInt(&expand
->equil_samples
);
409 serializer
->doInt(&expand
->equil_n_at_lam
);
410 serializer
->doReal(&expand
->equil_wl_delta
);
411 serializer
->doReal(&expand
->equil_ratio
);
415 static void do_simtempvals(gmx::ISerializer
*serializer
,
420 if (file_version
>= 79)
422 serializer
->doInt(&simtemp
->eSimTempScale
);
423 serializer
->doReal(&simtemp
->simtemp_high
);
424 serializer
->doReal(&simtemp
->simtemp_low
);
427 if (serializer
->reading())
429 snew(simtemp
->temperatures
, n_lambda
);
431 serializer
->doRealArray(simtemp
->temperatures
, n_lambda
);
436 static void do_imd(gmx::ISerializer
*serializer
, t_IMD
*imd
)
438 serializer
->doInt(&imd
->nat
);
439 if (serializer
->reading())
441 snew(imd
->ind
, imd
->nat
);
443 serializer
->doIntArray(imd
->ind
, imd
->nat
);
446 static void do_fepvals(gmx::ISerializer
*serializer
,
450 /* i is defined in the ndo_double macro; use g to iterate. */
454 /* free energy values */
456 if (file_version
>= 79)
458 serializer
->doInt(&fepvals
->init_fep_state
);
459 serializer
->doDouble(&fepvals
->init_lambda
);
460 serializer
->doDouble(&fepvals
->delta_lambda
);
462 else if (file_version
>= 59)
464 serializer
->doDouble(&fepvals
->init_lambda
);
465 serializer
->doDouble(&fepvals
->delta_lambda
);
469 serializer
->doReal(&rdum
);
470 fepvals
->init_lambda
= rdum
;
471 serializer
->doReal(&rdum
);
472 fepvals
->delta_lambda
= rdum
;
474 if (file_version
>= 79)
476 serializer
->doInt(&fepvals
->n_lambda
);
477 if (serializer
->reading())
479 snew(fepvals
->all_lambda
, efptNR
);
481 for (g
= 0; g
< efptNR
; g
++)
483 if (fepvals
->n_lambda
> 0)
485 if (serializer
->reading())
487 snew(fepvals
->all_lambda
[g
], fepvals
->n_lambda
);
489 serializer
->doDoubleArray(fepvals
->all_lambda
[g
], fepvals
->n_lambda
);
490 serializer
->doBoolArray(fepvals
->separate_dvdl
, efptNR
);
492 else if (fepvals
->init_lambda
>= 0)
494 fepvals
->separate_dvdl
[efptFEP
] = TRUE
;
498 else if (file_version
>= 64)
500 serializer
->doInt(&fepvals
->n_lambda
);
501 if (serializer
->reading())
505 snew(fepvals
->all_lambda
, efptNR
);
506 /* still allocate the all_lambda array's contents. */
507 for (g
= 0; g
< efptNR
; g
++)
509 if (fepvals
->n_lambda
> 0)
511 snew(fepvals
->all_lambda
[g
], fepvals
->n_lambda
);
515 serializer
->doDoubleArray(fepvals
->all_lambda
[efptFEP
], fepvals
->n_lambda
);
516 if (fepvals
->init_lambda
>= 0)
520 fepvals
->separate_dvdl
[efptFEP
] = TRUE
;
522 if (serializer
->reading())
524 /* copy the contents of the efptFEP lambda component to all
525 the other components */
526 for (g
= 0; g
< efptNR
; g
++)
528 for (h
= 0; h
< fepvals
->n_lambda
; h
++)
532 fepvals
->all_lambda
[g
][h
] =
533 fepvals
->all_lambda
[efptFEP
][h
];
542 fepvals
->n_lambda
= 0;
543 fepvals
->all_lambda
= nullptr;
544 if (fepvals
->init_lambda
>= 0)
546 fepvals
->separate_dvdl
[efptFEP
] = TRUE
;
549 serializer
->doReal(&fepvals
->sc_alpha
);
550 serializer
->doInt(&fepvals
->sc_power
);
551 if (file_version
>= 79)
553 serializer
->doReal(&fepvals
->sc_r_power
);
557 fepvals
->sc_r_power
= 6.0;
559 serializer
->doReal(&fepvals
->sc_sigma
);
560 if (serializer
->reading())
562 if (file_version
>= 71)
564 fepvals
->sc_sigma_min
= fepvals
->sc_sigma
;
568 fepvals
->sc_sigma_min
= 0;
571 if (file_version
>= 79)
573 serializer
->doBool(&fepvals
->bScCoul
);
577 fepvals
->bScCoul
= TRUE
;
579 if (file_version
>= 64)
581 serializer
->doInt(&fepvals
->nstdhdl
);
585 fepvals
->nstdhdl
= 1;
588 if (file_version
>= 73)
590 serializer
->doInt(&fepvals
->separate_dhdl_file
);
591 serializer
->doInt(&fepvals
->dhdl_derivatives
);
595 fepvals
->separate_dhdl_file
= esepdhdlfileYES
;
596 fepvals
->dhdl_derivatives
= edhdlderivativesYES
;
598 if (file_version
>= 71)
600 serializer
->doInt(&fepvals
->dh_hist_size
);
601 serializer
->doDouble(&fepvals
->dh_hist_spacing
);
605 fepvals
->dh_hist_size
= 0;
606 fepvals
->dh_hist_spacing
= 0.1;
608 if (file_version
>= 79)
610 serializer
->doInt(&fepvals
->edHdLPrintEnergy
);
614 fepvals
->edHdLPrintEnergy
= edHdLPrintEnergyNO
;
617 /* handle lambda_neighbors */
618 if ((file_version
>= 83 && file_version
< 90) || file_version
>= 92)
620 serializer
->doInt(&fepvals
->lambda_neighbors
);
621 if ( (fepvals
->lambda_neighbors
>= 0) && (fepvals
->init_fep_state
>= 0) &&
622 (fepvals
->init_lambda
< 0) )
624 fepvals
->lambda_start_n
= (fepvals
->init_fep_state
-
625 fepvals
->lambda_neighbors
);
626 fepvals
->lambda_stop_n
= (fepvals
->init_fep_state
+
627 fepvals
->lambda_neighbors
+ 1);
628 if (fepvals
->lambda_start_n
< 0)
630 fepvals
->lambda_start_n
= 0;;
632 if (fepvals
->lambda_stop_n
>= fepvals
->n_lambda
)
634 fepvals
->lambda_stop_n
= fepvals
->n_lambda
;
639 fepvals
->lambda_start_n
= 0;
640 fepvals
->lambda_stop_n
= fepvals
->n_lambda
;
645 fepvals
->lambda_start_n
= 0;
646 fepvals
->lambda_stop_n
= fepvals
->n_lambda
;
650 static void do_awhBias(gmx::ISerializer
*serializer
, gmx::AwhBiasParams
*awhBiasParams
,
651 int gmx_unused file_version
)
653 serializer
->doInt(&awhBiasParams
->eTarget
);
654 serializer
->doDouble(&awhBiasParams
->targetBetaScaling
);
655 serializer
->doDouble(&awhBiasParams
->targetCutoff
);
656 serializer
->doInt(&awhBiasParams
->eGrowth
);
657 serializer
->doInt(&awhBiasParams
->bUserData
);
658 serializer
->doDouble(&awhBiasParams
->errorInitial
);
659 serializer
->doInt(&awhBiasParams
->ndim
);
660 serializer
->doInt(&awhBiasParams
->shareGroup
);
661 serializer
->doBool(&awhBiasParams
->equilibrateHistogram
);
663 if (serializer
->reading())
665 snew(awhBiasParams
->dimParams
, awhBiasParams
->ndim
);
668 for (int d
= 0; d
< awhBiasParams
->ndim
; d
++)
670 gmx::AwhDimParams
*dimParams
= &awhBiasParams
->dimParams
[d
];
672 serializer
->doInt(&dimParams
->eCoordProvider
);
673 serializer
->doInt(&dimParams
->coordIndex
);
674 serializer
->doDouble(&dimParams
->origin
);
675 serializer
->doDouble(&dimParams
->end
);
676 serializer
->doDouble(&dimParams
->period
);
677 serializer
->doDouble(&dimParams
->forceConstant
);
678 serializer
->doDouble(&dimParams
->diffusion
);
679 serializer
->doDouble(&dimParams
->coordValueInit
);
680 serializer
->doDouble(&dimParams
->coverDiameter
);
684 static void do_awh(gmx::ISerializer
*serializer
,
685 gmx::AwhParams
*awhParams
,
686 int gmx_unused file_version
)
688 serializer
->doInt(&awhParams
->numBias
);
689 serializer
->doInt(&awhParams
->nstOut
);
690 serializer
->doInt64(&awhParams
->seed
);
691 serializer
->doInt(&awhParams
->nstSampleCoord
);
692 serializer
->doInt(&awhParams
->numSamplesUpdateFreeEnergy
);
693 serializer
->doInt(&awhParams
->ePotential
);
694 serializer
->doBool(&awhParams
->shareBiasMultisim
);
696 if (awhParams
->numBias
> 0)
698 if (serializer
->reading())
700 snew(awhParams
->awhBiasParams
, awhParams
->numBias
);
703 for (int k
= 0; k
< awhParams
->numBias
; k
++)
705 do_awhBias(serializer
, &awhParams
->awhBiasParams
[k
], file_version
);
710 static void do_pull(gmx::ISerializer
*serializer
,
719 if (file_version
>= 95)
721 serializer
->doInt(&pull
->ngroup
);
723 serializer
->doInt(&pull
->ncoord
);
724 if (file_version
< 95)
726 pull
->ngroup
= pull
->ncoord
+ 1;
728 if (file_version
< tpxv_PullCoordTypeGeom
)
732 serializer
->doInt(&eGeomOld
);
733 serializer
->doIvec(&dimOld
);
734 /* The inner cylinder radius, now removed */
735 serializer
->doReal(&dum
);
737 serializer
->doReal(&pull
->cylinder_r
);
738 serializer
->doReal(&pull
->constr_tol
);
739 if (file_version
>= 95)
741 serializer
->doBool(&pull
->bPrintCOM
);
742 /* With file_version < 95 this value is set below */
744 if (file_version
>= tpxv_ReplacePullPrintCOM12
)
746 serializer
->doBool(&pull
->bPrintRefValue
);
747 serializer
->doBool(&pull
->bPrintComp
);
749 else if (file_version
>= tpxv_PullCoordTypeGeom
)
752 serializer
->doInt(&idum
); /* used to be bPrintCOM2 */
753 serializer
->doBool(&pull
->bPrintRefValue
);
754 serializer
->doBool(&pull
->bPrintComp
);
758 pull
->bPrintRefValue
= FALSE
;
759 pull
->bPrintComp
= TRUE
;
761 serializer
->doInt(&pull
->nstxout
);
762 serializer
->doInt(&pull
->nstfout
);
763 if (file_version
>= tpxv_PullPrevStepCOMAsReference
)
765 serializer
->doBool(&pull
->bSetPbcRefToPrevStepCOM
);
769 pull
->bSetPbcRefToPrevStepCOM
= FALSE
;
771 if (serializer
->reading())
773 snew(pull
->group
, pull
->ngroup
);
774 snew(pull
->coord
, pull
->ncoord
);
776 if (file_version
< 95)
778 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
779 if (eGeomOld
== epullgDIRPBC
)
781 gmx_fatal(FARGS
, "pull-geometry=position is no longer supported");
783 if (eGeomOld
> epullgDIRPBC
)
788 for (g
= 0; g
< pull
->ngroup
; g
++)
790 /* We read and ignore a pull coordinate for group 0 */
791 do_pullgrp_tpx_pre95(serializer
, &pull
->group
[g
], &pull
->coord
[std::max(g
-1, 0)]);
794 pull
->coord
[g
-1].group
[0] = 0;
795 pull
->coord
[g
-1].group
[1] = g
;
799 pull
->bPrintCOM
= (pull
->group
[0].nat
> 0);
803 for (g
= 0; g
< pull
->ngroup
; g
++)
805 do_pull_group(serializer
, &pull
->group
[g
]);
807 for (g
= 0; g
< pull
->ncoord
; g
++)
809 do_pull_coord(serializer
, &pull
->coord
[g
],
810 file_version
, ePullOld
, eGeomOld
, dimOld
);
813 if (file_version
>= tpxv_PullAverage
)
817 v
= pull
->bXOutAverage
;
818 serializer
->doBool(&v
);
819 pull
->bXOutAverage
= v
;
820 v
= pull
->bFOutAverage
;
821 serializer
->doBool(&v
);
822 pull
->bFOutAverage
= v
;
827 static void do_rotgrp(gmx::ISerializer
*serializer
,
830 serializer
->doInt(&rotg
->eType
);
831 serializer
->doInt(&rotg
->bMassW
);
832 serializer
->doInt(&rotg
->nat
);
833 if (serializer
->reading())
835 snew(rotg
->ind
, rotg
->nat
);
837 serializer
->doIntArray(rotg
->ind
, rotg
->nat
);
838 if (serializer
->reading())
840 snew(rotg
->x_ref
, rotg
->nat
);
842 serializer
->doRvecArray(rotg
->x_ref
, rotg
->nat
);
843 serializer
->doRvec(&rotg
->inputVec
);
844 serializer
->doRvec(&rotg
->pivot
);
845 serializer
->doReal(&rotg
->rate
);
846 serializer
->doReal(&rotg
->k
);
847 serializer
->doReal(&rotg
->slab_dist
);
848 serializer
->doReal(&rotg
->min_gaussian
);
849 serializer
->doReal(&rotg
->eps
);
850 serializer
->doInt(&rotg
->eFittype
);
851 serializer
->doInt(&rotg
->PotAngle_nstep
);
852 serializer
->doReal(&rotg
->PotAngle_step
);
855 static void do_rot(gmx::ISerializer
*serializer
,
860 serializer
->doInt(&rot
->ngrp
);
861 serializer
->doInt(&rot
->nstrout
);
862 serializer
->doInt(&rot
->nstsout
);
863 if (serializer
->reading())
865 snew(rot
->grp
, rot
->ngrp
);
867 for (g
= 0; g
< rot
->ngrp
; g
++)
869 do_rotgrp(serializer
, &rot
->grp
[g
]);
874 static void do_swapgroup(gmx::ISerializer
*serializer
,
878 /* Name of the group or molecule */
880 if (serializer
->reading())
882 serializer
->doString(&buf
);
883 g
->molname
= gmx_strdup(buf
.c_str());
888 serializer
->doString(&buf
);
891 /* Number of atoms in the group */
892 serializer
->doInt(&g
->nat
);
894 /* The group's atom indices */
895 if (serializer
->reading())
897 snew(g
->ind
, g
->nat
);
899 serializer
->doIntArray(g
->ind
, g
->nat
);
901 /* Requested counts for compartments A and B */
902 serializer
->doIntArray(g
->nmolReq
, eCompNR
);
905 static void do_swapcoords_tpx(gmx::ISerializer
*serializer
,
909 /* Enums for better readability of the code */
914 eChannel0
= 0, eChannel1
918 if (file_version
>= tpxv_CompElPolyatomicIonsAndMultipleIonTypes
)
920 /* The total number of swap groups is the sum of the fixed groups
921 * (split0, split1, solvent), and the user-defined groups (2+ types of ions)
923 serializer
->doInt(&swap
->ngrp
);
924 if (serializer
->reading())
926 snew(swap
->grp
, swap
->ngrp
);
928 for (int ig
= 0; ig
< swap
->ngrp
; ig
++)
930 do_swapgroup(serializer
, &swap
->grp
[ig
]);
932 serializer
->doBool(&swap
->massw_split
[eChannel0
]);
933 serializer
->doBool(&swap
->massw_split
[eChannel1
]);
934 serializer
->doInt(&swap
->nstswap
);
935 serializer
->doInt(&swap
->nAverage
);
936 serializer
->doReal(&swap
->threshold
);
937 serializer
->doReal(&swap
->cyl0r
);
938 serializer
->doReal(&swap
->cyl0u
);
939 serializer
->doReal(&swap
->cyl0l
);
940 serializer
->doReal(&swap
->cyl1r
);
941 serializer
->doReal(&swap
->cyl1u
);
942 serializer
->doReal(&swap
->cyl1l
);
946 /*** Support reading older CompEl .tpr files ***/
948 /* In the original CompEl .tpr files, we always have 5 groups: */
950 snew(swap
->grp
, swap
->ngrp
);
952 swap
->grp
[eGrpSplit0
].molname
= gmx_strdup("split0" ); // group 0: split0
953 swap
->grp
[eGrpSplit1
].molname
= gmx_strdup("split1" ); // group 1: split1
954 swap
->grp
[eGrpSolvent
].molname
= gmx_strdup("solvent"); // group 2: solvent
955 swap
->grp
[3 ].molname
= gmx_strdup("anions" ); // group 3: anions
956 swap
->grp
[4 ].molname
= gmx_strdup("cations"); // group 4: cations
958 serializer
->doInt(&swap
->grp
[3].nat
);
959 serializer
->doInt(&swap
->grp
[eGrpSolvent
].nat
);
960 serializer
->doInt(&swap
->grp
[eGrpSplit0
].nat
);
961 serializer
->doBool(&swap
->massw_split
[eChannel0
]);
962 serializer
->doInt(&swap
->grp
[eGrpSplit1
].nat
);
963 serializer
->doBool(&swap
->massw_split
[eChannel1
]);
964 serializer
->doInt(&swap
->nstswap
);
965 serializer
->doInt(&swap
->nAverage
);
966 serializer
->doReal(&swap
->threshold
);
967 serializer
->doReal(&swap
->cyl0r
);
968 serializer
->doReal(&swap
->cyl0u
);
969 serializer
->doReal(&swap
->cyl0l
);
970 serializer
->doReal(&swap
->cyl1r
);
971 serializer
->doReal(&swap
->cyl1u
);
972 serializer
->doReal(&swap
->cyl1l
);
974 // The order[] array keeps compatibility with older .tpr files
975 // by reading in the groups in the classic order
977 const int order
[4] = {3, eGrpSolvent
, eGrpSplit0
, eGrpSplit1
};
979 for (int ig
= 0; ig
< 4; ig
++)
982 snew(swap
->grp
[g
].ind
, swap
->grp
[g
].nat
);
983 serializer
->doIntArray(swap
->grp
[g
].ind
, swap
->grp
[g
].nat
);
987 for (int j
= eCompA
; j
<= eCompB
; j
++)
989 serializer
->doInt(&swap
->grp
[3].nmolReq
[j
]); // group 3 = anions
990 serializer
->doInt(&swap
->grp
[4].nmolReq
[j
]); // group 4 = cations
992 } /* End support reading older CompEl .tpr files */
994 if (file_version
>= tpxv_CompElWithSwapLayerOffset
)
996 serializer
->doReal(&swap
->bulkOffset
[eCompA
]);
997 serializer
->doReal(&swap
->bulkOffset
[eCompB
]);
1002 static void do_legacy_efield(gmx::ISerializer
*serializer
, gmx::KeyValueTreeObjectBuilder
*root
)
1004 const char *const dimName
[] = { "x", "y", "z" };
1006 auto appliedForcesObj
= root
->addObject("applied-forces");
1007 auto efieldObj
= appliedForcesObj
.addObject("electric-field");
1008 // The content of the tpr file for this feature has
1009 // been the same since gromacs 4.0 that was used for
1011 for (int j
= 0; j
< DIM
; ++j
)
1014 serializer
->doInt(&n
);
1015 serializer
->doInt(&nt
);
1016 std::vector
<real
> aa(n
+1), phi(nt
+1), at(nt
+1), phit(nt
+1);
1017 serializer
->doRealArray(aa
.data(), n
);
1018 serializer
->doRealArray(phi
.data(), n
);
1019 serializer
->doRealArray(at
.data(), nt
);
1020 serializer
->doRealArray(phit
.data(), nt
);
1023 if (n
> 1 || nt
> 1)
1025 gmx_fatal(FARGS
, "Can not handle tpr files with more than one electric field term per direction.");
1027 auto dimObj
= efieldObj
.addObject(dimName
[j
]);
1028 dimObj
.addValue
<real
>("E0", aa
[0]);
1029 dimObj
.addValue
<real
>("omega", at
[0]);
1030 dimObj
.addValue
<real
>("t0", phi
[0]);
1031 dimObj
.addValue
<real
>("sigma", phit
[0]);
1037 static void do_inputrec(gmx::ISerializer
*serializer
,
1041 int i
, j
, k
, idum
= 0;
1043 gmx_bool bdum
= false;
1045 if (file_version
!= tpx_version
)
1047 /* Give a warning about features that are not accessible */
1048 fprintf(stderr
, "Note: file tpx version %d, software tpx version %d\n",
1049 file_version
, tpx_version
);
1052 if (file_version
== 0)
1057 gmx::KeyValueTreeBuilder paramsBuilder
;
1058 gmx::KeyValueTreeObjectBuilder paramsObj
= paramsBuilder
.rootObject();
1060 /* Basic inputrec stuff */
1061 serializer
->doInt(&ir
->eI
);
1062 if (file_version
>= 62)
1064 serializer
->doInt64(&ir
->nsteps
);
1068 serializer
->doInt(&idum
);
1072 if (file_version
>= 62)
1074 serializer
->doInt64(&ir
->init_step
);
1078 serializer
->doInt(&idum
);
1079 ir
->init_step
= idum
;
1082 serializer
->doInt(&ir
->simulation_part
);
1084 if (file_version
>= 67)
1086 serializer
->doInt(&ir
->nstcalcenergy
);
1090 ir
->nstcalcenergy
= 1;
1092 if (file_version
>= 81)
1094 serializer
->doInt(&ir
->cutoff_scheme
);
1095 if (file_version
< 94)
1097 ir
->cutoff_scheme
= 1 - ir
->cutoff_scheme
;
1102 ir
->cutoff_scheme
= ecutsGROUP
;
1104 serializer
->doInt(&ir
->ns_type
);
1105 serializer
->doInt(&ir
->nstlist
);
1106 serializer
->doInt(&idum
); /* used to be ndelta; not used anymore */
1108 serializer
->doReal(&ir
->rtpi
);
1110 serializer
->doInt(&ir
->nstcomm
);
1111 serializer
->doInt(&ir
->comm_mode
);
1113 /* ignore nstcheckpoint */
1114 if (file_version
< tpxv_RemoveObsoleteParameters1
)
1116 serializer
->doInt(&idum
);
1119 serializer
->doInt(&ir
->nstcgsteep
);
1121 serializer
->doInt(&ir
->nbfgscorr
);
1123 serializer
->doInt(&ir
->nstlog
);
1124 serializer
->doInt(&ir
->nstxout
);
1125 serializer
->doInt(&ir
->nstvout
);
1126 serializer
->doInt(&ir
->nstfout
);
1127 serializer
->doInt(&ir
->nstenergy
);
1128 serializer
->doInt(&ir
->nstxout_compressed
);
1129 if (file_version
>= 59)
1131 serializer
->doDouble(&ir
->init_t
);
1132 serializer
->doDouble(&ir
->delta_t
);
1136 serializer
->doReal(&rdum
);
1138 serializer
->doReal(&rdum
);
1141 serializer
->doReal(&ir
->x_compression_precision
);
1142 if (file_version
>= 81)
1144 serializer
->doReal(&ir
->verletbuf_tol
);
1148 ir
->verletbuf_tol
= 0;
1150 serializer
->doReal(&ir
->rlist
);
1151 if (file_version
>= 67 && file_version
< tpxv_RemoveTwinRange
)
1153 if (serializer
->reading())
1155 // Reading such a file version could invoke the twin-range
1156 // scheme, about which mdrun should give a fatal error.
1157 real dummy_rlistlong
= -1;
1158 serializer
->doReal(&dummy_rlistlong
);
1160 if (ir
->rlist
> 0 && (dummy_rlistlong
== 0 || dummy_rlistlong
> ir
->rlist
))
1162 // Get mdrun to issue an error (regardless of
1163 // ir->cutoff_scheme).
1164 ir
->useTwinRange
= true;
1168 // grompp used to set rlistlong actively. Users were
1169 // probably also confused and set rlistlong == rlist.
1170 // However, in all remaining cases, it is safe to let
1171 // mdrun proceed normally.
1172 ir
->useTwinRange
= false;
1178 // No need to read or write anything
1179 ir
->useTwinRange
= false;
1181 if (file_version
>= 82 && file_version
!= 90)
1183 // Multiple time-stepping is no longer enabled, but the old
1184 // support required the twin-range scheme, for which mdrun
1185 // already emits a fatal error.
1186 int dummy_nstcalclr
= -1;
1187 serializer
->doInt(&dummy_nstcalclr
);
1189 serializer
->doInt(&ir
->coulombtype
);
1190 if (file_version
>= 81)
1192 serializer
->doInt(&ir
->coulomb_modifier
);
1196 ir
->coulomb_modifier
= (ir
->cutoff_scheme
== ecutsVERLET
? eintmodPOTSHIFT
: eintmodNONE
);
1198 serializer
->doReal(&ir
->rcoulomb_switch
);
1199 serializer
->doReal(&ir
->rcoulomb
);
1200 serializer
->doInt(&ir
->vdwtype
);
1201 if (file_version
>= 81)
1203 serializer
->doInt(&ir
->vdw_modifier
);
1207 ir
->vdw_modifier
= (ir
->cutoff_scheme
== ecutsVERLET
? eintmodPOTSHIFT
: eintmodNONE
);
1209 serializer
->doReal(&ir
->rvdw_switch
);
1210 serializer
->doReal(&ir
->rvdw
);
1211 serializer
->doInt(&ir
->eDispCorr
);
1212 serializer
->doReal(&ir
->epsilon_r
);
1213 serializer
->doReal(&ir
->epsilon_rf
);
1214 serializer
->doReal(&ir
->tabext
);
1216 // This permits reading a .tpr file that used implicit solvent,
1217 // and later permitting mdrun to refuse to run it.
1218 if (serializer
->reading())
1220 if (file_version
< tpxv_RemoveImplicitSolvation
)
1222 serializer
->doInt(&idum
);
1223 serializer
->doInt(&idum
);
1224 serializer
->doReal(&rdum
);
1225 serializer
->doReal(&rdum
);
1226 serializer
->doInt(&idum
);
1227 ir
->implicit_solvent
= (idum
> 0);
1231 ir
->implicit_solvent
= false;
1233 if (file_version
< tpxv_RemoveImplicitSolvation
)
1235 serializer
->doReal(&rdum
);
1236 serializer
->doReal(&rdum
);
1237 serializer
->doReal(&rdum
);
1238 serializer
->doReal(&rdum
);
1239 if (file_version
>= 60)
1241 serializer
->doReal(&rdum
);
1242 serializer
->doInt(&idum
);
1244 serializer
->doReal(&rdum
);
1248 if (file_version
>= 81)
1250 serializer
->doReal(&ir
->fourier_spacing
);
1254 ir
->fourier_spacing
= 0.0;
1256 serializer
->doInt(&ir
->nkx
);
1257 serializer
->doInt(&ir
->nky
);
1258 serializer
->doInt(&ir
->nkz
);
1259 serializer
->doInt(&ir
->pme_order
);
1260 serializer
->doReal(&ir
->ewald_rtol
);
1262 if (file_version
>= 93)
1264 serializer
->doReal(&ir
->ewald_rtol_lj
);
1268 ir
->ewald_rtol_lj
= ir
->ewald_rtol
;
1270 serializer
->doInt(&ir
->ewald_geometry
);
1271 serializer
->doReal(&ir
->epsilon_surface
);
1273 /* ignore bOptFFT */
1274 if (file_version
< tpxv_RemoveObsoleteParameters1
)
1276 serializer
->doBool(&bdum
);
1279 if (file_version
>= 93)
1281 serializer
->doInt(&ir
->ljpme_combination_rule
);
1283 serializer
->doBool(&ir
->bContinuation
);
1284 serializer
->doInt(&ir
->etc
);
1285 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1286 * but the values 0 and 1 still mean no and
1287 * berendsen temperature coupling, respectively.
1289 if (file_version
>= 79)
1291 serializer
->doBool(&ir
->bPrintNHChains
);
1293 if (file_version
>= 71)
1295 serializer
->doInt(&ir
->nsttcouple
);
1299 ir
->nsttcouple
= ir
->nstcalcenergy
;
1301 serializer
->doInt(&ir
->epc
);
1302 serializer
->doInt(&ir
->epct
);
1303 if (file_version
>= 71)
1305 serializer
->doInt(&ir
->nstpcouple
);
1309 ir
->nstpcouple
= ir
->nstcalcenergy
;
1311 serializer
->doReal(&ir
->tau_p
);
1312 serializer
->doRvec(&ir
->ref_p
[XX
]);
1313 serializer
->doRvec(&ir
->ref_p
[YY
]);
1314 serializer
->doRvec(&ir
->ref_p
[ZZ
]);
1315 serializer
->doRvec(&ir
->compress
[XX
]);
1316 serializer
->doRvec(&ir
->compress
[YY
]);
1317 serializer
->doRvec(&ir
->compress
[ZZ
]);
1318 serializer
->doInt(&ir
->refcoord_scaling
);
1319 serializer
->doRvec(&ir
->posres_com
);
1320 serializer
->doRvec(&ir
->posres_comB
);
1322 if (file_version
< 79)
1324 serializer
->doInt(&ir
->andersen_seed
);
1328 ir
->andersen_seed
= 0;
1331 serializer
->doReal(&ir
->shake_tol
);
1333 serializer
->doInt(&ir
->efep
);
1334 do_fepvals(serializer
, ir
->fepvals
, file_version
);
1336 if (file_version
>= 79)
1338 serializer
->doBool(&ir
->bSimTemp
);
1341 ir
->bSimTemp
= TRUE
;
1346 ir
->bSimTemp
= FALSE
;
1350 do_simtempvals(serializer
, ir
->simtempvals
, ir
->fepvals
->n_lambda
, file_version
);
1353 if (file_version
>= 79)
1355 serializer
->doBool(&ir
->bExpanded
);
1358 ir
->bExpanded
= TRUE
;
1362 ir
->bExpanded
= FALSE
;
1367 do_expandedvals(serializer
, ir
->expandedvals
, ir
->fepvals
, file_version
);
1370 serializer
->doInt(&ir
->eDisre
);
1371 serializer
->doInt(&ir
->eDisreWeighting
);
1372 serializer
->doBool(&ir
->bDisreMixed
);
1373 serializer
->doReal(&ir
->dr_fc
);
1374 serializer
->doReal(&ir
->dr_tau
);
1375 serializer
->doInt(&ir
->nstdisreout
);
1376 serializer
->doReal(&ir
->orires_fc
);
1377 serializer
->doReal(&ir
->orires_tau
);
1378 serializer
->doInt(&ir
->nstorireout
);
1380 /* ignore dihre_fc */
1381 if (file_version
< 79)
1383 serializer
->doReal(&rdum
);
1386 serializer
->doReal(&ir
->em_stepsize
);
1387 serializer
->doReal(&ir
->em_tol
);
1388 serializer
->doBool(&ir
->bShakeSOR
);
1389 serializer
->doInt(&ir
->niter
);
1390 serializer
->doReal(&ir
->fc_stepsize
);
1391 serializer
->doInt(&ir
->eConstrAlg
);
1392 serializer
->doInt(&ir
->nProjOrder
);
1393 serializer
->doReal(&ir
->LincsWarnAngle
);
1394 serializer
->doInt(&ir
->nLincsIter
);
1395 serializer
->doReal(&ir
->bd_fric
);
1396 if (file_version
>= tpxv_Use64BitRandomSeed
)
1398 serializer
->doInt64(&ir
->ld_seed
);
1402 serializer
->doInt(&idum
);
1406 for (i
= 0; i
< DIM
; i
++)
1408 serializer
->doRvec(&ir
->deform
[i
]);
1410 serializer
->doReal(&ir
->cos_accel
);
1412 serializer
->doInt(&ir
->userint1
);
1413 serializer
->doInt(&ir
->userint2
);
1414 serializer
->doInt(&ir
->userint3
);
1415 serializer
->doInt(&ir
->userint4
);
1416 serializer
->doReal(&ir
->userreal1
);
1417 serializer
->doReal(&ir
->userreal2
);
1418 serializer
->doReal(&ir
->userreal3
);
1419 serializer
->doReal(&ir
->userreal4
);
1421 /* AdResS is removed, but we need to be able to read old files,
1422 and let mdrun refuse to run them */
1423 if (file_version
>= 77 && file_version
< tpxv_RemoveAdress
)
1425 serializer
->doBool(&ir
->bAdress
);
1428 int idum
, numThermoForceGroups
, numEnergyGroups
;
1431 serializer
->doInt(&idum
);
1432 serializer
->doReal(&rdum
);
1433 serializer
->doReal(&rdum
);
1434 serializer
->doReal(&rdum
);
1435 serializer
->doInt(&idum
);
1436 serializer
->doInt(&idum
);
1437 serializer
->doRvec(&rvecdum
);
1438 serializer
->doInt(&numThermoForceGroups
);
1439 serializer
->doReal(&rdum
);
1440 serializer
->doInt(&numEnergyGroups
);
1441 serializer
->doInt(&idum
);
1443 if (numThermoForceGroups
> 0)
1445 std::vector
<int> idumn(numThermoForceGroups
);
1446 serializer
->doIntArray(idumn
.data(), idumn
.size());
1448 if (numEnergyGroups
> 0)
1450 std::vector
<int> idumn(numEnergyGroups
);
1451 serializer
->doIntArray(idumn
.data(), idumn
.size());
1457 ir
->bAdress
= FALSE
;
1464 if (file_version
>= tpxv_PullCoordTypeGeom
)
1466 serializer
->doBool(&ir
->bPull
);
1470 serializer
->doInt(&ePullOld
);
1471 ir
->bPull
= (ePullOld
> 0);
1472 /* We removed the first ePull=ePullNo for the enum */
1477 if (serializer
->reading())
1481 do_pull(serializer
, ir
->pull
, file_version
, ePullOld
);
1485 if (file_version
>= tpxv_AcceleratedWeightHistogram
)
1487 serializer
->doBool(&ir
->bDoAwh
);
1491 if (serializer
->reading())
1493 snew(ir
->awhParams
, 1);
1495 do_awh(serializer
, ir
->awhParams
, file_version
);
1503 /* Enforced rotation */
1504 if (file_version
>= 74)
1506 serializer
->doBool(&ir
->bRot
);
1509 if (serializer
->reading())
1513 do_rot(serializer
, ir
->rot
);
1521 /* Interactive molecular dynamics */
1522 if (file_version
>= tpxv_InteractiveMolecularDynamics
)
1524 serializer
->doBool(&ir
->bIMD
);
1527 if (serializer
->reading())
1531 do_imd(serializer
, ir
->imd
);
1536 /* We don't support IMD sessions for old .tpr files */
1541 serializer
->doInt(&ir
->opts
.ngtc
);
1542 if (file_version
>= 69)
1544 serializer
->doInt(&ir
->opts
.nhchainlength
);
1548 ir
->opts
.nhchainlength
= 1;
1550 serializer
->doInt(&ir
->opts
.ngacc
);
1551 serializer
->doInt(&ir
->opts
.ngfrz
);
1552 serializer
->doInt(&ir
->opts
.ngener
);
1554 if (serializer
->reading())
1556 snew(ir
->opts
.nrdf
, ir
->opts
.ngtc
);
1557 snew(ir
->opts
.ref_t
, ir
->opts
.ngtc
);
1558 snew(ir
->opts
.annealing
, ir
->opts
.ngtc
);
1559 snew(ir
->opts
.anneal_npoints
, ir
->opts
.ngtc
);
1560 snew(ir
->opts
.anneal_time
, ir
->opts
.ngtc
);
1561 snew(ir
->opts
.anneal_temp
, ir
->opts
.ngtc
);
1562 snew(ir
->opts
.tau_t
, ir
->opts
.ngtc
);
1563 snew(ir
->opts
.nFreeze
, ir
->opts
.ngfrz
);
1564 snew(ir
->opts
.acc
, ir
->opts
.ngacc
);
1565 snew(ir
->opts
.egp_flags
, ir
->opts
.ngener
*ir
->opts
.ngener
);
1567 if (ir
->opts
.ngtc
> 0)
1569 serializer
->doRealArray(ir
->opts
.nrdf
, ir
->opts
.ngtc
);
1570 serializer
->doRealArray(ir
->opts
.ref_t
, ir
->opts
.ngtc
);
1571 serializer
->doRealArray(ir
->opts
.tau_t
, ir
->opts
.ngtc
);
1573 if (ir
->opts
.ngfrz
> 0)
1575 serializer
->doIvecArray(ir
->opts
.nFreeze
, ir
->opts
.ngfrz
);
1577 if (ir
->opts
.ngacc
> 0)
1579 serializer
->doRvecArray(ir
->opts
.acc
, ir
->opts
.ngacc
);
1581 serializer
->doIntArray(ir
->opts
.egp_flags
,
1582 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
);
1632 serializer
->doBool(&ir
->bQMMM
);
1633 serializer
->doInt(&ir
->QMMMscheme
);
1634 serializer
->doReal(&ir
->scalefactor
);
1635 serializer
->doInt(&ir
->opts
.ngQM
);
1636 if (serializer
->reading())
1638 snew(ir
->opts
.QMmethod
, ir
->opts
.ngQM
);
1639 snew(ir
->opts
.QMbasis
, ir
->opts
.ngQM
);
1640 snew(ir
->opts
.QMcharge
, ir
->opts
.ngQM
);
1641 snew(ir
->opts
.QMmult
, ir
->opts
.ngQM
);
1642 snew(ir
->opts
.bSH
, ir
->opts
.ngQM
);
1643 snew(ir
->opts
.CASorbitals
, ir
->opts
.ngQM
);
1644 snew(ir
->opts
.CASelectrons
, ir
->opts
.ngQM
);
1645 snew(ir
->opts
.SAon
, ir
->opts
.ngQM
);
1646 snew(ir
->opts
.SAoff
, ir
->opts
.ngQM
);
1647 snew(ir
->opts
.SAsteps
, ir
->opts
.ngQM
);
1649 if (ir
->opts
.ngQM
> 0 && ir
->bQMMM
)
1651 serializer
->doIntArray(ir
->opts
.QMmethod
, ir
->opts
.ngQM
);
1652 serializer
->doIntArray(ir
->opts
.QMbasis
, ir
->opts
.ngQM
);
1653 serializer
->doIntArray(ir
->opts
.QMcharge
, ir
->opts
.ngQM
);
1654 serializer
->doIntArray(ir
->opts
.QMmult
, ir
->opts
.ngQM
);
1655 serializer
->doBoolArray(ir
->opts
.bSH
, ir
->opts
.ngQM
);
1656 serializer
->doIntArray(ir
->opts
.CASorbitals
, ir
->opts
.ngQM
);
1657 serializer
->doIntArray(ir
->opts
.CASelectrons
, ir
->opts
.ngQM
);
1658 serializer
->doRealArray(ir
->opts
.SAon
, ir
->opts
.ngQM
);
1659 serializer
->doRealArray(ir
->opts
.SAoff
, ir
->opts
.ngQM
);
1660 serializer
->doIntArray(ir
->opts
.SAsteps
, ir
->opts
.ngQM
);
1661 /* We leave in dummy i/o for removed parameters to avoid
1662 * changing the tpr format for every QMMM change.
1664 std::vector
<int> dummy(ir
->opts
.ngQM
, 0);
1665 serializer
->doIntArray(dummy
.data(), ir
->opts
.ngQM
);
1666 serializer
->doIntArray(dummy
.data(), ir
->opts
.ngQM
);
1668 /* end of QMMM stuff */
1671 if (file_version
>= tpxv_GenericParamsForElectricField
)
1673 if (serializer
->reading())
1675 paramsObj
.mergeObject(
1676 gmx::deserializeKeyValueTree(serializer
));
1680 GMX_RELEASE_ASSERT(ir
->params
!= nullptr,
1681 "Parameters should be present when writing inputrec");
1682 gmx::serializeKeyValueTree(*ir
->params
, serializer
);
1685 if (serializer
->reading())
1687 ir
->params
= new gmx::KeyValueTreeObject(paramsBuilder
.build());
1688 // Initialize internal parameters to an empty kvt for all tpr versions
1689 ir
->internalParameters
= std::make_unique
<gmx::KeyValueTreeObject
>();
1692 if (file_version
>= tpxv_GenericInternalParameters
)
1694 if (serializer
->reading())
1696 ir
->internalParameters
= 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
,
1732 do_harm(serializer
, iparams
);
1733 if ((ftype
== F_ANGRES
|| ftype
== F_ANGRESZ
) && serializer
->reading())
1735 /* Correct incorrect storage of parameters */
1736 iparams
->pdihs
.phiB
= iparams
->pdihs
.phiA
;
1737 iparams
->pdihs
.cpB
= iparams
->pdihs
.cpA
;
1741 serializer
->doReal(&iparams
->harmonic
.rA
);
1742 serializer
->doReal(&iparams
->harmonic
.krA
);
1744 case F_LINEAR_ANGLES
:
1745 serializer
->doReal(&iparams
->linangle
.klinA
);
1746 serializer
->doReal(&iparams
->linangle
.aA
);
1747 serializer
->doReal(&iparams
->linangle
.klinB
);
1748 serializer
->doReal(&iparams
->linangle
.aB
);
1751 serializer
->doReal(&iparams
->fene
.bm
);
1752 serializer
->doReal(&iparams
->fene
.kb
);
1756 serializer
->doReal(&iparams
->restraint
.lowA
);
1757 serializer
->doReal(&iparams
->restraint
.up1A
);
1758 serializer
->doReal(&iparams
->restraint
.up2A
);
1759 serializer
->doReal(&iparams
->restraint
.kA
);
1760 serializer
->doReal(&iparams
->restraint
.lowB
);
1761 serializer
->doReal(&iparams
->restraint
.up1B
);
1762 serializer
->doReal(&iparams
->restraint
.up2B
);
1763 serializer
->doReal(&iparams
->restraint
.kB
);
1769 serializer
->doReal(&iparams
->tab
.kA
);
1770 serializer
->doInt(&iparams
->tab
.table
);
1771 serializer
->doReal(&iparams
->tab
.kB
);
1773 case F_CROSS_BOND_BONDS
:
1774 serializer
->doReal(&iparams
->cross_bb
.r1e
);
1775 serializer
->doReal(&iparams
->cross_bb
.r2e
);
1776 serializer
->doReal(&iparams
->cross_bb
.krr
);
1778 case F_CROSS_BOND_ANGLES
:
1779 serializer
->doReal(&iparams
->cross_ba
.r1e
);
1780 serializer
->doReal(&iparams
->cross_ba
.r2e
);
1781 serializer
->doReal(&iparams
->cross_ba
.r3e
);
1782 serializer
->doReal(&iparams
->cross_ba
.krt
);
1784 case F_UREY_BRADLEY
:
1785 serializer
->doReal(&iparams
->u_b
.thetaA
);
1786 serializer
->doReal(&iparams
->u_b
.kthetaA
);
1787 serializer
->doReal(&iparams
->u_b
.r13A
);
1788 serializer
->doReal(&iparams
->u_b
.kUBA
);
1789 if (file_version
>= 79)
1791 serializer
->doReal(&iparams
->u_b
.thetaB
);
1792 serializer
->doReal(&iparams
->u_b
.kthetaB
);
1793 serializer
->doReal(&iparams
->u_b
.r13B
);
1794 serializer
->doReal(&iparams
->u_b
.kUBB
);
1798 iparams
->u_b
.thetaB
= iparams
->u_b
.thetaA
;
1799 iparams
->u_b
.kthetaB
= iparams
->u_b
.kthetaA
;
1800 iparams
->u_b
.r13B
= iparams
->u_b
.r13A
;
1801 iparams
->u_b
.kUBB
= iparams
->u_b
.kUBA
;
1804 case F_QUARTIC_ANGLES
:
1805 serializer
->doReal(&iparams
->qangle
.theta
);
1806 serializer
->doRealArray(iparams
->qangle
.c
, 5);
1809 serializer
->doReal(&iparams
->bham
.a
);
1810 serializer
->doReal(&iparams
->bham
.b
);
1811 serializer
->doReal(&iparams
->bham
.c
);
1814 serializer
->doReal(&iparams
->morse
.b0A
);
1815 serializer
->doReal(&iparams
->morse
.cbA
);
1816 serializer
->doReal(&iparams
->morse
.betaA
);
1817 if (file_version
>= 79)
1819 serializer
->doReal(&iparams
->morse
.b0B
);
1820 serializer
->doReal(&iparams
->morse
.cbB
);
1821 serializer
->doReal(&iparams
->morse
.betaB
);
1825 iparams
->morse
.b0B
= iparams
->morse
.b0A
;
1826 iparams
->morse
.cbB
= iparams
->morse
.cbA
;
1827 iparams
->morse
.betaB
= iparams
->morse
.betaA
;
1831 serializer
->doReal(&iparams
->cubic
.b0
);
1832 serializer
->doReal(&iparams
->cubic
.kb
);
1833 serializer
->doReal(&iparams
->cubic
.kcub
);
1837 case F_POLARIZATION
:
1838 serializer
->doReal(&iparams
->polarize
.alpha
);
1841 serializer
->doReal(&iparams
->anharm_polarize
.alpha
);
1842 serializer
->doReal(&iparams
->anharm_polarize
.drcut
);
1843 serializer
->doReal(&iparams
->anharm_polarize
.khyp
);
1846 serializer
->doReal(&iparams
->wpol
.al_x
);
1847 serializer
->doReal(&iparams
->wpol
.al_y
);
1848 serializer
->doReal(&iparams
->wpol
.al_z
);
1849 serializer
->doReal(&iparams
->wpol
.rOH
);
1850 serializer
->doReal(&iparams
->wpol
.rHH
);
1851 serializer
->doReal(&iparams
->wpol
.rOD
);
1854 serializer
->doReal(&iparams
->thole
.a
);
1855 serializer
->doReal(&iparams
->thole
.alpha1
);
1856 serializer
->doReal(&iparams
->thole
.alpha2
);
1857 serializer
->doReal(&iparams
->thole
.rfac
);
1860 serializer
->doReal(&iparams
->lj
.c6
);
1861 serializer
->doReal(&iparams
->lj
.c12
);
1864 serializer
->doReal(&iparams
->lj14
.c6A
);
1865 serializer
->doReal(&iparams
->lj14
.c12A
);
1866 serializer
->doReal(&iparams
->lj14
.c6B
);
1867 serializer
->doReal(&iparams
->lj14
.c12B
);
1870 serializer
->doReal(&iparams
->ljc14
.fqq
);
1871 serializer
->doReal(&iparams
->ljc14
.qi
);
1872 serializer
->doReal(&iparams
->ljc14
.qj
);
1873 serializer
->doReal(&iparams
->ljc14
.c6
);
1874 serializer
->doReal(&iparams
->ljc14
.c12
);
1876 case F_LJC_PAIRS_NB
:
1877 serializer
->doReal(&iparams
->ljcnb
.qi
);
1878 serializer
->doReal(&iparams
->ljcnb
.qj
);
1879 serializer
->doReal(&iparams
->ljcnb
.c6
);
1880 serializer
->doReal(&iparams
->ljcnb
.c12
);
1886 serializer
->doReal(&iparams
->pdihs
.phiA
);
1887 serializer
->doReal(&iparams
->pdihs
.cpA
);
1888 serializer
->doReal(&iparams
->pdihs
.phiB
);
1889 serializer
->doReal(&iparams
->pdihs
.cpB
);
1890 serializer
->doInt(&iparams
->pdihs
.mult
);
1893 serializer
->doReal(&iparams
->pdihs
.phiA
);
1894 serializer
->doReal(&iparams
->pdihs
.cpA
);
1897 serializer
->doInt(&iparams
->disres
.label
);
1898 serializer
->doInt(&iparams
->disres
.type
);
1899 serializer
->doReal(&iparams
->disres
.low
);
1900 serializer
->doReal(&iparams
->disres
.up1
);
1901 serializer
->doReal(&iparams
->disres
.up2
);
1902 serializer
->doReal(&iparams
->disres
.kfac
);
1905 serializer
->doInt(&iparams
->orires
.ex
);
1906 serializer
->doInt(&iparams
->orires
.label
);
1907 serializer
->doInt(&iparams
->orires
.power
);
1908 serializer
->doReal(&iparams
->orires
.c
);
1909 serializer
->doReal(&iparams
->orires
.obs
);
1910 serializer
->doReal(&iparams
->orires
.kfac
);
1913 if (file_version
< 82)
1915 serializer
->doInt(&idum
);
1916 serializer
->doInt(&idum
);
1918 serializer
->doReal(&iparams
->dihres
.phiA
);
1919 serializer
->doReal(&iparams
->dihres
.dphiA
);
1920 serializer
->doReal(&iparams
->dihres
.kfacA
);
1921 if (file_version
>= 82)
1923 serializer
->doReal(&iparams
->dihres
.phiB
);
1924 serializer
->doReal(&iparams
->dihres
.dphiB
);
1925 serializer
->doReal(&iparams
->dihres
.kfacB
);
1929 iparams
->dihres
.phiB
= iparams
->dihres
.phiA
;
1930 iparams
->dihres
.dphiB
= iparams
->dihres
.dphiA
;
1931 iparams
->dihres
.kfacB
= iparams
->dihres
.kfacA
;
1935 serializer
->doRvec(&iparams
->posres
.pos0A
);
1936 serializer
->doRvec(&iparams
->posres
.fcA
);
1937 serializer
->doRvec(&iparams
->posres
.pos0B
);
1938 serializer
->doRvec(&iparams
->posres
.fcB
);
1941 serializer
->doInt(&iparams
->fbposres
.geom
);
1942 serializer
->doRvec(&iparams
->fbposres
.pos0
);
1943 serializer
->doReal(&iparams
->fbposres
.r
);
1944 serializer
->doReal(&iparams
->fbposres
.k
);
1947 serializer
->doRealArray(iparams
->cbtdihs
.cbtcA
, NR_CBTDIHS
);
1950 serializer
->doRealArray(iparams
->rbdihs
.rbcA
, NR_RBDIHS
);
1951 serializer
->doRealArray(iparams
->rbdihs
.rbcB
, NR_RBDIHS
);
1954 /* Fourier dihedrals are internally represented
1955 * as Ryckaert-Bellemans since those are faster to compute.
1957 serializer
->doRealArray(iparams
->rbdihs
.rbcA
, NR_RBDIHS
);
1958 serializer
->doRealArray(iparams
->rbdihs
.rbcB
, NR_RBDIHS
);
1962 serializer
->doReal(&iparams
->constr
.dA
);
1963 serializer
->doReal(&iparams
->constr
.dB
);
1966 serializer
->doReal(&iparams
->settle
.doh
);
1967 serializer
->doReal(&iparams
->settle
.dhh
);
1970 serializer
->doReal(&iparams
->vsite
.a
);
1975 serializer
->doReal(&iparams
->vsite
.a
);
1976 serializer
->doReal(&iparams
->vsite
.b
);
1981 serializer
->doReal(&iparams
->vsite
.a
);
1982 serializer
->doReal(&iparams
->vsite
.b
);
1983 serializer
->doReal(&iparams
->vsite
.c
);
1986 serializer
->doInt(&iparams
->vsiten
.n
);
1987 serializer
->doReal(&iparams
->vsiten
.a
);
1989 case F_GB12_NOLONGERUSED
:
1990 case F_GB13_NOLONGERUSED
:
1991 case F_GB14_NOLONGERUSED
:
1992 // Implicit solvent parameters can still be read, but never used
1993 if (serializer
->reading())
1995 if (file_version
< 68)
1997 serializer
->doReal(&rdum
);
1998 serializer
->doReal(&rdum
);
1999 serializer
->doReal(&rdum
);
2000 serializer
->doReal(&rdum
);
2002 if (file_version
< tpxv_RemoveImplicitSolvation
)
2004 serializer
->doReal(&rdum
);
2005 serializer
->doReal(&rdum
);
2006 serializer
->doReal(&rdum
);
2007 serializer
->doReal(&rdum
);
2008 serializer
->doReal(&rdum
);
2013 serializer
->doInt(&iparams
->cmap
.cmapA
);
2014 serializer
->doInt(&iparams
->cmap
.cmapB
);
2017 gmx_fatal(FARGS
, "unknown function type %d (%s) in %s line %d",
2018 ftype
, interaction_function
[ftype
].name
, __FILE__
, __LINE__
);
2022 static void do_ilist(gmx::ISerializer
*serializer
, InteractionList
*ilist
)
2024 int nr
= ilist
->size();
2025 serializer
->doInt(&nr
);
2026 if (serializer
->reading())
2028 ilist
->iatoms
.resize(nr
);
2030 serializer
->doIntArray(ilist
->iatoms
.data(), ilist
->size());
2033 static void do_ffparams(gmx::ISerializer
*serializer
,
2034 gmx_ffparams_t
*ffparams
,
2037 serializer
->doInt(&ffparams
->atnr
);
2038 int numTypes
= ffparams
->numTypes();
2039 serializer
->doInt(&numTypes
);
2040 if (serializer
->reading())
2042 ffparams
->functype
.resize(numTypes
);
2043 ffparams
->iparams
.resize(numTypes
);
2045 /* Read/write all the function types */
2046 serializer
->doIntArray(ffparams
->functype
.data(), ffparams
->functype
.size());
2048 if (file_version
>= 66)
2050 serializer
->doDouble(&ffparams
->reppow
);
2054 ffparams
->reppow
= 12.0;
2057 serializer
->doReal(&ffparams
->fudgeQQ
);
2059 /* Check whether all these function types are supported by the code.
2060 * In practice the code is backwards compatible, which means that the
2061 * numbering may have to be altered from old numbering to new numbering
2063 for (int i
= 0; i
< ffparams
->numTypes(); i
++)
2065 if (serializer
->reading())
2067 /* Loop over file versions */
2068 for (int k
= 0; k
< NFTUPD
; k
++)
2070 /* Compare the read file_version to the update table */
2071 if ((file_version
< ftupd
[k
].fvnr
) &&
2072 (ffparams
->functype
[i
] >= ftupd
[k
].ftype
))
2074 ffparams
->functype
[i
] += 1;
2079 do_iparams(serializer
, ffparams
->functype
[i
], &ffparams
->iparams
[i
],
2084 static void add_settle_atoms(InteractionList
*ilist
)
2088 /* Settle used to only store the first atom: add the other two */
2089 ilist
->iatoms
.resize(2*ilist
->size());
2090 for (i
= ilist
->size()/4 - 1; i
>= 0; i
--)
2092 ilist
->iatoms
[4*i
+0] = ilist
->iatoms
[2*i
+0];
2093 ilist
->iatoms
[4*i
+1] = ilist
->iatoms
[2*i
+1];
2094 ilist
->iatoms
[4*i
+2] = ilist
->iatoms
[2*i
+1] + 1;
2095 ilist
->iatoms
[4*i
+3] = ilist
->iatoms
[2*i
+1] + 2;
2099 static void do_ilists(gmx::ISerializer
*serializer
,
2100 InteractionLists
*ilists
,
2103 GMX_RELEASE_ASSERT(ilists
, "Need a valid ilists object");
2104 GMX_RELEASE_ASSERT(ilists
->size() == F_NRE
, "The code needs to be in sync with InteractionLists");
2106 for (int j
= 0; j
< F_NRE
; j
++)
2108 InteractionList
&ilist
= (*ilists
)[j
];
2109 gmx_bool bClear
= FALSE
;
2110 if (serializer
->reading())
2112 for (int k
= 0; k
< NFTUPD
; k
++)
2114 if ((file_version
< ftupd
[k
].fvnr
) && (j
== ftupd
[k
].ftype
))
2122 ilist
.iatoms
.clear();
2126 do_ilist(serializer
, &ilist
);
2127 if (file_version
< 78 && j
== F_SETTLE
&& ilist
.size() > 0)
2129 add_settle_atoms(&ilist
);
2135 static void do_block(gmx::ISerializer
*serializer
, t_block
*block
)
2137 serializer
->doInt(&block
->nr
);
2138 if (serializer
->reading())
2140 if ((block
->nalloc_index
> 0) && (nullptr != block
->index
))
2142 sfree(block
->index
);
2144 block
->nalloc_index
= block
->nr
+1;
2145 snew(block
->index
, block
->nalloc_index
);
2147 serializer
->doIntArray(block
->index
, block
->nr
+1);
2150 static void do_blocka(gmx::ISerializer
*serializer
, t_blocka
*block
)
2152 serializer
->doInt(&block
->nr
);
2153 serializer
->doInt(&block
->nra
);
2154 if (serializer
->reading())
2156 block
->nalloc_index
= block
->nr
+1;
2157 snew(block
->index
, block
->nalloc_index
);
2158 block
->nalloc_a
= block
->nra
;
2159 snew(block
->a
, block
->nalloc_a
);
2161 serializer
->doIntArray(block
->index
, block
->nr
+1);
2162 serializer
->doIntArray(block
->a
, block
->nra
);
2165 /* This is a primitive routine to make it possible to translate atomic numbers
2166 * to element names when reading TPR files, without making the Gromacs library
2167 * directory a dependency on mdrun (which is the case if we need elements.dat).
2170 atomicnumber_to_element(int atomicnumber
)
2174 /* This does not have to be complete, so we only include elements likely
2175 * to occur in PDB files.
2177 switch (atomicnumber
)
2179 case 1: p
= "H"; break;
2180 case 5: p
= "B"; break;
2181 case 6: p
= "C"; break;
2182 case 7: p
= "N"; break;
2183 case 8: p
= "O"; break;
2184 case 9: p
= "F"; break;
2185 case 11: p
= "Na"; break;
2186 case 12: p
= "Mg"; break;
2187 case 15: p
= "P"; break;
2188 case 16: p
= "S"; break;
2189 case 17: p
= "Cl"; break;
2190 case 18: p
= "Ar"; break;
2191 case 19: p
= "K"; break;
2192 case 20: p
= "Ca"; break;
2193 case 25: p
= "Mn"; break;
2194 case 26: p
= "Fe"; break;
2195 case 28: p
= "Ni"; break;
2196 case 29: p
= "Cu"; break;
2197 case 30: p
= "Zn"; break;
2198 case 35: p
= "Br"; break;
2199 case 47: p
= "Ag"; break;
2200 default: p
= ""; break;
2206 static void do_atom(gmx::ISerializer
*serializer
, t_atom
*atom
)
2208 serializer
->doReal(&atom
->m
);
2209 serializer
->doReal(&atom
->q
);
2210 serializer
->doReal(&atom
->mB
);
2211 serializer
->doReal(&atom
->qB
);
2212 serializer
->doUShort(&atom
->type
);
2213 serializer
->doUShort(&atom
->typeB
);
2214 serializer
->doInt(&atom
->ptype
);
2215 serializer
->doInt(&atom
->resind
);
2216 serializer
->doInt(&atom
->atomnumber
);
2217 if (serializer
->reading())
2219 /* Set element string from atomic number if present.
2220 * This routine returns an empty string if the name is not found.
2222 std::strncpy(atom
->elem
, atomicnumber_to_element(atom
->atomnumber
), 4);
2223 /* avoid warnings about potentially unterminated string */
2224 atom
->elem
[3] = '\0';
2228 static void do_grps(gmx::ISerializer
*serializer
,
2229 gmx::ArrayRef
<AtomGroupIndices
> grps
)
2231 for (auto &group
: grps
)
2233 int size
= group
.size();
2234 serializer
->doInt(&size
);
2235 if (serializer
->reading())
2239 serializer
->doIntArray(group
.data(), size
);
2243 static void do_symstr(gmx::ISerializer
*serializer
, char ***nm
, t_symtab
*symtab
)
2247 if (serializer
->reading())
2249 serializer
->doInt(&ls
);
2250 *nm
= get_symtab_handle(symtab
, ls
);
2254 ls
= lookup_symtab(symtab
, *nm
);
2255 serializer
->doInt(&ls
);
2259 static void do_strstr(gmx::ISerializer
*serializer
,
2266 for (j
= 0; (j
< nstr
); j
++)
2268 do_symstr(serializer
, &(nm
[j
]), symtab
);
2272 static void do_resinfo(gmx::ISerializer
*serializer
,
2280 for (j
= 0; (j
< n
); j
++)
2282 do_symstr(serializer
, &(ri
[j
].name
), symtab
);
2283 if (file_version
>= 63)
2285 serializer
->doInt(&ri
[j
].nr
);
2286 serializer
->doUChar(&ri
[j
].ic
);
2296 static void do_atoms(gmx::ISerializer
*serializer
,
2303 serializer
->doInt(&atoms
->nr
);
2304 serializer
->doInt(&atoms
->nres
);
2305 if (serializer
->reading())
2307 /* Since we have always written all t_atom properties in the tpr file
2308 * (at least for all backward compatible versions), we don't store
2309 * but simple set the booleans here.
2311 atoms
->haveMass
= TRUE
;
2312 atoms
->haveCharge
= TRUE
;
2313 atoms
->haveType
= TRUE
;
2314 atoms
->haveBState
= TRUE
;
2315 atoms
->havePdbInfo
= FALSE
;
2317 snew(atoms
->atom
, atoms
->nr
);
2318 snew(atoms
->atomname
, atoms
->nr
);
2319 snew(atoms
->atomtype
, atoms
->nr
);
2320 snew(atoms
->atomtypeB
, atoms
->nr
);
2321 snew(atoms
->resinfo
, atoms
->nres
);
2322 atoms
->pdbinfo
= nullptr;
2326 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");
2328 for (i
= 0; (i
< atoms
->nr
); i
++)
2330 do_atom(serializer
, &atoms
->atom
[i
]);
2332 do_strstr(serializer
, atoms
->nr
, atoms
->atomname
, symtab
);
2333 do_strstr(serializer
, atoms
->nr
, atoms
->atomtype
, symtab
);
2334 do_strstr(serializer
, atoms
->nr
, atoms
->atomtypeB
, symtab
);
2336 do_resinfo(serializer
, atoms
->nres
, atoms
->resinfo
, symtab
, file_version
);
2339 static void do_groups(gmx::ISerializer
*serializer
,
2340 SimulationGroups
*groups
,
2343 do_grps(serializer
, groups
->groups
);
2344 int numberOfGroupNames
= groups
->groupNames
.size();
2345 serializer
->doInt(&numberOfGroupNames
);
2346 if (serializer
->reading())
2348 groups
->groupNames
.resize(numberOfGroupNames
);
2350 do_strstr(serializer
, numberOfGroupNames
, groups
->groupNames
.data(), symtab
);
2351 for (auto group
: gmx::keysOf(groups
->groupNumbers
))
2353 int numberOfGroupNumbers
= groups
->numberOfGroupNumbers(group
);
2354 serializer
->doInt(&numberOfGroupNumbers
);
2355 if (numberOfGroupNumbers
!= 0)
2357 if (serializer
->reading())
2359 groups
->groupNumbers
[group
].resize(numberOfGroupNumbers
);
2361 serializer
->doUCharArray(groups
->groupNumbers
[group
].data(), numberOfGroupNumbers
);
2366 static void do_atomtypes(gmx::ISerializer
*serializer
, t_atomtypes
*atomtypes
,
2371 serializer
->doInt(&atomtypes
->nr
);
2373 if (serializer
->reading())
2375 snew(atomtypes
->atomnumber
, j
);
2377 if (serializer
->reading() && file_version
< tpxv_RemoveImplicitSolvation
)
2379 std::vector
<real
> dummy(atomtypes
->nr
, 0);
2380 serializer
->doRealArray(dummy
.data(), dummy
.size());
2381 serializer
->doRealArray(dummy
.data(), dummy
.size());
2382 serializer
->doRealArray(dummy
.data(), dummy
.size());
2384 serializer
->doIntArray(atomtypes
->atomnumber
, j
);
2386 if (serializer
->reading() && file_version
>= 60 && file_version
< tpxv_RemoveImplicitSolvation
)
2388 std::vector
<real
> dummy(atomtypes
->nr
, 0);
2389 serializer
->doRealArray(dummy
.data(), dummy
.size());
2390 serializer
->doRealArray(dummy
.data(), dummy
.size());
2394 static void do_symtab(gmx::ISerializer
*serializer
, t_symtab
*symtab
)
2399 serializer
->doInt(&symtab
->nr
);
2401 if (serializer
->reading())
2403 snew(symtab
->symbuf
, 1);
2404 symbuf
= symtab
->symbuf
;
2405 symbuf
->bufsize
= nr
;
2406 snew(symbuf
->buf
, nr
);
2407 for (i
= 0; (i
< nr
); i
++)
2410 serializer
->doString(&buf
);
2411 symbuf
->buf
[i
] = gmx_strdup(buf
.c_str());
2416 symbuf
= symtab
->symbuf
;
2417 while (symbuf
!= nullptr)
2419 for (i
= 0; (i
< symbuf
->bufsize
) && (i
< nr
); i
++)
2421 std::string buf
= symbuf
->buf
[i
];
2422 serializer
->doString(&buf
);
2425 symbuf
= symbuf
->next
;
2429 gmx_fatal(FARGS
, "nr of symtab strings left: %d", nr
);
2434 static void do_cmap(gmx::ISerializer
*serializer
, gmx_cmap_t
*cmap_grid
)
2437 int ngrid
= cmap_grid
->cmapdata
.size();
2438 serializer
->doInt(&ngrid
);
2439 serializer
->doInt(&cmap_grid
->grid_spacing
);
2441 int gs
= cmap_grid
->grid_spacing
;
2442 int nelem
= gs
* gs
;
2444 if (serializer
->reading())
2446 cmap_grid
->cmapdata
.resize(ngrid
);
2448 for (int i
= 0; i
< ngrid
; i
++)
2450 cmap_grid
->cmapdata
[i
].cmap
.resize(4*nelem
);
2454 for (int i
= 0; i
< ngrid
; i
++)
2456 for (int j
= 0; j
< nelem
; j
++)
2458 serializer
->doReal(&cmap_grid
->cmapdata
[i
].cmap
[j
*4]);
2459 serializer
->doReal(&cmap_grid
->cmapdata
[i
].cmap
[j
*4+1]);
2460 serializer
->doReal(&cmap_grid
->cmapdata
[i
].cmap
[j
*4+2]);
2461 serializer
->doReal(&cmap_grid
->cmapdata
[i
].cmap
[j
*4+3]);
2467 static void do_moltype(gmx::ISerializer
*serializer
,
2468 gmx_moltype_t
*molt
,
2472 do_symstr(serializer
, &(molt
->name
), symtab
);
2474 do_atoms(serializer
, &molt
->atoms
, symtab
, file_version
);
2476 do_ilists(serializer
, &molt
->ilist
, file_version
);
2478 do_block(serializer
, &molt
->cgs
);
2480 /* This used to be in the atoms struct */
2481 do_blocka(serializer
, &molt
->excls
);
2484 static void do_molblock(gmx::ISerializer
*serializer
,
2485 gmx_molblock_t
*molb
,
2486 int numAtomsPerMolecule
)
2488 serializer
->doInt(&molb
->type
);
2489 serializer
->doInt(&molb
->nmol
);
2490 /* To maintain forward topology reading compatibility, we store #atoms.
2491 * TODO: Change this to conditional reading of a dummy int when we
2492 * increase tpx_generation.
2494 serializer
->doInt(&numAtomsPerMolecule
);
2495 /* Position restraint coordinates */
2496 int numPosres_xA
= molb
->posres_xA
.size();
2497 serializer
->doInt(&numPosres_xA
);
2498 if (numPosres_xA
> 0)
2500 if (serializer
->reading())
2502 molb
->posres_xA
.resize(numPosres_xA
);
2504 serializer
->doRvecArray(as_rvec_array(molb
->posres_xA
.data()), numPosres_xA
);
2506 int numPosres_xB
= molb
->posres_xB
.size();
2507 serializer
->doInt(&numPosres_xB
);
2508 if (numPosres_xB
> 0)
2510 if (serializer
->reading())
2512 molb
->posres_xB
.resize(numPosres_xB
);
2514 serializer
->doRvecArray(as_rvec_array(molb
->posres_xB
.data()), numPosres_xB
);
2519 static void set_disres_npair(gmx_mtop_t
*mtop
)
2521 gmx_mtop_ilistloop_t iloop
;
2524 gmx::ArrayRef
<t_iparams
> ip
= mtop
->ffparams
.iparams
;
2526 iloop
= gmx_mtop_ilistloop_init(mtop
);
2527 while (const InteractionLists
*ilist
= gmx_mtop_ilistloop_next(iloop
, &nmol
))
2529 const InteractionList
&il
= (*ilist
)[F_DISRES
];
2533 gmx::ArrayRef
<const int> a
= il
.iatoms
;
2535 for (int i
= 0; i
< il
.size(); i
+= 3)
2538 if (i
+3 == il
.size() || ip
[a
[i
]].disres
.label
!= ip
[a
[i
+3]].disres
.label
)
2540 ip
[a
[i
]].disres
.npair
= npair
;
2548 static void do_mtop(gmx::ISerializer
*serializer
,
2552 do_symtab(serializer
, &(mtop
->symtab
));
2554 do_symstr(serializer
, &(mtop
->name
), &(mtop
->symtab
));
2556 do_ffparams(serializer
, &mtop
->ffparams
, file_version
);
2558 int nmoltype
= mtop
->moltype
.size();
2559 serializer
->doInt(&nmoltype
);
2560 if (serializer
->reading())
2562 mtop
->moltype
.resize(nmoltype
);
2564 for (gmx_moltype_t
&moltype
: mtop
->moltype
)
2566 do_moltype(serializer
, &moltype
, &mtop
->symtab
, file_version
);
2569 int nmolblock
= mtop
->molblock
.size();
2570 serializer
->doInt(&nmolblock
);
2571 if (serializer
->reading())
2573 mtop
->molblock
.resize(nmolblock
);
2575 for (gmx_molblock_t
&molblock
: mtop
->molblock
)
2577 int numAtomsPerMolecule
= (serializer
->reading() ? 0 : mtop
->moltype
[molblock
.type
].atoms
.nr
);
2578 do_molblock(serializer
, &molblock
, numAtomsPerMolecule
);
2580 serializer
->doInt(&mtop
->natoms
);
2582 if (file_version
>= tpxv_IntermolecularBondeds
)
2584 serializer
->doBool(&mtop
->bIntermolecularInteractions
);
2585 if (mtop
->bIntermolecularInteractions
)
2587 if (serializer
->reading())
2589 mtop
->intermolecular_ilist
= std::make_unique
<InteractionLists
>();
2591 do_ilists(serializer
, mtop
->intermolecular_ilist
.get(), file_version
);
2596 mtop
->bIntermolecularInteractions
= FALSE
;
2599 do_atomtypes (serializer
, &(mtop
->atomtypes
), file_version
);
2601 if (file_version
>= 65)
2603 do_cmap(serializer
, &mtop
->ffparams
.cmap_grid
);
2607 mtop
->ffparams
.cmap_grid
.grid_spacing
= 0;
2608 mtop
->ffparams
.cmap_grid
.cmapdata
.clear();
2611 do_groups(serializer
, &mtop
->groups
, &(mtop
->symtab
));
2613 mtop
->haveMoleculeIndices
= true;
2615 if (serializer
->reading())
2617 close_symtab(&(mtop
->symtab
));
2622 * Read the first part of the TPR file to find general system information.
2624 * If \p TopOnlyOK is true then we can read even future versions
2625 * of tpx files, provided the \p fileGeneration hasn't changed.
2626 * If it is false, we need the \p ir too, and bail out
2627 * if the file is newer than the program.
2629 * The version and generation of the topology (see top of this file)
2630 * are returned in the two last arguments, if those arguments are non-nullptr.
2632 * If possible, we will read the \p ir even when \p TopOnlyOK is true.
2634 * \param[in,out] serializer The serializer used to handle header processing.
2635 * \param[in,out] tpx File header datastructure.
2636 * \param[in] filename The name of the file being read/written
2637 * \param[in,out] fio File handle.
2638 * \param[in] TopOnlyOK If not reading \p ir is fine or not.
2640 static void do_tpxheader(gmx::FileIOXdrSerializer
*serializer
,
2642 const char *filename
,
2651 /* XDR binary topology file */
2652 precision
= sizeof(real
);
2654 std::string fileTag
;
2655 if (serializer
->reading())
2657 serializer
->doString(&buf
);
2658 if (std::strncmp(buf
.c_str(), "VERSION", 7) != 0)
2660 gmx_fatal(FARGS
, "Can not read file %s,\n"
2661 " this file is from a GROMACS version which is older than 2.0\n"
2662 " Make a new one with grompp or use a gro or pdb file, if possible",
2665 serializer
->doInt(&precision
);
2666 bDouble
= (precision
== sizeof(double));
2667 if ((precision
!= sizeof(float)) && !bDouble
)
2669 gmx_fatal(FARGS
, "Unknown precision in file %s: real is %d bytes "
2670 "instead of %zu or %zu",
2671 filename
, precision
, sizeof(float), sizeof(double));
2673 gmx_fio_setprecision(fio
, bDouble
);
2674 fprintf(stderr
, "Reading file %s, %s (%s precision)\n",
2675 filename
, buf
.c_str(), bDouble
? "double" : "single");
2679 buf
= gmx::formatString("VERSION %s", gmx_version());
2680 serializer
->doString(&buf
);
2681 bDouble
= (precision
== sizeof(double));
2682 gmx_fio_setprecision(fio
, bDouble
);
2683 serializer
->doInt(&precision
);
2684 fileTag
= gmx::formatString("%s", tpx_tag
);
2687 /* Check versions! */
2688 serializer
->doInt(&tpx
->fileVersion
);
2690 /* This is for backward compatibility with development versions 77-79
2691 * where the tag was, mistakenly, placed before the generation,
2692 * which would cause a segv instead of a proper error message
2693 * when reading the topology only from tpx with <77 code.
2695 if (tpx
->fileVersion
>= 77 && tpx
->fileVersion
<= 79)
2697 serializer
->doString(&fileTag
);
2700 serializer
->doInt(&tpx
->fileGeneration
);
2702 if (tpx
->fileVersion
>= 81)
2704 serializer
->doString(&fileTag
);
2706 if (serializer
->reading())
2708 if (tpx
->fileVersion
< 77)
2710 /* Versions before 77 don't have the tag, set it to release */
2711 fileTag
= gmx::formatString("%s", TPX_TAG_RELEASE
);
2714 if (fileTag
!= tpx_tag
)
2716 fprintf(stderr
, "Note: file tpx tag '%s', software tpx tag '%s'\n",
2717 fileTag
.c_str(), tpx_tag
);
2719 /* We only support reading tpx files with the same tag as the code
2720 * or tpx files with the release tag and with lower version number.
2722 if (fileTag
!= TPX_TAG_RELEASE
&& tpx
->fileVersion
< tpx_version
)
2724 gmx_fatal(FARGS
, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2725 filename
, tpx
->fileVersion
, fileTag
.c_str(),
2726 tpx_version
, tpx_tag
);
2731 if ((tpx
->fileVersion
<= tpx_incompatible_version
) ||
2732 ((tpx
->fileVersion
> tpx_version
) && !TopOnlyOK
) ||
2733 (tpx
->fileGeneration
> tpx_generation
) ||
2734 tpx_version
== 80) /*80 was used by both 5.0-dev and 4.6-dev*/
2736 gmx_fatal(FARGS
, "reading tpx file (%s) version %d with version %d program",
2737 filename
, tpx
->fileVersion
, tpx_version
);
2740 serializer
->doInt(&tpx
->natoms
);
2741 serializer
->doInt(&tpx
->ngtc
);
2743 if (tpx
->fileVersion
< 62)
2745 serializer
->doInt(&idum
);
2746 serializer
->doReal(&rdum
);
2748 if (tpx
->fileVersion
>= 79)
2750 serializer
->doInt(&tpx
->fep_state
);
2752 serializer
->doReal(&tpx
->lambda
);
2753 serializer
->doBool(&tpx
->bIr
);
2754 serializer
->doBool(&tpx
->bTop
);
2755 serializer
->doBool(&tpx
->bX
);
2756 serializer
->doBool(&tpx
->bV
);
2757 serializer
->doBool(&tpx
->bF
);
2758 serializer
->doBool(&tpx
->bBox
);
2760 if ((tpx
->fileGeneration
> tpx_generation
))
2762 /* This can only happen if TopOnlyOK=TRUE */
2767 static int do_tpx_body(gmx::ISerializer
*serializer
,
2776 gmx_bool bPeriodicMols
;
2778 if (!serializer
->reading())
2780 GMX_RELEASE_ASSERT(x
== nullptr && v
== nullptr, "Passing separate x and v pointers to do_tpx() is not supported when writing");
2784 GMX_RELEASE_ASSERT(!(x
== nullptr && v
!= nullptr), "Passing x==NULL and v!=NULL is not supported");
2787 if (serializer
->reading())
2790 init_gtc_state(state
, tpx
->ngtc
, 0, 0);
2793 // v is also nullptr by the above assertion, so we may
2794 // need to make memory in state for storing the contents
2798 state
->flags
|= (1 << estX
);
2802 state
->flags
|= (1 << estV
);
2804 state_change_natoms(state
, tpx
->natoms
);
2810 x
= state
->x
.rvec_array();
2811 v
= state
->v
.rvec_array();
2814 #define do_test(serializer, b, p) if ((serializer)->reading() && ((p) != nullptr) && !(b)) gmx_fatal(FARGS, "No %s in input file",#p)
2816 do_test(serializer
, tpx
->bBox
, state
->box
);
2819 serializer
->doRvecArray(state
->box
, DIM
);
2820 if (tpx
->fileVersion
>= 51)
2822 serializer
->doRvecArray(state
->box_rel
, DIM
);
2826 /* We initialize box_rel after reading the inputrec */
2827 clear_mat(state
->box_rel
);
2829 serializer
->doRvecArray(state
->boxv
, DIM
);
2830 if (tpx
->fileVersion
< 56)
2833 serializer
->doRvecArray(mdum
, DIM
);
2837 if (state
->ngtc
> 0)
2840 snew(dumv
, state
->ngtc
);
2841 if (tpx
->fileVersion
< 69)
2843 serializer
->doRealArray(dumv
, state
->ngtc
);
2845 /* These used to be the Berendsen tcoupl_lambda's */
2846 serializer
->doRealArray(dumv
, state
->ngtc
);
2850 do_test(serializer
, tpx
->bTop
, mtop
);
2855 do_mtop(serializer
, mtop
, tpx
->fileVersion
);
2860 do_mtop(serializer
, &dum_top
, tpx
->fileVersion
);
2863 do_test(serializer
, tpx
->bX
, x
);
2866 if (serializer
->reading())
2868 state
->flags
|= (1<<estX
);
2870 serializer
->doRvecArray(x
, tpx
->natoms
);
2873 do_test(serializer
, tpx
->bV
, v
);
2876 if (serializer
->reading())
2878 state
->flags
|= (1<<estV
);
2882 std::vector
<gmx::RVec
> dummyVelocities(tpx
->natoms
);
2883 serializer
->doRvecArray(as_rvec_array(dummyVelocities
.data()), tpx
->natoms
);
2887 serializer
->doRvecArray(v
, tpx
->natoms
);
2891 // No need to run do_test when the last argument is NULL
2894 std::vector
<gmx::RVec
> dummyForces(state
->natoms
);
2895 serializer
->doRvecArray(as_rvec_array(dummyForces
.data()), tpx
->natoms
);
2898 /* Starting with tpx version 26, we have the inputrec
2899 * at the end of the file, so we can ignore it
2900 * if the file is never than the software (but still the
2901 * same generation - see comments at the top of this file.
2906 bPeriodicMols
= FALSE
;
2908 do_test(serializer
, tpx
->bIr
, ir
);
2911 if (tpx
->fileVersion
>= 53)
2913 /* Removed the pbc info from do_inputrec, since we always want it */
2914 if (!serializer
->reading())
2917 bPeriodicMols
= ir
->bPeriodicMols
;
2919 serializer
->doInt(&ePBC
);
2920 serializer
->doBool(&bPeriodicMols
);
2922 if (tpx
->fileGeneration
<= tpx_generation
&& ir
)
2924 do_inputrec(serializer
, ir
, tpx
->fileVersion
);
2925 if (tpx
->fileVersion
< 51)
2927 set_box_rel(ir
, state
);
2929 if (tpx
->fileVersion
< 53)
2932 bPeriodicMols
= ir
->bPeriodicMols
;
2935 if (serializer
->reading() && ir
&& tpx
->fileVersion
>= 53)
2937 /* We need to do this after do_inputrec, since that initializes ir */
2939 ir
->bPeriodicMols
= bPeriodicMols
;
2943 if (serializer
->reading())
2947 if (state
->ngtc
== 0)
2949 /* Reading old version without tcoupl state data: set it */
2950 init_gtc_state(state
, ir
->opts
.ngtc
, 0, ir
->opts
.nhchainlength
);
2952 if (tpx
->bTop
&& mtop
)
2954 if (tpx
->fileVersion
< 57)
2956 if (mtop
->moltype
[0].ilist
[F_DISRES
].size() > 0)
2958 ir
->eDisre
= edrSimple
;
2962 ir
->eDisre
= edrNone
;
2965 set_disres_npair(mtop
);
2969 if (tpx
->bTop
&& mtop
)
2971 gmx_mtop_finalize(mtop
);
2978 static t_fileio
*open_tpx(const char *fn
, const char *mode
)
2980 return gmx_fio_open(fn
, mode
);
2983 static void close_tpx(t_fileio
*fio
)
2989 * Fill information into the header only from state before writing.
2991 * Populating the header needs to be independent from writing the information
2992 * to file to allow things like writing the raw byte stream.
2994 * \param[in] state The current simulation state. Can't write without it.
2995 * \param[in] ir Parameter and system information.
2996 * \param[in] mtop Global topology.
2997 * \returns Fully populated header.
2999 static TpxFileHeader
populateTpxHeader(const t_state
&state
,
3000 const t_inputrec
*ir
,
3001 const gmx_mtop_t
*mtop
)
3003 TpxFileHeader header
;
3004 header
.natoms
= state
.natoms
;
3005 header
.ngtc
= state
.ngtc
;
3006 header
.fep_state
= state
.fep_state
;
3007 header
.lambda
= state
.lambda
[efptFEP
];
3008 header
.bIr
= ir
!= nullptr;
3009 header
.bTop
= mtop
!= nullptr;
3010 header
.bX
= (state
.flags
& (1 << estX
)) != 0;
3011 header
.bV
= (state
.flags
& (1 << estV
)) != 0;
3014 header
.fileVersion
= tpx_version
;
3015 header
.fileGeneration
= tpx_generation
;
3020 /************************************************************
3022 * The following routines are the exported ones
3024 ************************************************************/
3026 TpxFileHeader
readTpxHeader(const char *fileName
, bool canReadTopologyOnly
)
3030 fio
= open_tpx(fileName
, "r");
3031 gmx::FileIOXdrSerializer
serializer(fio
);
3034 do_tpxheader(&serializer
, &tpx
, fileName
, fio
, canReadTopologyOnly
);
3039 void write_tpx_state(const char *fn
,
3040 const t_inputrec
*ir
, const t_state
*state
, const gmx_mtop_t
*mtop
)
3044 fio
= open_tpx(fn
, "w");
3046 TpxFileHeader tpx
= populateTpxHeader(*state
,
3050 gmx::FileIOXdrSerializer
serializer(fio
);
3051 do_tpxheader(&serializer
,
3056 do_tpx_body(&serializer
,
3058 const_cast<t_inputrec
*>(ir
),
3059 const_cast<t_state
*>(state
), nullptr, nullptr,
3060 const_cast<gmx_mtop_t
*>(mtop
));
3064 void read_tpx_state(const char *fn
,
3065 t_inputrec
*ir
, t_state
*state
, gmx_mtop_t
*mtop
)
3070 fio
= open_tpx(fn
, "r");
3071 gmx::FileIOXdrSerializer
serializer(fio
);
3072 do_tpxheader(&serializer
,
3077 do_tpx_body(&serializer
,
3087 int read_tpx(const char *fn
,
3088 t_inputrec
*ir
, matrix box
, int *natoms
,
3089 rvec
*x
, rvec
*v
, gmx_mtop_t
*mtop
)
3096 fio
= open_tpx(fn
, "r");
3097 gmx::FileIOXdrSerializer
serializer(fio
);
3098 do_tpxheader(&serializer
,
3103 ePBC
= do_tpx_body(&serializer
,
3105 ir
, &state
, x
, v
, mtop
);
3107 if (mtop
!= nullptr && natoms
!= nullptr)
3109 *natoms
= mtop
->natoms
;
3113 copy_mat(state
.box
, box
);
3119 int read_tpx_top(const char *fn
,
3120 t_inputrec
*ir
, matrix box
, int *natoms
,
3121 rvec
*x
, rvec
*v
, t_topology
*top
)
3126 ePBC
= read_tpx(fn
, ir
, box
, natoms
, x
, v
, &mtop
);
3128 *top
= gmx_mtop_t_to_t_topology(&mtop
, true);
3133 gmx_bool
fn2bTPX(const char *file
)
3135 return (efTPR
== fn2ftp(file
));
3138 void pr_tpxheader(FILE *fp
, int indent
, const char *title
, const TpxFileHeader
*sh
)
3140 if (available(fp
, sh
, indent
, title
))
3142 indent
= pr_title(fp
, indent
, title
);
3143 pr_indent(fp
, indent
);
3144 fprintf(fp
, "bIr = %spresent\n", sh
->bIr
? "" : "not ");
3145 pr_indent(fp
, indent
);
3146 fprintf(fp
, "bBox = %spresent\n", sh
->bBox
? "" : "not ");
3147 pr_indent(fp
, indent
);
3148 fprintf(fp
, "bTop = %spresent\n", sh
->bTop
? "" : "not ");
3149 pr_indent(fp
, indent
);
3150 fprintf(fp
, "bX = %spresent\n", sh
->bX
? "" : "not ");
3151 pr_indent(fp
, indent
);
3152 fprintf(fp
, "bV = %spresent\n", sh
->bV
? "" : "not ");
3153 pr_indent(fp
, indent
);
3154 fprintf(fp
, "bF = %spresent\n", sh
->bF
? "" : "not ");
3156 pr_indent(fp
, indent
);
3157 fprintf(fp
, "natoms = %d\n", sh
->natoms
);
3158 pr_indent(fp
, indent
);
3159 fprintf(fp
, "lambda = %e\n", sh
->lambda
);