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/inmemoryserializer.h"
76 #include "gromacs/utility/keyvaluetreebuilder.h"
77 #include "gromacs/utility/keyvaluetreeserializer.h"
78 #include "gromacs/utility/smalloc.h"
79 #include "gromacs/utility/snprintf.h"
80 #include "gromacs/utility/txtdump.h"
82 #define TPX_TAG_RELEASE "release"
84 /*! \brief Tag string for the file format written to run input files
85 * written by this version of the code.
87 * Change this if you want to change the run input format in a feature
88 * branch. This ensures that there will not be different run input
89 * formats around which cannot be distinguished, while not causing
90 * problems rebasing the feature branch onto upstream changes. When
91 * merging with mainstream GROMACS, set this tag string back to
92 * TPX_TAG_RELEASE, and instead add an element to tpxv.
94 static const char* tpx_tag
= TPX_TAG_RELEASE
;
96 /*! \brief Enum of values that describe the contents of a tpr file
97 * whose format matches a version number
99 * The enum helps the code be more self-documenting and ensure merges
100 * do not silently resolve when two patches make the same bump. When
101 * adding new functionality, add a new element just above tpxv_Count
102 * in this enumeration, and write code below that does the right thing
103 * according to the value of file_version.
107 tpxv_ComputationalElectrophysiology
=
108 96, /**< support for ion/water position swaps (computational electrophysiology) */
109 tpxv_Use64BitRandomSeed
, /**< change ld_seed from int to int64_t */
110 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials
, /**< potentials for supporting coarse-grained force fields */
111 tpxv_InteractiveMolecularDynamics
, /**< interactive molecular dynamics (IMD) */
112 tpxv_RemoveObsoleteParameters1
, /**< remove optimize_fft, dihre_fc, nstcheckpoint */
113 tpxv_PullCoordTypeGeom
, /**< add pull type and geometry per group and flat-bottom */
114 tpxv_PullGeomDirRel
, /**< add pull geometry direction-relative */
115 tpxv_IntermolecularBondeds
, /**< permit inter-molecular bonded interactions in the topology */
116 tpxv_CompElWithSwapLayerOffset
, /**< added parameters for improved CompEl setups */
117 tpxv_CompElPolyatomicIonsAndMultipleIonTypes
, /**< CompEl now can handle polyatomic ions and more than two types of ions */
118 tpxv_RemoveAdress
, /**< removed support for AdResS */
119 tpxv_PullCoordNGroup
, /**< add ngroup to pull coord */
120 tpxv_RemoveTwinRange
, /**< removed support for twin-range interactions */
121 tpxv_ReplacePullPrintCOM12
, /**< Replaced print-com-1, 2 with pull-print-com */
122 tpxv_PullExternalPotential
, /**< Added pull type external potential */
123 tpxv_GenericParamsForElectricField
, /**< Introduced KeyValueTree and moved electric field parameters */
124 tpxv_AcceleratedWeightHistogram
, /**< sampling with accelerated weight histogram method (AWH) */
125 tpxv_RemoveImplicitSolvation
, /**< removed support for implicit solvation */
126 tpxv_PullPrevStepCOMAsReference
, /**< Enabled using the COM of the pull group of the last frame as reference for PBC */
127 tpxv_MimicQMMM
, /**< Introduced support for MiMiC QM/MM interface */
128 tpxv_PullAverage
, /**< Added possibility to output average pull force and position */
129 tpxv_GenericInternalParameters
, /**< Added internal parameters for mdrun modules*/
130 tpxv_VSite2FD
, /**< Added 2FD type virtual site */
131 tpxv_AddSizeField
, /**< Added field with information about the size of the serialized tpr file in bytes, excluding the header */
132 tpxv_Count
/**< the total number of tpxv versions */
135 /*! \brief Version number of the file format written to run input
136 * files by this version of the code.
138 * The tpx_version increases whenever the file format in the main
139 * development branch changes, due to an extension of the tpxv enum above.
140 * Backward compatibility for reading old run input files is maintained
141 * by checking this version number against that of the file and then using
142 * the correct code path.
144 * When developing a feature branch that needs to change the run input
145 * file format, change tpx_tag instead. */
146 static const int tpx_version
= tpxv_Count
- 1;
149 /* This number should only be increased when you edit the TOPOLOGY section
150 * or the HEADER of the tpx format.
151 * This way we can maintain forward compatibility too for all analysis tools
152 * and/or external programs that only need to know the atom/residue names,
153 * charges, and bond connectivity.
155 * It first appeared in tpx version 26, when I also moved the inputrecord
156 * to the end of the tpx file, so we can just skip it if we only
159 * In particular, it must be increased when adding new elements to
160 * ftupd, so that old code can read new .tpr files.
162 * Updated for added field that contains the number of bytes of the tpr body, excluding the header.
164 static const int tpx_generation
= 27;
166 /* This number should be the most recent backwards incompatible version
167 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
169 static const int tpx_incompatible_version
= 57; // GMX4.0 has version 58
172 /* Struct used to maintain tpx compatibility when function types are added */
175 int fvnr
; /* file version number in which the function type first appeared */
176 int ftype
; /* function type */
180 * TODO The following three lines make little sense, please clarify if
181 * you've had to work out how ftupd works.
183 * The entries should be ordered in:
184 * 1. ascending function type number
185 * 2. ascending file version number
187 * Because we support reading of old .tpr file versions (even when
188 * mdrun can no longer run the simulation), we need to be able to read
189 * obsolete t_interaction_function types. Any data read from such
190 * fields is discarded. Their names have _NOLONGERUSED appended to
191 * them to make things clear.
193 static const t_ftupd ftupd
[] = {
194 { 70, F_RESTRBONDS
},
195 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials
, F_RESTRANGLES
},
196 { 76, F_LINEAR_ANGLES
},
197 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials
, F_RESTRDIHS
},
198 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials
, F_CBTDIHS
},
200 { 60, F_GB12_NOLONGERUSED
},
201 { 61, F_GB13_NOLONGERUSED
},
202 { 61, F_GB14_NOLONGERUSED
},
203 { 72, F_GBPOL_NOLONGERUSED
},
204 { 72, F_NPSOLVATION_NOLONGERUSED
},
207 { 69, F_VTEMP_NOLONGERUSED
},
209 { 76, F_ANHARM_POL
},
219 { 79, F_DVDL_RESTRAINT
},
220 { 79, F_DVDL_TEMPERATURE
},
221 { tpxv_GenericInternalParameters
, F_DENSITYFITTING
},
222 { tpxv_VSite2FD
, F_VSITE2FD
},
224 #define NFTUPD asize(ftupd)
226 /* Needed for backward compatibility */
229 /**************************************************************
231 * Now the higer level routines that do io of the structures and arrays
233 **************************************************************/
234 static void do_pullgrp_tpx_pre95(gmx::ISerializer
* serializer
, t_pull_group
* pgrp
, t_pull_coord
* pcrd
)
238 serializer
->doInt(&pgrp
->nat
);
239 if (serializer
->reading())
241 snew(pgrp
->ind
, pgrp
->nat
);
243 serializer
->doIntArray(pgrp
->ind
, pgrp
->nat
);
244 serializer
->doInt(&pgrp
->nweight
);
245 if (serializer
->reading())
247 snew(pgrp
->weight
, pgrp
->nweight
);
249 serializer
->doRealArray(pgrp
->weight
, pgrp
->nweight
);
250 serializer
->doInt(&pgrp
->pbcatom
);
251 serializer
->doRvec(&pcrd
->vec
);
252 clear_rvec(pcrd
->origin
);
253 serializer
->doRvec(&tmp
);
255 serializer
->doReal(&pcrd
->rate
);
256 serializer
->doReal(&pcrd
->k
);
257 serializer
->doReal(&pcrd
->kB
);
260 static void do_pull_group(gmx::ISerializer
* serializer
, t_pull_group
* pgrp
)
262 serializer
->doInt(&pgrp
->nat
);
263 if (serializer
->reading())
265 snew(pgrp
->ind
, pgrp
->nat
);
267 serializer
->doIntArray(pgrp
->ind
, pgrp
->nat
);
268 serializer
->doInt(&pgrp
->nweight
);
269 if (serializer
->reading())
271 snew(pgrp
->weight
, pgrp
->nweight
);
273 serializer
->doRealArray(pgrp
->weight
, pgrp
->nweight
);
274 serializer
->doInt(&pgrp
->pbcatom
);
277 static void do_pull_coord(gmx::ISerializer
* serializer
,
284 if (file_version
>= tpxv_PullCoordNGroup
)
286 serializer
->doInt(&pcrd
->eType
);
287 if (file_version
>= tpxv_PullExternalPotential
)
289 if (pcrd
->eType
== epullEXTERNAL
)
292 if (serializer
->reading())
294 serializer
->doString(&buf
);
295 pcrd
->externalPotentialProvider
= gmx_strdup(buf
.c_str());
299 buf
= pcrd
->externalPotentialProvider
;
300 serializer
->doString(&buf
);
305 pcrd
->externalPotentialProvider
= nullptr;
310 if (serializer
->reading())
312 pcrd
->externalPotentialProvider
= nullptr;
315 /* Note that we try to support adding new geometries without
316 * changing the tpx version. This requires checks when printing the
317 * geometry string and a check and fatal_error in init_pull.
319 serializer
->doInt(&pcrd
->eGeom
);
320 serializer
->doInt(&pcrd
->ngroup
);
321 if (pcrd
->ngroup
<= c_pullCoordNgroupMax
)
323 serializer
->doIntArray(pcrd
->group
, pcrd
->ngroup
);
327 /* More groups in file than supported, this must be a new geometry
328 * that is not supported by our current code. Since we will not
329 * use the groups for this coord (checks in the pull and WHAM code
330 * ensure this), we can ignore the groups and set ngroup=0.
333 snew(dum
, pcrd
->ngroup
);
334 serializer
->doIntArray(dum
, pcrd
->ngroup
);
339 serializer
->doIvec(&pcrd
->dim
);
344 serializer
->doInt(&pcrd
->group
[0]);
345 serializer
->doInt(&pcrd
->group
[1]);
346 if (file_version
>= tpxv_PullCoordTypeGeom
)
348 pcrd
->ngroup
= (pcrd
->eGeom
== epullgDIRRELATIVE
? 4 : 2);
349 serializer
->doInt(&pcrd
->eType
);
350 serializer
->doInt(&pcrd
->eGeom
);
351 if (pcrd
->ngroup
== 4)
353 serializer
->doInt(&pcrd
->group
[2]);
354 serializer
->doInt(&pcrd
->group
[3]);
356 serializer
->doIvec(&pcrd
->dim
);
360 pcrd
->eType
= ePullOld
;
361 pcrd
->eGeom
= eGeomOld
;
362 copy_ivec(dimOld
, pcrd
->dim
);
365 serializer
->doRvec(&pcrd
->origin
);
366 serializer
->doRvec(&pcrd
->vec
);
367 if (file_version
>= tpxv_PullCoordTypeGeom
)
369 serializer
->doBool(&pcrd
->bStart
);
373 /* This parameter is only printed, but not actually used by mdrun */
374 pcrd
->bStart
= FALSE
;
376 serializer
->doReal(&pcrd
->init
);
377 serializer
->doReal(&pcrd
->rate
);
378 serializer
->doReal(&pcrd
->k
);
379 serializer
->doReal(&pcrd
->kB
);
382 static void do_expandedvals(gmx::ISerializer
* serializer
, t_expanded
* expand
, t_lambda
* fepvals
, int file_version
)
384 int n_lambda
= fepvals
->n_lambda
;
386 /* reset the lambda calculation window */
387 fepvals
->lambda_start_n
= 0;
388 fepvals
->lambda_stop_n
= n_lambda
;
389 if (file_version
>= 79)
393 if (serializer
->reading())
395 snew(expand
->init_lambda_weights
, n_lambda
);
397 serializer
->doRealArray(expand
->init_lambda_weights
, n_lambda
);
398 serializer
->doBool(&expand
->bInit_weights
);
401 serializer
->doInt(&expand
->nstexpanded
);
402 serializer
->doInt(&expand
->elmcmove
);
403 serializer
->doInt(&expand
->elamstats
);
404 serializer
->doInt(&expand
->lmc_repeats
);
405 serializer
->doInt(&expand
->gibbsdeltalam
);
406 serializer
->doInt(&expand
->lmc_forced_nstart
);
407 serializer
->doInt(&expand
->lmc_seed
);
408 serializer
->doReal(&expand
->mc_temp
);
409 serializer
->doBool(&expand
->bSymmetrizedTMatrix
);
410 serializer
->doInt(&expand
->nstTij
);
411 serializer
->doInt(&expand
->minvarmin
);
412 serializer
->doInt(&expand
->c_range
);
413 serializer
->doReal(&expand
->wl_scale
);
414 serializer
->doReal(&expand
->wl_ratio
);
415 serializer
->doReal(&expand
->init_wl_delta
);
416 serializer
->doBool(&expand
->bWLoneovert
);
417 serializer
->doInt(&expand
->elmceq
);
418 serializer
->doInt(&expand
->equil_steps
);
419 serializer
->doInt(&expand
->equil_samples
);
420 serializer
->doInt(&expand
->equil_n_at_lam
);
421 serializer
->doReal(&expand
->equil_wl_delta
);
422 serializer
->doReal(&expand
->equil_ratio
);
426 static void do_simtempvals(gmx::ISerializer
* serializer
, t_simtemp
* simtemp
, int n_lambda
, int file_version
)
428 if (file_version
>= 79)
430 serializer
->doInt(&simtemp
->eSimTempScale
);
431 serializer
->doReal(&simtemp
->simtemp_high
);
432 serializer
->doReal(&simtemp
->simtemp_low
);
435 if (serializer
->reading())
437 snew(simtemp
->temperatures
, n_lambda
);
439 serializer
->doRealArray(simtemp
->temperatures
, n_lambda
);
444 static void do_imd(gmx::ISerializer
* serializer
, t_IMD
* imd
)
446 serializer
->doInt(&imd
->nat
);
447 if (serializer
->reading())
449 snew(imd
->ind
, imd
->nat
);
451 serializer
->doIntArray(imd
->ind
, imd
->nat
);
454 static void do_fepvals(gmx::ISerializer
* serializer
, t_lambda
* fepvals
, int file_version
)
456 /* i is defined in the ndo_double macro; use g to iterate. */
460 /* free energy values */
462 if (file_version
>= 79)
464 serializer
->doInt(&fepvals
->init_fep_state
);
465 serializer
->doDouble(&fepvals
->init_lambda
);
466 serializer
->doDouble(&fepvals
->delta_lambda
);
468 else if (file_version
>= 59)
470 serializer
->doDouble(&fepvals
->init_lambda
);
471 serializer
->doDouble(&fepvals
->delta_lambda
);
475 serializer
->doReal(&rdum
);
476 fepvals
->init_lambda
= rdum
;
477 serializer
->doReal(&rdum
);
478 fepvals
->delta_lambda
= rdum
;
480 if (file_version
>= 79)
482 serializer
->doInt(&fepvals
->n_lambda
);
483 if (serializer
->reading())
485 snew(fepvals
->all_lambda
, efptNR
);
487 for (g
= 0; g
< efptNR
; g
++)
489 if (fepvals
->n_lambda
> 0)
491 if (serializer
->reading())
493 snew(fepvals
->all_lambda
[g
], fepvals
->n_lambda
);
495 serializer
->doDoubleArray(fepvals
->all_lambda
[g
], fepvals
->n_lambda
);
496 serializer
->doBoolArray(fepvals
->separate_dvdl
, efptNR
);
498 else if (fepvals
->init_lambda
>= 0)
500 fepvals
->separate_dvdl
[efptFEP
] = TRUE
;
504 else if (file_version
>= 64)
506 serializer
->doInt(&fepvals
->n_lambda
);
507 if (serializer
->reading())
511 snew(fepvals
->all_lambda
, efptNR
);
512 /* still allocate the all_lambda array's contents. */
513 for (g
= 0; g
< efptNR
; g
++)
515 if (fepvals
->n_lambda
> 0)
517 snew(fepvals
->all_lambda
[g
], fepvals
->n_lambda
);
521 serializer
->doDoubleArray(fepvals
->all_lambda
[efptFEP
], fepvals
->n_lambda
);
522 if (fepvals
->init_lambda
>= 0)
526 fepvals
->separate_dvdl
[efptFEP
] = TRUE
;
528 if (serializer
->reading())
530 /* copy the contents of the efptFEP lambda component to all
531 the other components */
532 for (g
= 0; g
< efptNR
; g
++)
534 for (h
= 0; h
< fepvals
->n_lambda
; h
++)
538 fepvals
->all_lambda
[g
][h
] = fepvals
->all_lambda
[efptFEP
][h
];
547 fepvals
->n_lambda
= 0;
548 fepvals
->all_lambda
= nullptr;
549 if (fepvals
->init_lambda
>= 0)
551 fepvals
->separate_dvdl
[efptFEP
] = TRUE
;
554 serializer
->doReal(&fepvals
->sc_alpha
);
555 serializer
->doInt(&fepvals
->sc_power
);
556 if (file_version
>= 79)
558 serializer
->doReal(&fepvals
->sc_r_power
);
562 fepvals
->sc_r_power
= 6.0;
564 serializer
->doReal(&fepvals
->sc_sigma
);
565 if (serializer
->reading())
567 if (file_version
>= 71)
569 fepvals
->sc_sigma_min
= fepvals
->sc_sigma
;
573 fepvals
->sc_sigma_min
= 0;
576 if (file_version
>= 79)
578 serializer
->doBool(&fepvals
->bScCoul
);
582 fepvals
->bScCoul
= TRUE
;
584 if (file_version
>= 64)
586 serializer
->doInt(&fepvals
->nstdhdl
);
590 fepvals
->nstdhdl
= 1;
593 if (file_version
>= 73)
595 serializer
->doInt(&fepvals
->separate_dhdl_file
);
596 serializer
->doInt(&fepvals
->dhdl_derivatives
);
600 fepvals
->separate_dhdl_file
= esepdhdlfileYES
;
601 fepvals
->dhdl_derivatives
= edhdlderivativesYES
;
603 if (file_version
>= 71)
605 serializer
->doInt(&fepvals
->dh_hist_size
);
606 serializer
->doDouble(&fepvals
->dh_hist_spacing
);
610 fepvals
->dh_hist_size
= 0;
611 fepvals
->dh_hist_spacing
= 0.1;
613 if (file_version
>= 79)
615 serializer
->doInt(&fepvals
->edHdLPrintEnergy
);
619 fepvals
->edHdLPrintEnergy
= edHdLPrintEnergyNO
;
622 /* handle lambda_neighbors */
623 if ((file_version
>= 83 && file_version
< 90) || file_version
>= 92)
625 serializer
->doInt(&fepvals
->lambda_neighbors
);
626 if ((fepvals
->lambda_neighbors
>= 0) && (fepvals
->init_fep_state
>= 0)
627 && (fepvals
->init_lambda
< 0))
629 fepvals
->lambda_start_n
= (fepvals
->init_fep_state
- fepvals
->lambda_neighbors
);
630 fepvals
->lambda_stop_n
= (fepvals
->init_fep_state
+ fepvals
->lambda_neighbors
+ 1);
631 if (fepvals
->lambda_start_n
< 0)
633 fepvals
->lambda_start_n
= 0;
635 if (fepvals
->lambda_stop_n
>= fepvals
->n_lambda
)
637 fepvals
->lambda_stop_n
= fepvals
->n_lambda
;
642 fepvals
->lambda_start_n
= 0;
643 fepvals
->lambda_stop_n
= fepvals
->n_lambda
;
648 fepvals
->lambda_start_n
= 0;
649 fepvals
->lambda_stop_n
= fepvals
->n_lambda
;
653 static void do_awhBias(gmx::ISerializer
* serializer
, gmx::AwhBiasParams
* awhBiasParams
, int gmx_unused file_version
)
655 serializer
->doInt(&awhBiasParams
->eTarget
);
656 serializer
->doDouble(&awhBiasParams
->targetBetaScaling
);
657 serializer
->doDouble(&awhBiasParams
->targetCutoff
);
658 serializer
->doInt(&awhBiasParams
->eGrowth
);
659 serializer
->doInt(&awhBiasParams
->bUserData
);
660 serializer
->doDouble(&awhBiasParams
->errorInitial
);
661 serializer
->doInt(&awhBiasParams
->ndim
);
662 serializer
->doInt(&awhBiasParams
->shareGroup
);
663 serializer
->doBool(&awhBiasParams
->equilibrateHistogram
);
665 if (serializer
->reading())
667 snew(awhBiasParams
->dimParams
, awhBiasParams
->ndim
);
670 for (int d
= 0; d
< awhBiasParams
->ndim
; d
++)
672 gmx::AwhDimParams
* dimParams
= &awhBiasParams
->dimParams
[d
];
674 serializer
->doInt(&dimParams
->eCoordProvider
);
675 serializer
->doInt(&dimParams
->coordIndex
);
676 serializer
->doDouble(&dimParams
->origin
);
677 serializer
->doDouble(&dimParams
->end
);
678 serializer
->doDouble(&dimParams
->period
);
679 serializer
->doDouble(&dimParams
->forceConstant
);
680 serializer
->doDouble(&dimParams
->diffusion
);
681 serializer
->doDouble(&dimParams
->coordValueInit
);
682 serializer
->doDouble(&dimParams
->coverDiameter
);
686 static void do_awh(gmx::ISerializer
* serializer
, gmx::AwhParams
* awhParams
, 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
, pull_params_t
* pull
, int file_version
, int ePullOld
)
716 if (file_version
>= 95)
718 serializer
->doInt(&pull
->ngroup
);
720 serializer
->doInt(&pull
->ncoord
);
721 if (file_version
< 95)
723 pull
->ngroup
= pull
->ncoord
+ 1;
725 if (file_version
< tpxv_PullCoordTypeGeom
)
729 serializer
->doInt(&eGeomOld
);
730 serializer
->doIvec(&dimOld
);
731 /* The inner cylinder radius, now removed */
732 serializer
->doReal(&dum
);
734 serializer
->doReal(&pull
->cylinder_r
);
735 serializer
->doReal(&pull
->constr_tol
);
736 if (file_version
>= 95)
738 serializer
->doBool(&pull
->bPrintCOM
);
739 /* With file_version < 95 this value is set below */
741 if (file_version
>= tpxv_ReplacePullPrintCOM12
)
743 serializer
->doBool(&pull
->bPrintRefValue
);
744 serializer
->doBool(&pull
->bPrintComp
);
746 else if (file_version
>= tpxv_PullCoordTypeGeom
)
749 serializer
->doInt(&idum
); /* used to be bPrintCOM2 */
750 serializer
->doBool(&pull
->bPrintRefValue
);
751 serializer
->doBool(&pull
->bPrintComp
);
755 pull
->bPrintRefValue
= FALSE
;
756 pull
->bPrintComp
= TRUE
;
758 serializer
->doInt(&pull
->nstxout
);
759 serializer
->doInt(&pull
->nstfout
);
760 if (file_version
>= tpxv_PullPrevStepCOMAsReference
)
762 serializer
->doBool(&pull
->bSetPbcRefToPrevStepCOM
);
766 pull
->bSetPbcRefToPrevStepCOM
= FALSE
;
768 if (serializer
->reading())
770 snew(pull
->group
, pull
->ngroup
);
771 snew(pull
->coord
, pull
->ncoord
);
773 if (file_version
< 95)
775 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
776 if (eGeomOld
== epullgDIRPBC
)
778 gmx_fatal(FARGS
, "pull-geometry=position is no longer supported");
780 if (eGeomOld
> epullgDIRPBC
)
785 for (g
= 0; g
< pull
->ngroup
; g
++)
787 /* We read and ignore a pull coordinate for group 0 */
788 do_pullgrp_tpx_pre95(serializer
, &pull
->group
[g
], &pull
->coord
[std::max(g
- 1, 0)]);
791 pull
->coord
[g
- 1].group
[0] = 0;
792 pull
->coord
[g
- 1].group
[1] = g
;
796 pull
->bPrintCOM
= (pull
->group
[0].nat
> 0);
800 for (g
= 0; g
< pull
->ngroup
; g
++)
802 do_pull_group(serializer
, &pull
->group
[g
]);
804 for (g
= 0; g
< pull
->ncoord
; g
++)
806 do_pull_coord(serializer
, &pull
->coord
[g
], file_version
, ePullOld
, eGeomOld
, dimOld
);
809 if (file_version
>= tpxv_PullAverage
)
813 v
= pull
->bXOutAverage
;
814 serializer
->doBool(&v
);
815 pull
->bXOutAverage
= v
;
816 v
= pull
->bFOutAverage
;
817 serializer
->doBool(&v
);
818 pull
->bFOutAverage
= v
;
823 static void do_rotgrp(gmx::ISerializer
* serializer
, t_rotgrp
* rotg
)
825 serializer
->doInt(&rotg
->eType
);
826 serializer
->doInt(&rotg
->bMassW
);
827 serializer
->doInt(&rotg
->nat
);
828 if (serializer
->reading())
830 snew(rotg
->ind
, rotg
->nat
);
832 serializer
->doIntArray(rotg
->ind
, rotg
->nat
);
833 if (serializer
->reading())
835 snew(rotg
->x_ref
, rotg
->nat
);
837 serializer
->doRvecArray(rotg
->x_ref
, rotg
->nat
);
838 serializer
->doRvec(&rotg
->inputVec
);
839 serializer
->doRvec(&rotg
->pivot
);
840 serializer
->doReal(&rotg
->rate
);
841 serializer
->doReal(&rotg
->k
);
842 serializer
->doReal(&rotg
->slab_dist
);
843 serializer
->doReal(&rotg
->min_gaussian
);
844 serializer
->doReal(&rotg
->eps
);
845 serializer
->doInt(&rotg
->eFittype
);
846 serializer
->doInt(&rotg
->PotAngle_nstep
);
847 serializer
->doReal(&rotg
->PotAngle_step
);
850 static void do_rot(gmx::ISerializer
* serializer
, t_rot
* rot
)
854 serializer
->doInt(&rot
->ngrp
);
855 serializer
->doInt(&rot
->nstrout
);
856 serializer
->doInt(&rot
->nstsout
);
857 if (serializer
->reading())
859 snew(rot
->grp
, rot
->ngrp
);
861 for (g
= 0; g
< rot
->ngrp
; g
++)
863 do_rotgrp(serializer
, &rot
->grp
[g
]);
868 static void do_swapgroup(gmx::ISerializer
* serializer
, t_swapGroup
* g
)
871 /* Name of the group or molecule */
873 if (serializer
->reading())
875 serializer
->doString(&buf
);
876 g
->molname
= gmx_strdup(buf
.c_str());
881 serializer
->doString(&buf
);
884 /* Number of atoms in the group */
885 serializer
->doInt(&g
->nat
);
887 /* The group's atom indices */
888 if (serializer
->reading())
890 snew(g
->ind
, g
->nat
);
892 serializer
->doIntArray(g
->ind
, g
->nat
);
894 /* Requested counts for compartments A and B */
895 serializer
->doIntArray(g
->nmolReq
, eCompNR
);
898 static void do_swapcoords_tpx(gmx::ISerializer
* serializer
, t_swapcoords
* swap
, int file_version
)
900 /* Enums for better readability of the code */
913 if (file_version
>= tpxv_CompElPolyatomicIonsAndMultipleIonTypes
)
915 /* The total number of swap groups is the sum of the fixed groups
916 * (split0, split1, solvent), and the user-defined groups (2+ types of ions)
918 serializer
->doInt(&swap
->ngrp
);
919 if (serializer
->reading())
921 snew(swap
->grp
, swap
->ngrp
);
923 for (int ig
= 0; ig
< swap
->ngrp
; ig
++)
925 do_swapgroup(serializer
, &swap
->grp
[ig
]);
927 serializer
->doBool(&swap
->massw_split
[eChannel0
]);
928 serializer
->doBool(&swap
->massw_split
[eChannel1
]);
929 serializer
->doInt(&swap
->nstswap
);
930 serializer
->doInt(&swap
->nAverage
);
931 serializer
->doReal(&swap
->threshold
);
932 serializer
->doReal(&swap
->cyl0r
);
933 serializer
->doReal(&swap
->cyl0u
);
934 serializer
->doReal(&swap
->cyl0l
);
935 serializer
->doReal(&swap
->cyl1r
);
936 serializer
->doReal(&swap
->cyl1u
);
937 serializer
->doReal(&swap
->cyl1l
);
941 /*** Support reading older CompEl .tpr files ***/
943 /* In the original CompEl .tpr files, we always have 5 groups: */
945 snew(swap
->grp
, swap
->ngrp
);
947 swap
->grp
[eGrpSplit0
].molname
= gmx_strdup("split0"); // group 0: split0
948 swap
->grp
[eGrpSplit1
].molname
= gmx_strdup("split1"); // group 1: split1
949 swap
->grp
[eGrpSolvent
].molname
= gmx_strdup("solvent"); // group 2: solvent
950 swap
->grp
[3].molname
= gmx_strdup("anions"); // group 3: anions
951 swap
->grp
[4].molname
= gmx_strdup("cations"); // group 4: cations
953 serializer
->doInt(&swap
->grp
[3].nat
);
954 serializer
->doInt(&swap
->grp
[eGrpSolvent
].nat
);
955 serializer
->doInt(&swap
->grp
[eGrpSplit0
].nat
);
956 serializer
->doBool(&swap
->massw_split
[eChannel0
]);
957 serializer
->doInt(&swap
->grp
[eGrpSplit1
].nat
);
958 serializer
->doBool(&swap
->massw_split
[eChannel1
]);
959 serializer
->doInt(&swap
->nstswap
);
960 serializer
->doInt(&swap
->nAverage
);
961 serializer
->doReal(&swap
->threshold
);
962 serializer
->doReal(&swap
->cyl0r
);
963 serializer
->doReal(&swap
->cyl0u
);
964 serializer
->doReal(&swap
->cyl0l
);
965 serializer
->doReal(&swap
->cyl1r
);
966 serializer
->doReal(&swap
->cyl1u
);
967 serializer
->doReal(&swap
->cyl1l
);
969 // The order[] array keeps compatibility with older .tpr files
970 // by reading in the groups in the classic order
972 const int order
[4] = { 3, eGrpSolvent
, eGrpSplit0
, eGrpSplit1
};
974 for (int ig
= 0; ig
< 4; ig
++)
977 snew(swap
->grp
[g
].ind
, swap
->grp
[g
].nat
);
978 serializer
->doIntArray(swap
->grp
[g
].ind
, swap
->grp
[g
].nat
);
982 for (int j
= eCompA
; j
<= eCompB
; j
++)
984 serializer
->doInt(&swap
->grp
[3].nmolReq
[j
]); // group 3 = anions
985 serializer
->doInt(&swap
->grp
[4].nmolReq
[j
]); // group 4 = cations
987 } /* End support reading older CompEl .tpr files */
989 if (file_version
>= tpxv_CompElWithSwapLayerOffset
)
991 serializer
->doReal(&swap
->bulkOffset
[eCompA
]);
992 serializer
->doReal(&swap
->bulkOffset
[eCompB
]);
996 static void do_legacy_efield(gmx::ISerializer
* serializer
, gmx::KeyValueTreeObjectBuilder
* root
)
998 const char* const dimName
[] = { "x", "y", "z" };
1000 auto appliedForcesObj
= root
->addObject("applied-forces");
1001 auto efieldObj
= appliedForcesObj
.addObject("electric-field");
1002 // The content of the tpr file for this feature has
1003 // been the same since gromacs 4.0 that was used for
1005 for (int j
= 0; j
< DIM
; ++j
)
1008 serializer
->doInt(&n
);
1009 serializer
->doInt(&nt
);
1010 std::vector
<real
> aa(n
+ 1), phi(nt
+ 1), at(nt
+ 1), phit(nt
+ 1);
1011 serializer
->doRealArray(aa
.data(), n
);
1012 serializer
->doRealArray(phi
.data(), n
);
1013 serializer
->doRealArray(at
.data(), nt
);
1014 serializer
->doRealArray(phit
.data(), nt
);
1017 if (n
> 1 || nt
> 1)
1020 "Can not handle tpr files with more than one electric field term per "
1023 auto dimObj
= efieldObj
.addObject(dimName
[j
]);
1024 dimObj
.addValue
<real
>("E0", aa
[0]);
1025 dimObj
.addValue
<real
>("omega", at
[0]);
1026 dimObj
.addValue
<real
>("t0", phi
[0]);
1027 dimObj
.addValue
<real
>("sigma", phit
[0]);
1033 static void do_inputrec(gmx::ISerializer
* serializer
, t_inputrec
* ir
, int file_version
)
1035 int i
, j
, k
, idum
= 0;
1037 gmx_bool bdum
= false;
1039 if (file_version
!= tpx_version
)
1041 /* Give a warning about features that are not accessible */
1042 fprintf(stderr
, "Note: file tpx version %d, software tpx version %d\n", file_version
, tpx_version
);
1045 if (file_version
== 0)
1050 gmx::KeyValueTreeBuilder paramsBuilder
;
1051 gmx::KeyValueTreeObjectBuilder paramsObj
= paramsBuilder
.rootObject();
1053 /* Basic inputrec stuff */
1054 serializer
->doInt(&ir
->eI
);
1055 if (file_version
>= 62)
1057 serializer
->doInt64(&ir
->nsteps
);
1061 serializer
->doInt(&idum
);
1065 if (file_version
>= 62)
1067 serializer
->doInt64(&ir
->init_step
);
1071 serializer
->doInt(&idum
);
1072 ir
->init_step
= idum
;
1075 serializer
->doInt(&ir
->simulation_part
);
1077 if (file_version
>= 67)
1079 serializer
->doInt(&ir
->nstcalcenergy
);
1083 ir
->nstcalcenergy
= 1;
1085 if (file_version
>= 81)
1087 serializer
->doInt(&ir
->cutoff_scheme
);
1088 if (file_version
< 94)
1090 ir
->cutoff_scheme
= 1 - ir
->cutoff_scheme
;
1095 ir
->cutoff_scheme
= ecutsGROUP
;
1097 serializer
->doInt(&idum
); /* used to be ns_type; not used anymore */
1098 serializer
->doInt(&ir
->nstlist
);
1099 serializer
->doInt(&idum
); /* used to be ndelta; not used anymore */
1101 serializer
->doReal(&ir
->rtpi
);
1103 serializer
->doInt(&ir
->nstcomm
);
1104 serializer
->doInt(&ir
->comm_mode
);
1106 /* ignore nstcheckpoint */
1107 if (file_version
< tpxv_RemoveObsoleteParameters1
)
1109 serializer
->doInt(&idum
);
1112 serializer
->doInt(&ir
->nstcgsteep
);
1114 serializer
->doInt(&ir
->nbfgscorr
);
1116 serializer
->doInt(&ir
->nstlog
);
1117 serializer
->doInt(&ir
->nstxout
);
1118 serializer
->doInt(&ir
->nstvout
);
1119 serializer
->doInt(&ir
->nstfout
);
1120 serializer
->doInt(&ir
->nstenergy
);
1121 serializer
->doInt(&ir
->nstxout_compressed
);
1122 if (file_version
>= 59)
1124 serializer
->doDouble(&ir
->init_t
);
1125 serializer
->doDouble(&ir
->delta_t
);
1129 serializer
->doReal(&rdum
);
1131 serializer
->doReal(&rdum
);
1134 serializer
->doReal(&ir
->x_compression_precision
);
1135 if (file_version
>= 81)
1137 serializer
->doReal(&ir
->verletbuf_tol
);
1141 ir
->verletbuf_tol
= 0;
1143 serializer
->doReal(&ir
->rlist
);
1144 if (file_version
>= 67 && file_version
< tpxv_RemoveTwinRange
)
1146 if (serializer
->reading())
1148 // Reading such a file version could invoke the twin-range
1149 // scheme, about which mdrun should give a fatal error.
1150 real dummy_rlistlong
= -1;
1151 serializer
->doReal(&dummy_rlistlong
);
1153 if (ir
->rlist
> 0 && (dummy_rlistlong
== 0 || dummy_rlistlong
> ir
->rlist
))
1155 // Get mdrun to issue an error (regardless of
1156 // ir->cutoff_scheme).
1157 ir
->useTwinRange
= true;
1161 // grompp used to set rlistlong actively. Users were
1162 // probably also confused and set rlistlong == rlist.
1163 // However, in all remaining cases, it is safe to let
1164 // mdrun proceed normally.
1165 ir
->useTwinRange
= false;
1171 // No need to read or write anything
1172 ir
->useTwinRange
= false;
1174 if (file_version
>= 82 && file_version
!= 90)
1176 // Multiple time-stepping is no longer enabled, but the old
1177 // support required the twin-range scheme, for which mdrun
1178 // already emits a fatal error.
1179 int dummy_nstcalclr
= -1;
1180 serializer
->doInt(&dummy_nstcalclr
);
1182 serializer
->doInt(&ir
->coulombtype
);
1183 if (file_version
>= 81)
1185 serializer
->doInt(&ir
->coulomb_modifier
);
1189 ir
->coulomb_modifier
= (ir
->cutoff_scheme
== ecutsVERLET
? eintmodPOTSHIFT
: eintmodNONE
);
1191 serializer
->doReal(&ir
->rcoulomb_switch
);
1192 serializer
->doReal(&ir
->rcoulomb
);
1193 serializer
->doInt(&ir
->vdwtype
);
1194 if (file_version
>= 81)
1196 serializer
->doInt(&ir
->vdw_modifier
);
1200 ir
->vdw_modifier
= (ir
->cutoff_scheme
== ecutsVERLET
? eintmodPOTSHIFT
: eintmodNONE
);
1202 serializer
->doReal(&ir
->rvdw_switch
);
1203 serializer
->doReal(&ir
->rvdw
);
1204 serializer
->doInt(&ir
->eDispCorr
);
1205 serializer
->doReal(&ir
->epsilon_r
);
1206 serializer
->doReal(&ir
->epsilon_rf
);
1207 serializer
->doReal(&ir
->tabext
);
1209 // This permits reading a .tpr file that used implicit solvent,
1210 // and later permitting mdrun to refuse to run it.
1211 if (serializer
->reading())
1213 if (file_version
< tpxv_RemoveImplicitSolvation
)
1215 serializer
->doInt(&idum
);
1216 serializer
->doInt(&idum
);
1217 serializer
->doReal(&rdum
);
1218 serializer
->doReal(&rdum
);
1219 serializer
->doInt(&idum
);
1220 ir
->implicit_solvent
= (idum
> 0);
1224 ir
->implicit_solvent
= false;
1226 if (file_version
< tpxv_RemoveImplicitSolvation
)
1228 serializer
->doReal(&rdum
);
1229 serializer
->doReal(&rdum
);
1230 serializer
->doReal(&rdum
);
1231 serializer
->doReal(&rdum
);
1232 if (file_version
>= 60)
1234 serializer
->doReal(&rdum
);
1235 serializer
->doInt(&idum
);
1237 serializer
->doReal(&rdum
);
1241 if (file_version
>= 81)
1243 serializer
->doReal(&ir
->fourier_spacing
);
1247 ir
->fourier_spacing
= 0.0;
1249 serializer
->doInt(&ir
->nkx
);
1250 serializer
->doInt(&ir
->nky
);
1251 serializer
->doInt(&ir
->nkz
);
1252 serializer
->doInt(&ir
->pme_order
);
1253 serializer
->doReal(&ir
->ewald_rtol
);
1255 if (file_version
>= 93)
1257 serializer
->doReal(&ir
->ewald_rtol_lj
);
1261 ir
->ewald_rtol_lj
= ir
->ewald_rtol
;
1263 serializer
->doInt(&ir
->ewald_geometry
);
1264 serializer
->doReal(&ir
->epsilon_surface
);
1266 /* ignore bOptFFT */
1267 if (file_version
< tpxv_RemoveObsoleteParameters1
)
1269 serializer
->doBool(&bdum
);
1272 if (file_version
>= 93)
1274 serializer
->doInt(&ir
->ljpme_combination_rule
);
1276 serializer
->doBool(&ir
->bContinuation
);
1277 serializer
->doInt(&ir
->etc
);
1278 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1279 * but the values 0 and 1 still mean no and
1280 * berendsen temperature coupling, respectively.
1282 if (file_version
>= 79)
1284 serializer
->doBool(&ir
->bPrintNHChains
);
1286 if (file_version
>= 71)
1288 serializer
->doInt(&ir
->nsttcouple
);
1292 ir
->nsttcouple
= ir
->nstcalcenergy
;
1294 serializer
->doInt(&ir
->epc
);
1295 serializer
->doInt(&ir
->epct
);
1296 if (file_version
>= 71)
1298 serializer
->doInt(&ir
->nstpcouple
);
1302 ir
->nstpcouple
= ir
->nstcalcenergy
;
1304 serializer
->doReal(&ir
->tau_p
);
1305 serializer
->doRvec(&ir
->ref_p
[XX
]);
1306 serializer
->doRvec(&ir
->ref_p
[YY
]);
1307 serializer
->doRvec(&ir
->ref_p
[ZZ
]);
1308 serializer
->doRvec(&ir
->compress
[XX
]);
1309 serializer
->doRvec(&ir
->compress
[YY
]);
1310 serializer
->doRvec(&ir
->compress
[ZZ
]);
1311 serializer
->doInt(&ir
->refcoord_scaling
);
1312 serializer
->doRvec(&ir
->posres_com
);
1313 serializer
->doRvec(&ir
->posres_comB
);
1315 if (file_version
< 79)
1317 serializer
->doInt(&ir
->andersen_seed
);
1321 ir
->andersen_seed
= 0;
1324 serializer
->doReal(&ir
->shake_tol
);
1326 serializer
->doInt(&ir
->efep
);
1327 do_fepvals(serializer
, ir
->fepvals
, file_version
);
1329 if (file_version
>= 79)
1331 serializer
->doBool(&ir
->bSimTemp
);
1334 ir
->bSimTemp
= TRUE
;
1339 ir
->bSimTemp
= FALSE
;
1343 do_simtempvals(serializer
, ir
->simtempvals
, ir
->fepvals
->n_lambda
, file_version
);
1346 if (file_version
>= 79)
1348 serializer
->doBool(&ir
->bExpanded
);
1351 ir
->bExpanded
= TRUE
;
1355 ir
->bExpanded
= FALSE
;
1360 do_expandedvals(serializer
, ir
->expandedvals
, ir
->fepvals
, file_version
);
1363 serializer
->doInt(&ir
->eDisre
);
1364 serializer
->doInt(&ir
->eDisreWeighting
);
1365 serializer
->doBool(&ir
->bDisreMixed
);
1366 serializer
->doReal(&ir
->dr_fc
);
1367 serializer
->doReal(&ir
->dr_tau
);
1368 serializer
->doInt(&ir
->nstdisreout
);
1369 serializer
->doReal(&ir
->orires_fc
);
1370 serializer
->doReal(&ir
->orires_tau
);
1371 serializer
->doInt(&ir
->nstorireout
);
1373 /* ignore dihre_fc */
1374 if (file_version
< 79)
1376 serializer
->doReal(&rdum
);
1379 serializer
->doReal(&ir
->em_stepsize
);
1380 serializer
->doReal(&ir
->em_tol
);
1381 serializer
->doBool(&ir
->bShakeSOR
);
1382 serializer
->doInt(&ir
->niter
);
1383 serializer
->doReal(&ir
->fc_stepsize
);
1384 serializer
->doInt(&ir
->eConstrAlg
);
1385 serializer
->doInt(&ir
->nProjOrder
);
1386 serializer
->doReal(&ir
->LincsWarnAngle
);
1387 serializer
->doInt(&ir
->nLincsIter
);
1388 serializer
->doReal(&ir
->bd_fric
);
1389 if (file_version
>= tpxv_Use64BitRandomSeed
)
1391 serializer
->doInt64(&ir
->ld_seed
);
1395 serializer
->doInt(&idum
);
1399 for (i
= 0; i
< DIM
; i
++)
1401 serializer
->doRvec(&ir
->deform
[i
]);
1403 serializer
->doReal(&ir
->cos_accel
);
1405 serializer
->doInt(&ir
->userint1
);
1406 serializer
->doInt(&ir
->userint2
);
1407 serializer
->doInt(&ir
->userint3
);
1408 serializer
->doInt(&ir
->userint4
);
1409 serializer
->doReal(&ir
->userreal1
);
1410 serializer
->doReal(&ir
->userreal2
);
1411 serializer
->doReal(&ir
->userreal3
);
1412 serializer
->doReal(&ir
->userreal4
);
1414 /* AdResS is removed, but we need to be able to read old files,
1415 and let mdrun refuse to run them */
1416 if (file_version
>= 77 && file_version
< tpxv_RemoveAdress
)
1418 serializer
->doBool(&ir
->bAdress
);
1421 int idum
, numThermoForceGroups
, numEnergyGroups
;
1424 serializer
->doInt(&idum
);
1425 serializer
->doReal(&rdum
);
1426 serializer
->doReal(&rdum
);
1427 serializer
->doReal(&rdum
);
1428 serializer
->doInt(&idum
);
1429 serializer
->doInt(&idum
);
1430 serializer
->doRvec(&rvecdum
);
1431 serializer
->doInt(&numThermoForceGroups
);
1432 serializer
->doReal(&rdum
);
1433 serializer
->doInt(&numEnergyGroups
);
1434 serializer
->doInt(&idum
);
1436 if (numThermoForceGroups
> 0)
1438 std::vector
<int> idumn(numThermoForceGroups
);
1439 serializer
->doIntArray(idumn
.data(), idumn
.size());
1441 if (numEnergyGroups
> 0)
1443 std::vector
<int> idumn(numEnergyGroups
);
1444 serializer
->doIntArray(idumn
.data(), idumn
.size());
1450 ir
->bAdress
= FALSE
;
1457 if (file_version
>= tpxv_PullCoordTypeGeom
)
1459 serializer
->doBool(&ir
->bPull
);
1463 serializer
->doInt(&ePullOld
);
1464 ir
->bPull
= (ePullOld
> 0);
1465 /* We removed the first ePull=ePullNo for the enum */
1470 if (serializer
->reading())
1474 do_pull(serializer
, ir
->pull
, file_version
, ePullOld
);
1478 if (file_version
>= tpxv_AcceleratedWeightHistogram
)
1480 serializer
->doBool(&ir
->bDoAwh
);
1484 if (serializer
->reading())
1486 snew(ir
->awhParams
, 1);
1488 do_awh(serializer
, ir
->awhParams
, file_version
);
1496 /* Enforced rotation */
1497 if (file_version
>= 74)
1499 serializer
->doBool(&ir
->bRot
);
1502 if (serializer
->reading())
1506 do_rot(serializer
, ir
->rot
);
1514 /* Interactive molecular dynamics */
1515 if (file_version
>= tpxv_InteractiveMolecularDynamics
)
1517 serializer
->doBool(&ir
->bIMD
);
1520 if (serializer
->reading())
1524 do_imd(serializer
, ir
->imd
);
1529 /* We don't support IMD sessions for old .tpr files */
1534 serializer
->doInt(&ir
->opts
.ngtc
);
1535 if (file_version
>= 69)
1537 serializer
->doInt(&ir
->opts
.nhchainlength
);
1541 ir
->opts
.nhchainlength
= 1;
1543 serializer
->doInt(&ir
->opts
.ngacc
);
1544 serializer
->doInt(&ir
->opts
.ngfrz
);
1545 serializer
->doInt(&ir
->opts
.ngener
);
1547 if (serializer
->reading())
1549 snew(ir
->opts
.nrdf
, ir
->opts
.ngtc
);
1550 snew(ir
->opts
.ref_t
, ir
->opts
.ngtc
);
1551 snew(ir
->opts
.annealing
, ir
->opts
.ngtc
);
1552 snew(ir
->opts
.anneal_npoints
, ir
->opts
.ngtc
);
1553 snew(ir
->opts
.anneal_time
, ir
->opts
.ngtc
);
1554 snew(ir
->opts
.anneal_temp
, ir
->opts
.ngtc
);
1555 snew(ir
->opts
.tau_t
, ir
->opts
.ngtc
);
1556 snew(ir
->opts
.nFreeze
, ir
->opts
.ngfrz
);
1557 snew(ir
->opts
.acc
, ir
->opts
.ngacc
);
1558 snew(ir
->opts
.egp_flags
, ir
->opts
.ngener
* ir
->opts
.ngener
);
1560 if (ir
->opts
.ngtc
> 0)
1562 serializer
->doRealArray(ir
->opts
.nrdf
, ir
->opts
.ngtc
);
1563 serializer
->doRealArray(ir
->opts
.ref_t
, ir
->opts
.ngtc
);
1564 serializer
->doRealArray(ir
->opts
.tau_t
, ir
->opts
.ngtc
);
1566 if (ir
->opts
.ngfrz
> 0)
1568 serializer
->doIvecArray(ir
->opts
.nFreeze
, ir
->opts
.ngfrz
);
1570 if (ir
->opts
.ngacc
> 0)
1572 serializer
->doRvecArray(ir
->opts
.acc
, ir
->opts
.ngacc
);
1574 serializer
->doIntArray(ir
->opts
.egp_flags
, ir
->opts
.ngener
* ir
->opts
.ngener
);
1576 /* First read the lists with annealing and npoints for each group */
1577 serializer
->doIntArray(ir
->opts
.annealing
, ir
->opts
.ngtc
);
1578 serializer
->doIntArray(ir
->opts
.anneal_npoints
, ir
->opts
.ngtc
);
1579 for (j
= 0; j
< (ir
->opts
.ngtc
); j
++)
1581 k
= ir
->opts
.anneal_npoints
[j
];
1582 if (serializer
->reading())
1584 snew(ir
->opts
.anneal_time
[j
], k
);
1585 snew(ir
->opts
.anneal_temp
[j
], k
);
1587 serializer
->doRealArray(ir
->opts
.anneal_time
[j
], k
);
1588 serializer
->doRealArray(ir
->opts
.anneal_temp
[j
], k
);
1592 serializer
->doInt(&ir
->nwall
);
1593 serializer
->doInt(&ir
->wall_type
);
1594 serializer
->doReal(&ir
->wall_r_linpot
);
1595 serializer
->doInt(&ir
->wall_atomtype
[0]);
1596 serializer
->doInt(&ir
->wall_atomtype
[1]);
1597 serializer
->doReal(&ir
->wall_density
[0]);
1598 serializer
->doReal(&ir
->wall_density
[1]);
1599 serializer
->doReal(&ir
->wall_ewald_zfac
);
1602 /* Cosine stuff for electric fields */
1603 if (file_version
< tpxv_GenericParamsForElectricField
)
1605 do_legacy_efield(serializer
, ¶msObj
);
1609 if (file_version
>= tpxv_ComputationalElectrophysiology
)
1611 serializer
->doInt(&ir
->eSwapCoords
);
1612 if (ir
->eSwapCoords
!= eswapNO
)
1614 if (serializer
->reading())
1618 do_swapcoords_tpx(serializer
, ir
->swap
, file_version
);
1624 serializer
->doBool(&ir
->bQMMM
);
1625 serializer
->doInt(&ir
->QMMMscheme
);
1626 serializer
->doReal(&ir
->scalefactor
);
1627 serializer
->doInt(&ir
->opts
.ngQM
);
1628 if (serializer
->reading())
1630 snew(ir
->opts
.QMmethod
, ir
->opts
.ngQM
);
1631 snew(ir
->opts
.QMbasis
, ir
->opts
.ngQM
);
1632 snew(ir
->opts
.QMcharge
, ir
->opts
.ngQM
);
1633 snew(ir
->opts
.QMmult
, ir
->opts
.ngQM
);
1634 snew(ir
->opts
.bSH
, ir
->opts
.ngQM
);
1635 snew(ir
->opts
.CASorbitals
, ir
->opts
.ngQM
);
1636 snew(ir
->opts
.CASelectrons
, ir
->opts
.ngQM
);
1637 snew(ir
->opts
.SAon
, ir
->opts
.ngQM
);
1638 snew(ir
->opts
.SAoff
, ir
->opts
.ngQM
);
1639 snew(ir
->opts
.SAsteps
, ir
->opts
.ngQM
);
1641 if (ir
->opts
.ngQM
> 0 && ir
->bQMMM
)
1643 serializer
->doIntArray(ir
->opts
.QMmethod
, ir
->opts
.ngQM
);
1644 serializer
->doIntArray(ir
->opts
.QMbasis
, ir
->opts
.ngQM
);
1645 serializer
->doIntArray(ir
->opts
.QMcharge
, ir
->opts
.ngQM
);
1646 serializer
->doIntArray(ir
->opts
.QMmult
, ir
->opts
.ngQM
);
1647 serializer
->doBoolArray(ir
->opts
.bSH
, ir
->opts
.ngQM
);
1648 serializer
->doIntArray(ir
->opts
.CASorbitals
, ir
->opts
.ngQM
);
1649 serializer
->doIntArray(ir
->opts
.CASelectrons
, ir
->opts
.ngQM
);
1650 serializer
->doRealArray(ir
->opts
.SAon
, ir
->opts
.ngQM
);
1651 serializer
->doRealArray(ir
->opts
.SAoff
, ir
->opts
.ngQM
);
1652 serializer
->doIntArray(ir
->opts
.SAsteps
, ir
->opts
.ngQM
);
1653 /* We leave in dummy i/o for removed parameters to avoid
1654 * changing the tpr format for every QMMM change.
1656 std::vector
<int> dummy(ir
->opts
.ngQM
, 0);
1657 serializer
->doIntArray(dummy
.data(), ir
->opts
.ngQM
);
1658 serializer
->doIntArray(dummy
.data(), ir
->opts
.ngQM
);
1660 /* end of QMMM stuff */
1663 if (file_version
>= tpxv_GenericParamsForElectricField
)
1665 if (serializer
->reading())
1667 paramsObj
.mergeObject(gmx::deserializeKeyValueTree(serializer
));
1671 GMX_RELEASE_ASSERT(ir
->params
!= nullptr,
1672 "Parameters should be present when writing inputrec");
1673 gmx::serializeKeyValueTree(*ir
->params
, serializer
);
1676 if (serializer
->reading())
1678 ir
->params
= new gmx::KeyValueTreeObject(paramsBuilder
.build());
1679 // Initialize internal parameters to an empty kvt for all tpr versions
1680 ir
->internalParameters
= std::make_unique
<gmx::KeyValueTreeObject
>();
1683 if (file_version
>= tpxv_GenericInternalParameters
)
1685 if (serializer
->reading())
1687 ir
->internalParameters
=
1688 std::make_unique
<gmx::KeyValueTreeObject
>(gmx::deserializeKeyValueTree(serializer
));
1692 GMX_RELEASE_ASSERT(ir
->internalParameters
!= nullptr,
1693 "Parameters should be present when writing inputrec");
1694 gmx::serializeKeyValueTree(*ir
->internalParameters
, serializer
);
1700 static void do_harm(gmx::ISerializer
* serializer
, t_iparams
* iparams
)
1702 serializer
->doReal(&iparams
->harmonic
.rA
);
1703 serializer
->doReal(&iparams
->harmonic
.krA
);
1704 serializer
->doReal(&iparams
->harmonic
.rB
);
1705 serializer
->doReal(&iparams
->harmonic
.krB
);
1708 static void do_iparams(gmx::ISerializer
* serializer
, t_functype ftype
, t_iparams
* iparams
, int file_version
)
1721 do_harm(serializer
, iparams
);
1722 if ((ftype
== F_ANGRES
|| ftype
== F_ANGRESZ
) && serializer
->reading())
1724 /* Correct incorrect storage of parameters */
1725 iparams
->pdihs
.phiB
= iparams
->pdihs
.phiA
;
1726 iparams
->pdihs
.cpB
= iparams
->pdihs
.cpA
;
1730 serializer
->doReal(&iparams
->harmonic
.rA
);
1731 serializer
->doReal(&iparams
->harmonic
.krA
);
1733 case F_LINEAR_ANGLES
:
1734 serializer
->doReal(&iparams
->linangle
.klinA
);
1735 serializer
->doReal(&iparams
->linangle
.aA
);
1736 serializer
->doReal(&iparams
->linangle
.klinB
);
1737 serializer
->doReal(&iparams
->linangle
.aB
);
1740 serializer
->doReal(&iparams
->fene
.bm
);
1741 serializer
->doReal(&iparams
->fene
.kb
);
1745 serializer
->doReal(&iparams
->restraint
.lowA
);
1746 serializer
->doReal(&iparams
->restraint
.up1A
);
1747 serializer
->doReal(&iparams
->restraint
.up2A
);
1748 serializer
->doReal(&iparams
->restraint
.kA
);
1749 serializer
->doReal(&iparams
->restraint
.lowB
);
1750 serializer
->doReal(&iparams
->restraint
.up1B
);
1751 serializer
->doReal(&iparams
->restraint
.up2B
);
1752 serializer
->doReal(&iparams
->restraint
.kB
);
1758 serializer
->doReal(&iparams
->tab
.kA
);
1759 serializer
->doInt(&iparams
->tab
.table
);
1760 serializer
->doReal(&iparams
->tab
.kB
);
1762 case F_CROSS_BOND_BONDS
:
1763 serializer
->doReal(&iparams
->cross_bb
.r1e
);
1764 serializer
->doReal(&iparams
->cross_bb
.r2e
);
1765 serializer
->doReal(&iparams
->cross_bb
.krr
);
1767 case F_CROSS_BOND_ANGLES
:
1768 serializer
->doReal(&iparams
->cross_ba
.r1e
);
1769 serializer
->doReal(&iparams
->cross_ba
.r2e
);
1770 serializer
->doReal(&iparams
->cross_ba
.r3e
);
1771 serializer
->doReal(&iparams
->cross_ba
.krt
);
1773 case F_UREY_BRADLEY
:
1774 serializer
->doReal(&iparams
->u_b
.thetaA
);
1775 serializer
->doReal(&iparams
->u_b
.kthetaA
);
1776 serializer
->doReal(&iparams
->u_b
.r13A
);
1777 serializer
->doReal(&iparams
->u_b
.kUBA
);
1778 if (file_version
>= 79)
1780 serializer
->doReal(&iparams
->u_b
.thetaB
);
1781 serializer
->doReal(&iparams
->u_b
.kthetaB
);
1782 serializer
->doReal(&iparams
->u_b
.r13B
);
1783 serializer
->doReal(&iparams
->u_b
.kUBB
);
1787 iparams
->u_b
.thetaB
= iparams
->u_b
.thetaA
;
1788 iparams
->u_b
.kthetaB
= iparams
->u_b
.kthetaA
;
1789 iparams
->u_b
.r13B
= iparams
->u_b
.r13A
;
1790 iparams
->u_b
.kUBB
= iparams
->u_b
.kUBA
;
1793 case F_QUARTIC_ANGLES
:
1794 serializer
->doReal(&iparams
->qangle
.theta
);
1795 serializer
->doRealArray(iparams
->qangle
.c
, 5);
1798 serializer
->doReal(&iparams
->bham
.a
);
1799 serializer
->doReal(&iparams
->bham
.b
);
1800 serializer
->doReal(&iparams
->bham
.c
);
1803 serializer
->doReal(&iparams
->morse
.b0A
);
1804 serializer
->doReal(&iparams
->morse
.cbA
);
1805 serializer
->doReal(&iparams
->morse
.betaA
);
1806 if (file_version
>= 79)
1808 serializer
->doReal(&iparams
->morse
.b0B
);
1809 serializer
->doReal(&iparams
->morse
.cbB
);
1810 serializer
->doReal(&iparams
->morse
.betaB
);
1814 iparams
->morse
.b0B
= iparams
->morse
.b0A
;
1815 iparams
->morse
.cbB
= iparams
->morse
.cbA
;
1816 iparams
->morse
.betaB
= iparams
->morse
.betaA
;
1820 serializer
->doReal(&iparams
->cubic
.b0
);
1821 serializer
->doReal(&iparams
->cubic
.kb
);
1822 serializer
->doReal(&iparams
->cubic
.kcub
);
1824 case F_CONNBONDS
: break;
1825 case F_POLARIZATION
: serializer
->doReal(&iparams
->polarize
.alpha
); break;
1827 serializer
->doReal(&iparams
->anharm_polarize
.alpha
);
1828 serializer
->doReal(&iparams
->anharm_polarize
.drcut
);
1829 serializer
->doReal(&iparams
->anharm_polarize
.khyp
);
1832 serializer
->doReal(&iparams
->wpol
.al_x
);
1833 serializer
->doReal(&iparams
->wpol
.al_y
);
1834 serializer
->doReal(&iparams
->wpol
.al_z
);
1835 serializer
->doReal(&iparams
->wpol
.rOH
);
1836 serializer
->doReal(&iparams
->wpol
.rHH
);
1837 serializer
->doReal(&iparams
->wpol
.rOD
);
1840 serializer
->doReal(&iparams
->thole
.a
);
1841 serializer
->doReal(&iparams
->thole
.alpha1
);
1842 serializer
->doReal(&iparams
->thole
.alpha2
);
1843 serializer
->doReal(&iparams
->thole
.rfac
);
1846 serializer
->doReal(&iparams
->lj
.c6
);
1847 serializer
->doReal(&iparams
->lj
.c12
);
1850 serializer
->doReal(&iparams
->lj14
.c6A
);
1851 serializer
->doReal(&iparams
->lj14
.c12A
);
1852 serializer
->doReal(&iparams
->lj14
.c6B
);
1853 serializer
->doReal(&iparams
->lj14
.c12B
);
1856 serializer
->doReal(&iparams
->ljc14
.fqq
);
1857 serializer
->doReal(&iparams
->ljc14
.qi
);
1858 serializer
->doReal(&iparams
->ljc14
.qj
);
1859 serializer
->doReal(&iparams
->ljc14
.c6
);
1860 serializer
->doReal(&iparams
->ljc14
.c12
);
1862 case F_LJC_PAIRS_NB
:
1863 serializer
->doReal(&iparams
->ljcnb
.qi
);
1864 serializer
->doReal(&iparams
->ljcnb
.qj
);
1865 serializer
->doReal(&iparams
->ljcnb
.c6
);
1866 serializer
->doReal(&iparams
->ljcnb
.c12
);
1872 serializer
->doReal(&iparams
->pdihs
.phiA
);
1873 serializer
->doReal(&iparams
->pdihs
.cpA
);
1874 serializer
->doReal(&iparams
->pdihs
.phiB
);
1875 serializer
->doReal(&iparams
->pdihs
.cpB
);
1876 serializer
->doInt(&iparams
->pdihs
.mult
);
1879 serializer
->doReal(&iparams
->pdihs
.phiA
);
1880 serializer
->doReal(&iparams
->pdihs
.cpA
);
1883 serializer
->doInt(&iparams
->disres
.label
);
1884 serializer
->doInt(&iparams
->disres
.type
);
1885 serializer
->doReal(&iparams
->disres
.low
);
1886 serializer
->doReal(&iparams
->disres
.up1
);
1887 serializer
->doReal(&iparams
->disres
.up2
);
1888 serializer
->doReal(&iparams
->disres
.kfac
);
1891 serializer
->doInt(&iparams
->orires
.ex
);
1892 serializer
->doInt(&iparams
->orires
.label
);
1893 serializer
->doInt(&iparams
->orires
.power
);
1894 serializer
->doReal(&iparams
->orires
.c
);
1895 serializer
->doReal(&iparams
->orires
.obs
);
1896 serializer
->doReal(&iparams
->orires
.kfac
);
1899 if (file_version
< 82)
1901 serializer
->doInt(&idum
);
1902 serializer
->doInt(&idum
);
1904 serializer
->doReal(&iparams
->dihres
.phiA
);
1905 serializer
->doReal(&iparams
->dihres
.dphiA
);
1906 serializer
->doReal(&iparams
->dihres
.kfacA
);
1907 if (file_version
>= 82)
1909 serializer
->doReal(&iparams
->dihres
.phiB
);
1910 serializer
->doReal(&iparams
->dihres
.dphiB
);
1911 serializer
->doReal(&iparams
->dihres
.kfacB
);
1915 iparams
->dihres
.phiB
= iparams
->dihres
.phiA
;
1916 iparams
->dihres
.dphiB
= iparams
->dihres
.dphiA
;
1917 iparams
->dihres
.kfacB
= iparams
->dihres
.kfacA
;
1921 serializer
->doRvec(&iparams
->posres
.pos0A
);
1922 serializer
->doRvec(&iparams
->posres
.fcA
);
1923 serializer
->doRvec(&iparams
->posres
.pos0B
);
1924 serializer
->doRvec(&iparams
->posres
.fcB
);
1927 serializer
->doInt(&iparams
->fbposres
.geom
);
1928 serializer
->doRvec(&iparams
->fbposres
.pos0
);
1929 serializer
->doReal(&iparams
->fbposres
.r
);
1930 serializer
->doReal(&iparams
->fbposres
.k
);
1932 case F_CBTDIHS
: serializer
->doRealArray(iparams
->cbtdihs
.cbtcA
, NR_CBTDIHS
); break;
1934 serializer
->doRealArray(iparams
->rbdihs
.rbcA
, NR_RBDIHS
);
1935 serializer
->doRealArray(iparams
->rbdihs
.rbcB
, NR_RBDIHS
);
1938 /* Fourier dihedrals are internally represented
1939 * as Ryckaert-Bellemans since those are faster to compute.
1941 serializer
->doRealArray(iparams
->rbdihs
.rbcA
, NR_RBDIHS
);
1942 serializer
->doRealArray(iparams
->rbdihs
.rbcB
, NR_RBDIHS
);
1946 serializer
->doReal(&iparams
->constr
.dA
);
1947 serializer
->doReal(&iparams
->constr
.dB
);
1950 serializer
->doReal(&iparams
->settle
.doh
);
1951 serializer
->doReal(&iparams
->settle
.dhh
);
1954 case F_VSITE2FD
: serializer
->doReal(&iparams
->vsite
.a
); break;
1958 serializer
->doReal(&iparams
->vsite
.a
);
1959 serializer
->doReal(&iparams
->vsite
.b
);
1964 serializer
->doReal(&iparams
->vsite
.a
);
1965 serializer
->doReal(&iparams
->vsite
.b
);
1966 serializer
->doReal(&iparams
->vsite
.c
);
1969 serializer
->doInt(&iparams
->vsiten
.n
);
1970 serializer
->doReal(&iparams
->vsiten
.a
);
1972 case F_GB12_NOLONGERUSED
:
1973 case F_GB13_NOLONGERUSED
:
1974 case F_GB14_NOLONGERUSED
:
1975 // Implicit solvent parameters can still be read, but never used
1976 if (serializer
->reading())
1978 if (file_version
< 68)
1980 serializer
->doReal(&rdum
);
1981 serializer
->doReal(&rdum
);
1982 serializer
->doReal(&rdum
);
1983 serializer
->doReal(&rdum
);
1985 if (file_version
< tpxv_RemoveImplicitSolvation
)
1987 serializer
->doReal(&rdum
);
1988 serializer
->doReal(&rdum
);
1989 serializer
->doReal(&rdum
);
1990 serializer
->doReal(&rdum
);
1991 serializer
->doReal(&rdum
);
1996 serializer
->doInt(&iparams
->cmap
.cmapA
);
1997 serializer
->doInt(&iparams
->cmap
.cmapB
);
2000 gmx_fatal(FARGS
, "unknown function type %d (%s) in %s line %d", ftype
,
2001 interaction_function
[ftype
].name
, __FILE__
, __LINE__
);
2005 static void do_ilist(gmx::ISerializer
* serializer
, InteractionList
* ilist
)
2007 int nr
= ilist
->size();
2008 serializer
->doInt(&nr
);
2009 if (serializer
->reading())
2011 ilist
->iatoms
.resize(nr
);
2013 serializer
->doIntArray(ilist
->iatoms
.data(), ilist
->size());
2016 static void do_ffparams(gmx::ISerializer
* serializer
, gmx_ffparams_t
* ffparams
, int file_version
)
2018 serializer
->doInt(&ffparams
->atnr
);
2019 int numTypes
= ffparams
->numTypes();
2020 serializer
->doInt(&numTypes
);
2021 if (serializer
->reading())
2023 ffparams
->functype
.resize(numTypes
);
2024 ffparams
->iparams
.resize(numTypes
);
2026 /* Read/write all the function types */
2027 serializer
->doIntArray(ffparams
->functype
.data(), ffparams
->functype
.size());
2029 if (file_version
>= 66)
2031 serializer
->doDouble(&ffparams
->reppow
);
2035 ffparams
->reppow
= 12.0;
2038 serializer
->doReal(&ffparams
->fudgeQQ
);
2040 /* Check whether all these function types are supported by the code.
2041 * In practice the code is backwards compatible, which means that the
2042 * numbering may have to be altered from old numbering to new numbering
2044 for (int i
= 0; i
< ffparams
->numTypes(); i
++)
2046 if (serializer
->reading())
2048 /* Loop over file versions */
2049 for (int k
= 0; k
< NFTUPD
; k
++)
2051 /* Compare the read file_version to the update table */
2052 if ((file_version
< ftupd
[k
].fvnr
) && (ffparams
->functype
[i
] >= ftupd
[k
].ftype
))
2054 ffparams
->functype
[i
] += 1;
2059 do_iparams(serializer
, ffparams
->functype
[i
], &ffparams
->iparams
[i
], file_version
);
2063 static void add_settle_atoms(InteractionList
* ilist
)
2067 /* Settle used to only store the first atom: add the other two */
2068 ilist
->iatoms
.resize(2 * ilist
->size());
2069 for (i
= ilist
->size() / 4 - 1; i
>= 0; i
--)
2071 ilist
->iatoms
[4 * i
+ 0] = ilist
->iatoms
[2 * i
+ 0];
2072 ilist
->iatoms
[4 * i
+ 1] = ilist
->iatoms
[2 * i
+ 1];
2073 ilist
->iatoms
[4 * i
+ 2] = ilist
->iatoms
[2 * i
+ 1] + 1;
2074 ilist
->iatoms
[4 * i
+ 3] = ilist
->iatoms
[2 * i
+ 1] + 2;
2078 static void do_ilists(gmx::ISerializer
* serializer
, InteractionLists
* ilists
, int file_version
)
2080 GMX_RELEASE_ASSERT(ilists
, "Need a valid ilists object");
2081 GMX_RELEASE_ASSERT(ilists
->size() == F_NRE
,
2082 "The code needs to be in sync with InteractionLists");
2084 for (int j
= 0; j
< F_NRE
; j
++)
2086 InteractionList
& ilist
= (*ilists
)[j
];
2087 gmx_bool bClear
= FALSE
;
2088 if (serializer
->reading())
2090 for (int k
= 0; k
< NFTUPD
; k
++)
2092 if ((file_version
< ftupd
[k
].fvnr
) && (j
== ftupd
[k
].ftype
))
2100 ilist
.iatoms
.clear();
2104 do_ilist(serializer
, &ilist
);
2105 if (file_version
< 78 && j
== F_SETTLE
&& ilist
.size() > 0)
2107 add_settle_atoms(&ilist
);
2113 static void do_block(gmx::ISerializer
* serializer
, t_block
* block
)
2115 serializer
->doInt(&block
->nr
);
2116 if (serializer
->reading())
2118 if ((block
->nalloc_index
> 0) && (nullptr != block
->index
))
2120 sfree(block
->index
);
2122 block
->nalloc_index
= block
->nr
+ 1;
2123 snew(block
->index
, block
->nalloc_index
);
2125 serializer
->doIntArray(block
->index
, block
->nr
+ 1);
2128 static void do_blocka(gmx::ISerializer
* serializer
, t_blocka
* block
)
2130 serializer
->doInt(&block
->nr
);
2131 serializer
->doInt(&block
->nra
);
2132 if (serializer
->reading())
2134 block
->nalloc_index
= block
->nr
+ 1;
2135 snew(block
->index
, block
->nalloc_index
);
2136 block
->nalloc_a
= block
->nra
;
2137 snew(block
->a
, block
->nalloc_a
);
2139 serializer
->doIntArray(block
->index
, block
->nr
+ 1);
2140 serializer
->doIntArray(block
->a
, block
->nra
);
2143 /* This is a primitive routine to make it possible to translate atomic numbers
2144 * to element names when reading TPR files, without making the Gromacs library
2145 * directory a dependency on mdrun (which is the case if we need elements.dat).
2147 static const char* atomicnumber_to_element(int atomicnumber
)
2151 /* This does not have to be complete, so we only include elements likely
2152 * to occur in PDB files.
2154 switch (atomicnumber
)
2156 case 1: p
= "H"; break;
2157 case 5: p
= "B"; break;
2158 case 6: p
= "C"; break;
2159 case 7: p
= "N"; break;
2160 case 8: p
= "O"; break;
2161 case 9: p
= "F"; break;
2162 case 11: p
= "Na"; break;
2163 case 12: p
= "Mg"; break;
2164 case 15: p
= "P"; break;
2165 case 16: p
= "S"; break;
2166 case 17: p
= "Cl"; break;
2167 case 18: p
= "Ar"; break;
2168 case 19: p
= "K"; break;
2169 case 20: p
= "Ca"; break;
2170 case 25: p
= "Mn"; break;
2171 case 26: p
= "Fe"; break;
2172 case 28: p
= "Ni"; break;
2173 case 29: p
= "Cu"; break;
2174 case 30: p
= "Zn"; break;
2175 case 35: p
= "Br"; break;
2176 case 47: p
= "Ag"; break;
2177 default: p
= ""; break;
2183 static void do_atom(gmx::ISerializer
* serializer
, t_atom
* atom
)
2185 serializer
->doReal(&atom
->m
);
2186 serializer
->doReal(&atom
->q
);
2187 serializer
->doReal(&atom
->mB
);
2188 serializer
->doReal(&atom
->qB
);
2189 serializer
->doUShort(&atom
->type
);
2190 serializer
->doUShort(&atom
->typeB
);
2191 serializer
->doInt(&atom
->ptype
);
2192 serializer
->doInt(&atom
->resind
);
2193 serializer
->doInt(&atom
->atomnumber
);
2194 if (serializer
->reading())
2196 /* Set element string from atomic number if present.
2197 * This routine returns an empty string if the name is not found.
2199 std::strncpy(atom
->elem
, atomicnumber_to_element(atom
->atomnumber
), 4);
2200 /* avoid warnings about potentially unterminated string */
2201 atom
->elem
[3] = '\0';
2205 static void do_grps(gmx::ISerializer
* serializer
, gmx::ArrayRef
<AtomGroupIndices
> grps
)
2207 for (auto& group
: grps
)
2209 int size
= group
.size();
2210 serializer
->doInt(&size
);
2211 if (serializer
->reading())
2215 serializer
->doIntArray(group
.data(), size
);
2219 static void do_symstr(gmx::ISerializer
* serializer
, char*** nm
, t_symtab
* symtab
)
2223 if (serializer
->reading())
2225 serializer
->doInt(&ls
);
2226 *nm
= get_symtab_handle(symtab
, ls
);
2230 ls
= lookup_symtab(symtab
, *nm
);
2231 serializer
->doInt(&ls
);
2235 static void do_strstr(gmx::ISerializer
* serializer
, int nstr
, char*** nm
, t_symtab
* symtab
)
2239 for (j
= 0; (j
< nstr
); j
++)
2241 do_symstr(serializer
, &(nm
[j
]), symtab
);
2245 static void do_resinfo(gmx::ISerializer
* serializer
, int n
, t_resinfo
* ri
, t_symtab
* symtab
, int file_version
)
2249 for (j
= 0; (j
< n
); j
++)
2251 do_symstr(serializer
, &(ri
[j
].name
), symtab
);
2252 if (file_version
>= 63)
2254 serializer
->doInt(&ri
[j
].nr
);
2255 serializer
->doUChar(&ri
[j
].ic
);
2265 static void do_atoms(gmx::ISerializer
* serializer
, t_atoms
* atoms
, t_symtab
* symtab
, int file_version
)
2269 serializer
->doInt(&atoms
->nr
);
2270 serializer
->doInt(&atoms
->nres
);
2271 if (serializer
->reading())
2273 /* Since we have always written all t_atom properties in the tpr file
2274 * (at least for all backward compatible versions), we don't store
2275 * but simple set the booleans here.
2277 atoms
->haveMass
= TRUE
;
2278 atoms
->haveCharge
= TRUE
;
2279 atoms
->haveType
= TRUE
;
2280 atoms
->haveBState
= TRUE
;
2281 atoms
->havePdbInfo
= FALSE
;
2283 snew(atoms
->atom
, atoms
->nr
);
2284 snew(atoms
->atomname
, atoms
->nr
);
2285 snew(atoms
->atomtype
, atoms
->nr
);
2286 snew(atoms
->atomtypeB
, atoms
->nr
);
2287 snew(atoms
->resinfo
, atoms
->nres
);
2288 atoms
->pdbinfo
= nullptr;
2292 GMX_RELEASE_ASSERT(atoms
->haveMass
&& atoms
->haveCharge
&& atoms
->haveType
&& atoms
->haveBState
,
2293 "Mass, charge, atomtype and B-state parameters should be present in "
2294 "t_atoms when writing a tpr file");
2296 for (i
= 0; (i
< atoms
->nr
); i
++)
2298 do_atom(serializer
, &atoms
->atom
[i
]);
2300 do_strstr(serializer
, atoms
->nr
, atoms
->atomname
, symtab
);
2301 do_strstr(serializer
, atoms
->nr
, atoms
->atomtype
, symtab
);
2302 do_strstr(serializer
, atoms
->nr
, atoms
->atomtypeB
, symtab
);
2304 do_resinfo(serializer
, atoms
->nres
, atoms
->resinfo
, symtab
, file_version
);
2307 static void do_groups(gmx::ISerializer
* serializer
, SimulationGroups
* groups
, t_symtab
* symtab
)
2309 do_grps(serializer
, groups
->groups
);
2310 int numberOfGroupNames
= groups
->groupNames
.size();
2311 serializer
->doInt(&numberOfGroupNames
);
2312 if (serializer
->reading())
2314 groups
->groupNames
.resize(numberOfGroupNames
);
2316 do_strstr(serializer
, numberOfGroupNames
, groups
->groupNames
.data(), symtab
);
2317 for (auto group
: gmx::keysOf(groups
->groupNumbers
))
2319 int numberOfGroupNumbers
= groups
->numberOfGroupNumbers(group
);
2320 serializer
->doInt(&numberOfGroupNumbers
);
2321 if (numberOfGroupNumbers
!= 0)
2323 if (serializer
->reading())
2325 groups
->groupNumbers
[group
].resize(numberOfGroupNumbers
);
2327 serializer
->doUCharArray(groups
->groupNumbers
[group
].data(), numberOfGroupNumbers
);
2332 static void do_atomtypes(gmx::ISerializer
* serializer
, t_atomtypes
* atomtypes
, int file_version
)
2336 serializer
->doInt(&atomtypes
->nr
);
2338 if (serializer
->reading())
2340 snew(atomtypes
->atomnumber
, j
);
2342 if (serializer
->reading() && file_version
< tpxv_RemoveImplicitSolvation
)
2344 std::vector
<real
> dummy(atomtypes
->nr
, 0);
2345 serializer
->doRealArray(dummy
.data(), dummy
.size());
2346 serializer
->doRealArray(dummy
.data(), dummy
.size());
2347 serializer
->doRealArray(dummy
.data(), dummy
.size());
2349 serializer
->doIntArray(atomtypes
->atomnumber
, j
);
2351 if (serializer
->reading() && file_version
>= 60 && file_version
< tpxv_RemoveImplicitSolvation
)
2353 std::vector
<real
> dummy(atomtypes
->nr
, 0);
2354 serializer
->doRealArray(dummy
.data(), dummy
.size());
2355 serializer
->doRealArray(dummy
.data(), dummy
.size());
2359 static void do_symtab(gmx::ISerializer
* serializer
, t_symtab
* symtab
)
2364 serializer
->doInt(&symtab
->nr
);
2366 if (serializer
->reading())
2368 snew(symtab
->symbuf
, 1);
2369 symbuf
= symtab
->symbuf
;
2370 symbuf
->bufsize
= nr
;
2371 snew(symbuf
->buf
, nr
);
2372 for (i
= 0; (i
< nr
); i
++)
2375 serializer
->doString(&buf
);
2376 symbuf
->buf
[i
] = gmx_strdup(buf
.c_str());
2381 symbuf
= symtab
->symbuf
;
2382 while (symbuf
!= nullptr)
2384 for (i
= 0; (i
< symbuf
->bufsize
) && (i
< nr
); i
++)
2386 std::string buf
= symbuf
->buf
[i
];
2387 serializer
->doString(&buf
);
2390 symbuf
= symbuf
->next
;
2394 gmx_fatal(FARGS
, "nr of symtab strings left: %d", nr
);
2399 static void do_cmap(gmx::ISerializer
* serializer
, gmx_cmap_t
* cmap_grid
)
2402 int ngrid
= cmap_grid
->cmapdata
.size();
2403 serializer
->doInt(&ngrid
);
2404 serializer
->doInt(&cmap_grid
->grid_spacing
);
2406 int gs
= cmap_grid
->grid_spacing
;
2407 int nelem
= gs
* gs
;
2409 if (serializer
->reading())
2411 cmap_grid
->cmapdata
.resize(ngrid
);
2413 for (int i
= 0; i
< ngrid
; i
++)
2415 cmap_grid
->cmapdata
[i
].cmap
.resize(4 * nelem
);
2419 for (int i
= 0; i
< ngrid
; i
++)
2421 for (int j
= 0; j
< nelem
; j
++)
2423 serializer
->doReal(&cmap_grid
->cmapdata
[i
].cmap
[j
* 4]);
2424 serializer
->doReal(&cmap_grid
->cmapdata
[i
].cmap
[j
* 4 + 1]);
2425 serializer
->doReal(&cmap_grid
->cmapdata
[i
].cmap
[j
* 4 + 2]);
2426 serializer
->doReal(&cmap_grid
->cmapdata
[i
].cmap
[j
* 4 + 3]);
2432 static void do_moltype(gmx::ISerializer
* serializer
, gmx_moltype_t
* molt
, t_symtab
* symtab
, int file_version
)
2434 do_symstr(serializer
, &(molt
->name
), symtab
);
2436 do_atoms(serializer
, &molt
->atoms
, symtab
, file_version
);
2438 do_ilists(serializer
, &molt
->ilist
, file_version
);
2440 /* TODO: Remove the obsolete charge group index from the file */
2442 cgs
.nr
= molt
->atoms
.nr
;
2443 cgs
.nalloc_index
= cgs
.nr
+ 1;
2444 snew(cgs
.index
, cgs
.nalloc_index
);
2445 for (int i
= 0; i
< cgs
.nr
+ 1; i
++)
2449 do_block(serializer
, &cgs
);
2452 /* This used to be in the atoms struct */
2453 do_blocka(serializer
, &molt
->excls
);
2456 static void do_molblock(gmx::ISerializer
* serializer
, gmx_molblock_t
* molb
, int numAtomsPerMolecule
)
2458 serializer
->doInt(&molb
->type
);
2459 serializer
->doInt(&molb
->nmol
);
2460 /* To maintain forward topology reading compatibility, we store #atoms.
2461 * TODO: Change this to conditional reading of a dummy int when we
2462 * increase tpx_generation.
2464 serializer
->doInt(&numAtomsPerMolecule
);
2465 /* Position restraint coordinates */
2466 int numPosres_xA
= molb
->posres_xA
.size();
2467 serializer
->doInt(&numPosres_xA
);
2468 if (numPosres_xA
> 0)
2470 if (serializer
->reading())
2472 molb
->posres_xA
.resize(numPosres_xA
);
2474 serializer
->doRvecArray(as_rvec_array(molb
->posres_xA
.data()), numPosres_xA
);
2476 int numPosres_xB
= molb
->posres_xB
.size();
2477 serializer
->doInt(&numPosres_xB
);
2478 if (numPosres_xB
> 0)
2480 if (serializer
->reading())
2482 molb
->posres_xB
.resize(numPosres_xB
);
2484 serializer
->doRvecArray(as_rvec_array(molb
->posres_xB
.data()), numPosres_xB
);
2488 static void set_disres_npair(gmx_mtop_t
* mtop
)
2490 gmx_mtop_ilistloop_t iloop
;
2493 gmx::ArrayRef
<t_iparams
> ip
= mtop
->ffparams
.iparams
;
2495 iloop
= gmx_mtop_ilistloop_init(mtop
);
2496 while (const InteractionLists
* ilist
= gmx_mtop_ilistloop_next(iloop
, &nmol
))
2498 const InteractionList
& il
= (*ilist
)[F_DISRES
];
2502 gmx::ArrayRef
<const int> a
= il
.iatoms
;
2504 for (int i
= 0; i
< il
.size(); i
+= 3)
2507 if (i
+ 3 == il
.size() || ip
[a
[i
]].disres
.label
!= ip
[a
[i
+ 3]].disres
.label
)
2509 ip
[a
[i
]].disres
.npair
= npair
;
2517 static void do_mtop(gmx::ISerializer
* serializer
, gmx_mtop_t
* mtop
, int file_version
)
2519 do_symtab(serializer
, &(mtop
->symtab
));
2521 do_symstr(serializer
, &(mtop
->name
), &(mtop
->symtab
));
2523 do_ffparams(serializer
, &mtop
->ffparams
, file_version
);
2525 int nmoltype
= mtop
->moltype
.size();
2526 serializer
->doInt(&nmoltype
);
2527 if (serializer
->reading())
2529 mtop
->moltype
.resize(nmoltype
);
2531 for (gmx_moltype_t
& moltype
: mtop
->moltype
)
2533 do_moltype(serializer
, &moltype
, &mtop
->symtab
, file_version
);
2536 int nmolblock
= mtop
->molblock
.size();
2537 serializer
->doInt(&nmolblock
);
2538 if (serializer
->reading())
2540 mtop
->molblock
.resize(nmolblock
);
2542 for (gmx_molblock_t
& molblock
: mtop
->molblock
)
2544 int numAtomsPerMolecule
= (serializer
->reading() ? 0 : mtop
->moltype
[molblock
.type
].atoms
.nr
);
2545 do_molblock(serializer
, &molblock
, numAtomsPerMolecule
);
2547 serializer
->doInt(&mtop
->natoms
);
2549 if (file_version
>= tpxv_IntermolecularBondeds
)
2551 serializer
->doBool(&mtop
->bIntermolecularInteractions
);
2552 if (mtop
->bIntermolecularInteractions
)
2554 if (serializer
->reading())
2556 mtop
->intermolecular_ilist
= std::make_unique
<InteractionLists
>();
2558 do_ilists(serializer
, mtop
->intermolecular_ilist
.get(), file_version
);
2563 mtop
->bIntermolecularInteractions
= FALSE
;
2566 do_atomtypes(serializer
, &(mtop
->atomtypes
), file_version
);
2568 if (file_version
>= 65)
2570 do_cmap(serializer
, &mtop
->ffparams
.cmap_grid
);
2574 mtop
->ffparams
.cmap_grid
.grid_spacing
= 0;
2575 mtop
->ffparams
.cmap_grid
.cmapdata
.clear();
2578 do_groups(serializer
, &mtop
->groups
, &(mtop
->symtab
));
2580 mtop
->haveMoleculeIndices
= true;
2582 if (serializer
->reading())
2584 close_symtab(&(mtop
->symtab
));
2589 * Read the first part of the TPR file to find general system information.
2591 * If \p TopOnlyOK is true then we can read even future versions
2592 * of tpx files, provided the \p fileGeneration hasn't changed.
2593 * If it is false, we need the \p ir too, and bail out
2594 * if the file is newer than the program.
2596 * The version and generation of the topology (see top of this file)
2597 * are returned in the two last arguments, if those arguments are non-nullptr.
2599 * If possible, we will read the \p ir even when \p TopOnlyOK is true.
2601 * \param[in,out] serializer The serializer used to handle header processing.
2602 * \param[in,out] tpx File header datastructure.
2603 * \param[in] filename The name of the file being read/written
2604 * \param[in,out] fio File handle.
2605 * \param[in] TopOnlyOK If not reading \p ir is fine or not.
2607 static void do_tpxheader(gmx::FileIOXdrSerializer
* serializer
,
2609 const char* filename
,
2617 /* XDR binary topology file */
2618 precision
= sizeof(real
);
2620 std::string fileTag
;
2621 if (serializer
->reading())
2623 serializer
->doString(&buf
);
2624 if (std::strncmp(buf
.c_str(), "VERSION", 7) != 0)
2628 "Can not read file %s,\n"
2629 " this file is from a GROMACS version which is older than 2.0\n"
2630 " Make a new one with grompp or use a gro or pdb file, if possible",
2633 // We need to know the precision used to write the TPR file, to match it
2634 // to the precision of the currently running binary. If the precisions match
2635 // there is no problem, but mismatching precision needs to be accounted for
2636 // by reading into temporary variables of the correct precision instead
2637 // of the desired target datastructures.
2638 serializer
->doInt(&precision
);
2639 tpx
->isDouble
= (precision
== sizeof(double));
2640 if ((precision
!= sizeof(float)) && !tpx
->isDouble
)
2643 "Unknown precision in file %s: real is %d bytes "
2644 "instead of %zu or %zu",
2645 filename
, precision
, sizeof(float), sizeof(double));
2647 gmx_fio_setprecision(fio
, tpx
->isDouble
);
2648 fprintf(stderr
, "Reading file %s, %s (%s precision)\n", filename
, buf
.c_str(),
2649 tpx
->isDouble
? "double" : "single");
2653 buf
= gmx::formatString("VERSION %s", gmx_version());
2654 serializer
->doString(&buf
);
2655 gmx_fio_setprecision(fio
, tpx
->isDouble
);
2656 serializer
->doInt(&precision
);
2657 fileTag
= gmx::formatString("%s", tpx_tag
);
2660 /* Check versions! */
2661 serializer
->doInt(&tpx
->fileVersion
);
2663 /* This is for backward compatibility with development versions 77-79
2664 * where the tag was, mistakenly, placed before the generation,
2665 * which would cause a segv instead of a proper error message
2666 * when reading the topology only from tpx with <77 code.
2668 if (tpx
->fileVersion
>= 77 && tpx
->fileVersion
<= 79)
2670 serializer
->doString(&fileTag
);
2673 serializer
->doInt(&tpx
->fileGeneration
);
2675 if (tpx
->fileVersion
>= 81)
2677 serializer
->doString(&fileTag
);
2679 if (serializer
->reading())
2681 if (tpx
->fileVersion
< 77)
2683 /* Versions before 77 don't have the tag, set it to release */
2684 fileTag
= gmx::formatString("%s", TPX_TAG_RELEASE
);
2687 if (fileTag
!= tpx_tag
)
2689 fprintf(stderr
, "Note: file tpx tag '%s', software tpx tag '%s'\n", fileTag
.c_str(), tpx_tag
);
2691 /* We only support reading tpx files with the same tag as the code
2692 * or tpx files with the release tag and with lower version number.
2694 if (fileTag
!= TPX_TAG_RELEASE
&& tpx
->fileVersion
< tpx_version
)
2697 "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' "
2698 "with program for tpx version %d, tag '%s'",
2699 filename
, tpx
->fileVersion
, fileTag
.c_str(), tpx_version
, tpx_tag
);
2704 if ((tpx
->fileVersion
<= tpx_incompatible_version
)
2705 || ((tpx
->fileVersion
> tpx_version
) && !TopOnlyOK
) || (tpx
->fileGeneration
> tpx_generation
)
2706 || tpx_version
== 80) /*80 was used by both 5.0-dev and 4.6-dev*/
2708 gmx_fatal(FARGS
, "reading tpx file (%s) version %d with version %d program", filename
,
2709 tpx
->fileVersion
, tpx_version
);
2712 serializer
->doInt(&tpx
->natoms
);
2713 serializer
->doInt(&tpx
->ngtc
);
2715 if (tpx
->fileVersion
< 62)
2717 serializer
->doInt(&idum
);
2718 serializer
->doReal(&rdum
);
2720 if (tpx
->fileVersion
>= 79)
2722 serializer
->doInt(&tpx
->fep_state
);
2724 serializer
->doReal(&tpx
->lambda
);
2725 serializer
->doBool(&tpx
->bIr
);
2726 serializer
->doBool(&tpx
->bTop
);
2727 serializer
->doBool(&tpx
->bX
);
2728 serializer
->doBool(&tpx
->bV
);
2729 serializer
->doBool(&tpx
->bF
);
2730 serializer
->doBool(&tpx
->bBox
);
2732 if (tpx
->fileVersion
>= tpxv_AddSizeField
&& tpx
->fileGeneration
>= 27)
2734 if (!serializer
->reading())
2736 GMX_RELEASE_ASSERT(tpx
->sizeOfTprBody
!= 0,
2737 "Not possible to write new file with zero TPR body size");
2739 serializer
->doInt64(&tpx
->sizeOfTprBody
);
2742 if ((tpx
->fileGeneration
> tpx_generation
))
2744 /* This can only happen if TopOnlyOK=TRUE */
2749 #define do_test(serializer, b, p) \
2750 if ((serializer)->reading() && ((p) != nullptr) && !(b)) \
2751 gmx_fatal(FARGS, "No %s in input file", #p)
2754 * Process the first part of the TPR into the state datastructure.
2756 * Due to the structure of the legacy code, it is necessary
2757 * to split up the state reading into two parts, with the
2758 * box and legacy temperature coupling processed before the
2759 * topology datastructures.
2761 * See the documentation for do_tpx_body for the correct order of
2762 * the operations for reading a tpr file.
2764 * \param[in] serializer Abstract serializer used to read/write data.
2765 * \param[in] tpx The file header data.
2766 * \param[in, out] state Global state data.
2768 static void do_tpx_state_first(gmx::ISerializer
* serializer
, TpxFileHeader
* tpx
, t_state
* state
)
2770 if (serializer
->reading())
2773 init_gtc_state(state
, tpx
->ngtc
, 0, 0);
2775 do_test(serializer
, tpx
->bBox
, state
->box
);
2778 serializer
->doRvecArray(state
->box
, DIM
);
2779 if (tpx
->fileVersion
>= 51)
2781 serializer
->doRvecArray(state
->box_rel
, DIM
);
2785 /* We initialize box_rel after reading the inputrec */
2786 clear_mat(state
->box_rel
);
2788 serializer
->doRvecArray(state
->boxv
, DIM
);
2789 if (tpx
->fileVersion
< 56)
2792 serializer
->doRvecArray(mdum
, DIM
);
2796 if (state
->ngtc
> 0)
2799 snew(dumv
, state
->ngtc
);
2800 if (tpx
->fileVersion
< 69)
2802 serializer
->doRealArray(dumv
, state
->ngtc
);
2804 /* These used to be the Berendsen tcoupl_lambda's */
2805 serializer
->doRealArray(dumv
, state
->ngtc
);
2811 * Process global topology data.
2813 * See the documentation for do_tpx_body for the correct order of
2814 * the operations for reading a tpr file.
2816 * \param[in] serializer Abstract serializer used to read/write data.
2817 * \param[in] tpx The file header data.
2818 * \param[in,out] mtop Global topology.
2820 static void do_tpx_mtop(gmx::ISerializer
* serializer
, TpxFileHeader
* tpx
, gmx_mtop_t
* mtop
)
2822 do_test(serializer
, tpx
->bTop
, mtop
);
2827 do_mtop(serializer
, mtop
, tpx
->fileVersion
);
2828 set_disres_npair(mtop
);
2829 gmx_mtop_finalize(mtop
);
2834 do_mtop(serializer
, &dum_top
, tpx
->fileVersion
);
2839 * Process coordinate vectors for state data.
2841 * Main part of state gets processed here.
2843 * See the documentation for do_tpx_body for the correct order of
2844 * the operations for reading a tpr file.
2846 * \param[in] serializer Abstract serializer used to read/write data.
2847 * \param[in] tpx The file header data.
2848 * \param[in,out] state Global state data.
2849 * \param[in,out] x Individual coordinates for processing, deprecated.
2850 * \param[in,out] v Individual velocities for processing, deprecated.
2852 static void do_tpx_state_second(gmx::ISerializer
* serializer
, TpxFileHeader
* tpx
, t_state
* state
, rvec
* x
, rvec
* v
)
2854 if (!serializer
->reading())
2857 x
== nullptr && v
== nullptr,
2858 "Passing separate x and v pointers to do_tpx() is not supported when writing");
2862 GMX_RELEASE_ASSERT(!(x
== nullptr && v
!= nullptr),
2863 "Passing x==NULL and v!=NULL is not supported");
2866 if (serializer
->reading())
2870 // v is also nullptr by the above assertion, so we may
2871 // need to make memory in state for storing the contents
2875 state
->flags
|= (1 << estX
);
2879 state
->flags
|= (1 << estV
);
2881 state_change_natoms(state
, tpx
->natoms
);
2887 x
= state
->x
.rvec_array();
2888 v
= state
->v
.rvec_array();
2890 do_test(serializer
, tpx
->bX
, x
);
2893 if (serializer
->reading())
2895 state
->flags
|= (1 << estX
);
2897 serializer
->doRvecArray(x
, tpx
->natoms
);
2900 do_test(serializer
, tpx
->bV
, v
);
2903 if (serializer
->reading())
2905 state
->flags
|= (1 << estV
);
2909 std::vector
<gmx::RVec
> dummyVelocities(tpx
->natoms
);
2910 serializer
->doRvecArray(as_rvec_array(dummyVelocities
.data()), tpx
->natoms
);
2914 serializer
->doRvecArray(v
, tpx
->natoms
);
2918 // No need to run do_test when the last argument is NULL
2921 std::vector
<gmx::RVec
> dummyForces(state
->natoms
);
2922 serializer
->doRvecArray(as_rvec_array(dummyForces
.data()), tpx
->natoms
);
2926 * Process simulation parameters.
2928 * See the documentation for do_tpx_body for the correct order of
2929 * the operations for reading a tpr file.
2931 * \param[in] serializer Abstract serializer used to read/write data.
2932 * \param[in] tpx The file header data.
2933 * \param[in,out] ir Datastructure with simulation parameters.
2935 static int do_tpx_ir(gmx::ISerializer
* serializer
, TpxFileHeader
* tpx
, t_inputrec
* ir
)
2938 gmx_bool bPeriodicMols
;
2940 /* Starting with tpx version 26, we have the inputrec
2941 * at the end of the file, so we can ignore it
2942 * if the file is never than the software (but still the
2943 * same generation - see comments at the top of this file.
2948 bPeriodicMols
= FALSE
;
2950 do_test(serializer
, tpx
->bIr
, ir
);
2953 if (tpx
->fileVersion
>= 53)
2955 /* Removed the pbc info from do_inputrec, since we always want it */
2956 if (!serializer
->reading())
2959 bPeriodicMols
= ir
->bPeriodicMols
;
2961 serializer
->doInt(&ePBC
);
2962 serializer
->doBool(&bPeriodicMols
);
2964 if (tpx
->fileGeneration
<= tpx_generation
&& ir
)
2966 do_inputrec(serializer
, ir
, tpx
->fileVersion
);
2967 if (tpx
->fileVersion
< 53)
2970 bPeriodicMols
= ir
->bPeriodicMols
;
2973 if (serializer
->reading() && ir
&& tpx
->fileVersion
>= 53)
2975 /* We need to do this after do_inputrec, since that initializes ir */
2977 ir
->bPeriodicMols
= bPeriodicMols
;
2984 * Correct and finalize read information.
2986 * If \p state is nullptr, skip the parts dependent on it.
2988 * See the documentation for do_tpx_body for the correct order of
2989 * the operations for reading a tpr file.
2991 * \param[in] tpx The file header used to check version numbers.
2992 * \param[out] ir Input rec that needs correction.
2993 * \param[out] state State needing correction.
2994 * \param[out] mtop Topology to finalize.
2996 static void do_tpx_finalize(TpxFileHeader
* tpx
, t_inputrec
* ir
, t_state
* state
, gmx_mtop_t
* mtop
)
2998 if (tpx
->fileVersion
< 51 && state
)
3000 set_box_rel(ir
, state
);
3004 if (state
&& state
->ngtc
== 0)
3006 /* Reading old version without tcoupl state data: set it */
3007 init_gtc_state(state
, ir
->opts
.ngtc
, 0, ir
->opts
.nhchainlength
);
3009 if (tpx
->bTop
&& mtop
)
3011 if (tpx
->fileVersion
< 57)
3013 if (mtop
->moltype
[0].ilist
[F_DISRES
].size() > 0)
3015 ir
->eDisre
= edrSimple
;
3019 ir
->eDisre
= edrNone
;
3027 * Process TPR data for file reading/writing.
3029 * The TPR file gets processed in in four stages due to the organization
3030 * of the data within it.
3032 * First, state data for the box is processed in do_tpx_state_first.
3033 * This is followed by processing the topology in do_tpx_mtop.
3034 * Coordinate and velocity vectors are handled next in do_tpx_state_second.
3035 * The last file information processed is the collection of simulation parameters in do_tpx_ir.
3036 * When reading, a final processing step is undertaken at the end.
3038 * \param[in] serializer Abstract serializer used to read/write data.
3039 * \param[in] tpx The file header data.
3040 * \param[in,out] ir Datastructures with simulation parameters.
3041 * \param[in,out] state Global state data.
3042 * \param[in,out] x Individual coordinates for processing, deprecated.
3043 * \param[in,out] v Individual velocities for processing, deprecated.
3044 * \param[in,out] mtop Global topology.
3046 static int do_tpx_body(gmx::ISerializer
* serializer
,
3056 do_tpx_state_first(serializer
, tpx
, state
);
3058 do_tpx_mtop(serializer
, tpx
, mtop
);
3061 do_tpx_state_second(serializer
, tpx
, state
, x
, v
);
3063 int ePBC
= do_tpx_ir(serializer
, tpx
, ir
);
3064 if (serializer
->reading())
3066 do_tpx_finalize(tpx
, ir
, state
, mtop
);
3072 * Overload for do_tpx_body that defaults to state vectors being nullptr.
3074 * \param[in] serializer Abstract serializer used to read/write data.
3075 * \param[in] tpx The file header data.
3076 * \param[in,out] ir Datastructures with simulation parameters.
3077 * \param[in,out] mtop Global topology.
3079 static int do_tpx_body(gmx::ISerializer
* serializer
, TpxFileHeader
* tpx
, t_inputrec
* ir
, gmx_mtop_t
* mtop
)
3081 return do_tpx_body(serializer
, tpx
, ir
, nullptr, nullptr, nullptr, mtop
);
3084 static t_fileio
* open_tpx(const char* fn
, const char* mode
)
3086 return gmx_fio_open(fn
, mode
);
3089 static void close_tpx(t_fileio
* fio
)
3095 * Fill information into the header only from state before writing.
3097 * Populating the header needs to be independent from writing the information
3098 * to file to allow things like writing the raw byte stream.
3100 * \param[in] state The current simulation state. Can't write without it.
3101 * \param[in] ir Parameter and system information.
3102 * \param[in] mtop Global topology.
3103 * \returns Fully populated header.
3105 static TpxFileHeader
populateTpxHeader(const t_state
& state
, const t_inputrec
* ir
, const gmx_mtop_t
* mtop
)
3107 TpxFileHeader header
;
3108 header
.natoms
= state
.natoms
;
3109 header
.ngtc
= state
.ngtc
;
3110 header
.fep_state
= state
.fep_state
;
3111 header
.lambda
= state
.lambda
[efptFEP
];
3112 header
.bIr
= ir
!= nullptr;
3113 header
.bTop
= mtop
!= nullptr;
3114 header
.bX
= (state
.flags
& (1 << estX
)) != 0;
3115 header
.bV
= (state
.flags
& (1 << estV
)) != 0;
3118 header
.fileVersion
= tpx_version
;
3119 header
.fileGeneration
= tpx_generation
;
3120 header
.isDouble
= (sizeof(real
) == sizeof(double));
3126 * Process the body of a TPR file as an opaque data buffer.
3128 * Reads/writes the information in \p buffer from/to the \p serializer
3129 * provided to the function. Does not interact with the actual
3130 * TPR datastructures but with an in memory representation of the
3131 * data, so that this data can be efficiently read or written from/to
3132 * an original source.
3134 * \param[in] serializer The abstract serializer used for reading or writing
3135 * the information in \p buffer.
3136 * \param[in,out] buffer Information from TPR file as char buffer.
3138 static void doTpxBodyBuffer(gmx::ISerializer
* serializer
, gmx::ArrayRef
<char> buffer
)
3140 serializer
->doOpaque(buffer
.data(), buffer
.size());
3144 * Populates simulation datastructures.
3146 * Here the information from the serialization interface \p serializer
3147 * is used to first populate the datastructures containing the simulation
3148 * information. Depending on the version found in the header \p tpx,
3149 * this is done using the new reading of the data as one block from disk,
3150 * followed by complete deserialization of the information read from there.
3151 * Otherwise, the datastructures are populated as before one by one from disk.
3152 * The second version is the default for the legacy tools that read the
3153 * coordinates and velocities separate from the state.
3155 * After reading in the data, a separate buffer is populated from them
3156 * containing only \p ir and \p mtop that can be communicated directly
3157 * to nodes needing the information to set up a simulation.
3159 * \param[in] tpx The file header.
3160 * \param[in] serializer The Serialization interface used to read the TPR.
3161 * \param[out] ir Input rec to populate.
3162 * \param[out] state State vectors to populate.
3163 * \param[out] x Coordinates to populate if needed.
3164 * \param[out] v Velocities to populate if needed.
3165 * \param[out] mtop Global topology to populate.
3167 * \returns Partial de-serialized TPR used for communication to nodes.
3169 static PartialDeserializedTprFile
readTpxBody(TpxFileHeader
* tpx
,
3170 gmx::ISerializer
* serializer
,
3177 PartialDeserializedTprFile partialDeserializedTpr
;
3178 if (tpx
->fileVersion
>= tpxv_AddSizeField
&& tpx
->fileGeneration
>= 27)
3180 partialDeserializedTpr
.body
.resize(tpx
->sizeOfTprBody
);
3181 partialDeserializedTpr
.header
= *tpx
;
3182 doTpxBodyBuffer(serializer
, partialDeserializedTpr
.body
);
3184 partialDeserializedTpr
.ePBC
=
3185 completeTprDeserialization(&partialDeserializedTpr
, ir
, state
, x
, v
, mtop
);
3189 partialDeserializedTpr
.ePBC
= do_tpx_body(serializer
, tpx
, ir
, state
, x
, v
, mtop
);
3191 // Update header to system info for communication to nodes.
3192 // As we only need to communicate the inputrec and mtop to other nodes,
3193 // we prepare a new char buffer with the information we have already read
3195 partialDeserializedTpr
.header
= populateTpxHeader(*state
, ir
, mtop
);
3196 // Long-term we should move to use little endian in files to avoid extra byte swapping,
3197 // but since we just used the default XDR format (which is big endian) for the TPR
3198 // header it would cause third-party libraries reading our raw data to tear their hair
3199 // if we swap the endian in the middle of the file, so we stick to big endian in the
3200 // TPR file for now - and thus we ask the serializer to swap if this host is little endian.
3201 gmx::InMemorySerializer
tprBodySerializer(gmx::EndianSwapBehavior::SwapIfHostIsLittleEndian
);
3202 do_tpx_body(&tprBodySerializer
, &partialDeserializedTpr
.header
, ir
, mtop
);
3203 partialDeserializedTpr
.body
= tprBodySerializer
.finishAndGetBuffer();
3205 return partialDeserializedTpr
;
3208 /************************************************************
3210 * The following routines are the exported ones
3212 ************************************************************/
3214 TpxFileHeader
readTpxHeader(const char* fileName
, bool canReadTopologyOnly
)
3218 fio
= open_tpx(fileName
, "r");
3219 gmx::FileIOXdrSerializer
serializer(fio
);
3222 do_tpxheader(&serializer
, &tpx
, fileName
, fio
, canReadTopologyOnly
);
3227 void write_tpx_state(const char* fn
, const t_inputrec
* ir
, const t_state
* state
, const gmx_mtop_t
* mtop
)
3229 /* To write a state, we first need to write the state information to a buffer before
3230 * we append the raw bytes to the file. For this, the header information needs to be
3231 * populated before we write the main body because it has some information that is
3232 * otherwise not available.
3237 TpxFileHeader tpx
= populateTpxHeader(*state
, ir
, mtop
);
3238 // Long-term we should move to use little endian in files to avoid extra byte swapping,
3239 // but since we just used the default XDR format (which is big endian) for the TPR
3240 // header it would cause third-party libraries reading our raw data to tear their hair
3241 // if we swap the endian in the middle of the file, so we stick to big endian in the
3242 // TPR file for now - and thus we ask the serializer to swap if this host is little endian.
3243 gmx::InMemorySerializer
tprBodySerializer(gmx::EndianSwapBehavior::SwapIfHostIsLittleEndian
);
3245 do_tpx_body(&tprBodySerializer
, &tpx
, const_cast<t_inputrec
*>(ir
), const_cast<t_state
*>(state
),
3246 nullptr, nullptr, const_cast<gmx_mtop_t
*>(mtop
));
3248 std::vector
<char> tprBody
= tprBodySerializer
.finishAndGetBuffer();
3249 tpx
.sizeOfTprBody
= tprBody
.size();
3251 fio
= open_tpx(fn
, "w");
3252 gmx::FileIOXdrSerializer
serializer(fio
);
3253 do_tpxheader(&serializer
, &tpx
, fn
, fio
, ir
== nullptr);
3254 doTpxBodyBuffer(&serializer
, tprBody
);
3259 int completeTprDeserialization(PartialDeserializedTprFile
* partialDeserializedTpr
,
3266 // Long-term we should move to use little endian in files to avoid extra byte swapping,
3267 // but since we just used the default XDR format (which is big endian) for the TPR
3268 // header it would cause third-party libraries reading our raw data to tear their hair
3269 // if we swap the endian in the middle of the file, so we stick to big endian in the
3270 // TPR file for now - and thus we ask the serializer to swap if this host is little endian.
3271 gmx::InMemoryDeserializer
tprBodyDeserializer(partialDeserializedTpr
->body
,
3272 partialDeserializedTpr
->header
.isDouble
,
3273 gmx::EndianSwapBehavior::SwapIfHostIsLittleEndian
);
3274 return do_tpx_body(&tprBodyDeserializer
, &partialDeserializedTpr
->header
, ir
, state
, x
, v
, mtop
);
3277 int completeTprDeserialization(PartialDeserializedTprFile
* partialDeserializedTpr
,
3281 return completeTprDeserialization(partialDeserializedTpr
, ir
, nullptr, nullptr, nullptr, mtop
);
3284 PartialDeserializedTprFile
read_tpx_state(const char* fn
, t_inputrec
* ir
, t_state
* state
, gmx_mtop_t
* mtop
)
3287 fio
= open_tpx(fn
, "r");
3288 gmx::FileIOXdrSerializer
serializer(fio
);
3289 PartialDeserializedTprFile partialDeserializedTpr
;
3290 do_tpxheader(&serializer
, &partialDeserializedTpr
.header
, fn
, fio
, ir
== nullptr);
3291 partialDeserializedTpr
=
3292 readTpxBody(&partialDeserializedTpr
.header
, &serializer
, ir
, state
, nullptr, nullptr, mtop
);
3294 return partialDeserializedTpr
;
3297 int read_tpx(const char* fn
, t_inputrec
* ir
, matrix box
, int* natoms
, rvec
* x
, rvec
* v
, gmx_mtop_t
* mtop
)
3303 fio
= open_tpx(fn
, "r");
3304 gmx::FileIOXdrSerializer
serializer(fio
);
3305 do_tpxheader(&serializer
, &tpx
, fn
, fio
, ir
== nullptr);
3306 PartialDeserializedTprFile partialDeserializedTpr
=
3307 readTpxBody(&tpx
, &serializer
, ir
, &state
, x
, v
, mtop
);
3309 if (mtop
!= nullptr && natoms
!= nullptr)
3311 *natoms
= mtop
->natoms
;
3315 copy_mat(state
.box
, box
);
3317 return partialDeserializedTpr
.ePBC
;
3320 int read_tpx_top(const char* fn
, t_inputrec
* ir
, matrix box
, int* natoms
, rvec
* x
, rvec
* v
, t_topology
* top
)
3325 ePBC
= read_tpx(fn
, ir
, box
, natoms
, x
, v
, &mtop
);
3327 *top
= gmx_mtop_t_to_t_topology(&mtop
, true);
3332 gmx_bool
fn2bTPX(const char* file
)
3334 return (efTPR
== fn2ftp(file
));
3337 void pr_tpxheader(FILE* fp
, int indent
, const char* title
, const TpxFileHeader
* sh
)
3339 if (available(fp
, sh
, indent
, title
))
3341 indent
= pr_title(fp
, indent
, title
);
3342 pr_indent(fp
, indent
);
3343 fprintf(fp
, "bIr = %spresent\n", sh
->bIr
? "" : "not ");
3344 pr_indent(fp
, indent
);
3345 fprintf(fp
, "bBox = %spresent\n", sh
->bBox
? "" : "not ");
3346 pr_indent(fp
, indent
);
3347 fprintf(fp
, "bTop = %spresent\n", sh
->bTop
? "" : "not ");
3348 pr_indent(fp
, indent
);
3349 fprintf(fp
, "bX = %spresent\n", sh
->bX
? "" : "not ");
3350 pr_indent(fp
, indent
);
3351 fprintf(fp
, "bV = %spresent\n", sh
->bV
? "" : "not ");
3352 pr_indent(fp
, indent
);
3353 fprintf(fp
, "bF = %spresent\n", sh
->bF
? "" : "not ");
3355 pr_indent(fp
, indent
);
3356 fprintf(fp
, "natoms = %d\n", sh
->natoms
);
3357 pr_indent(fp
, indent
);
3358 fprintf(fp
, "lambda = %e\n", sh
->lambda
);
3359 pr_indent(fp
, indent
);
3360 fprintf(fp
, "buffer size = %" PRId64
"\n", sh
->sizeOfTprBody
);