Move natoms_mol out of gmx_molblock_t
[gromacs.git] / src / gromacs / fileio / tpxio.cpp
blob89c90e38a070dd1eedcca99fbd3b05d5773348c0
1 /*
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, 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.
37 #include "gmxpre.h"
39 /* This file is completely threadsafe - keep it that way! */
41 #include "tpxio.h"
43 #include <cstdio>
44 #include <cstdlib>
45 #include <cstring>
47 #include <algorithm>
48 #include <vector>
50 #include "gromacs/fileio/filetypes.h"
51 #include "gromacs/fileio/gmxfio.h"
52 #include "gromacs/fileio/gmxfio-xdr.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdtypes/awh-history.h"
56 #include "gromacs/mdtypes/awh-params.h"
57 #include "gromacs/mdtypes/inputrec.h"
58 #include "gromacs/mdtypes/md_enums.h"
59 #include "gromacs/mdtypes/pull-params.h"
60 #include "gromacs/mdtypes/state.h"
61 #include "gromacs/pbcutil/boxutilities.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/topology/block.h"
64 #include "gromacs/topology/ifunc.h"
65 #include "gromacs/topology/mtop_util.h"
66 #include "gromacs/topology/symtab.h"
67 #include "gromacs/topology/topology.h"
68 #include "gromacs/utility/arraysize.h"
69 #include "gromacs/utility/baseversion.h"
70 #include "gromacs/utility/cstringutil.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/futil.h"
73 #include "gromacs/utility/gmxassert.h"
74 #include "gromacs/utility/keyvaluetreebuilder.h"
75 #include "gromacs/utility/keyvaluetreeserializer.h"
76 #include "gromacs/utility/smalloc.h"
77 #include "gromacs/utility/snprintf.h"
78 #include "gromacs/utility/txtdump.h"
80 #define TPX_TAG_RELEASE "release"
82 /*! \brief Tag string for the file format written to run input files
83 * written by this version of the code.
85 * Change this if you want to change the run input format in a feature
86 * branch. This ensures that there will not be different run input
87 * formats around which cannot be distinguished, while not causing
88 * problems rebasing the feature branch onto upstream changes. When
89 * merging with mainstream GROMACS, set this tag string back to
90 * TPX_TAG_RELEASE, and instead add an element to tpxv.
92 static const char *tpx_tag = TPX_TAG_RELEASE;
94 /*! \brief Enum of values that describe the contents of a tpr file
95 * whose format matches a version number
97 * The enum helps the code be more self-documenting and ensure merges
98 * do not silently resolve when two patches make the same bump. When
99 * adding new functionality, add a new element just above tpxv_Count
100 * in this enumeration, and write code below that does the right thing
101 * according to the value of file_version.
103 enum tpxv {
104 tpxv_ComputationalElectrophysiology = 96, /**< support for ion/water position swaps (computational electrophysiology) */
105 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to gmx_int64_t */
106 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
107 tpxv_InteractiveMolecularDynamics, /**< interactive molecular dynamics (IMD) */
108 tpxv_RemoveObsoleteParameters1, /**< remove optimize_fft, dihre_fc, nstcheckpoint */
109 tpxv_PullCoordTypeGeom, /**< add pull type and geometry per group and flat-bottom */
110 tpxv_PullGeomDirRel, /**< add pull geometry direction-relative */
111 tpxv_IntermolecularBondeds, /**< permit inter-molecular bonded interactions in the topology */
112 tpxv_CompElWithSwapLayerOffset, /**< added parameters for improved CompEl setups */
113 tpxv_CompElPolyatomicIonsAndMultipleIonTypes, /**< CompEl now can handle polyatomic ions and more than two types of ions */
114 tpxv_RemoveAdress, /**< removed support for AdResS */
115 tpxv_PullCoordNGroup, /**< add ngroup to pull coord */
116 tpxv_RemoveTwinRange, /**< removed support for twin-range interactions */
117 tpxv_ReplacePullPrintCOM12, /**< Replaced print-com-1, 2 with pull-print-com */
118 tpxv_PullExternalPotential, /**< Added pull type external potential */
119 tpxv_GenericParamsForElectricField, /**< Introduced KeyValueTree and moved electric field parameters */
120 tpxv_AcceleratedWeightHistogram, /**< sampling with accelerated weight histogram method (AWH) */
121 tpxv_RemoveImplicitSolvation, /**< removed support for implicit solvation */
122 tpxv_Count /**< the total number of tpxv versions */
125 /*! \brief Version number of the file format written to run input
126 * files by this version of the code.
128 * The tpx_version increases whenever the file format in the main
129 * development branch changes, due to an extension of the tpxv enum above.
130 * Backward compatibility for reading old run input files is maintained
131 * by checking this version number against that of the file and then using
132 * the correct code path.
134 * When developing a feature branch that needs to change the run input
135 * file format, change tpx_tag instead. */
136 static const int tpx_version = tpxv_Count - 1;
139 /* This number should only be increased when you edit the TOPOLOGY section
140 * or the HEADER of the tpx format.
141 * This way we can maintain forward compatibility too for all analysis tools
142 * and/or external programs that only need to know the atom/residue names,
143 * charges, and bond connectivity.
145 * It first appeared in tpx version 26, when I also moved the inputrecord
146 * to the end of the tpx file, so we can just skip it if we only
147 * want the topology.
149 * In particular, it must be increased when adding new elements to
150 * ftupd, so that old code can read new .tpr files.
152 static const int tpx_generation = 26;
154 /* This number should be the most recent backwards incompatible version
155 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
157 static const int tpx_incompatible_version = 57; // GMX4.0 has version 58
161 /* Struct used to maintain tpx compatibility when function types are added */
162 typedef struct {
163 int fvnr; /* file version number in which the function type first appeared */
164 int ftype; /* function type */
165 } t_ftupd;
168 * TODO The following three lines make little sense, please clarify if
169 * you've had to work out how ftupd works.
171 * The entries should be ordered in:
172 * 1. ascending function type number
173 * 2. ascending file version number
175 * Because we support reading of old .tpr file versions (even when
176 * mdrun can no longer run the simulation), we need to be able to read
177 * obsolete t_interaction_function types. Any data read from such
178 * fields is discarded. Their names have _NOLONGERUSED appended to
179 * them to make things clear.
181 static const t_ftupd ftupd[] = {
182 { 70, F_RESTRBONDS },
183 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
184 { 76, F_LINEAR_ANGLES },
185 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
186 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
187 { 65, F_CMAP },
188 { 60, F_GB12_NOLONGERUSED },
189 { 61, F_GB13_NOLONGERUSED },
190 { 61, F_GB14_NOLONGERUSED },
191 { 72, F_GBPOL_NOLONGERUSED },
192 { 72, F_NPSOLVATION_NOLONGERUSED },
193 { 93, F_LJ_RECIP },
194 { 90, F_FBPOSRES },
195 { 69, F_VTEMP_NOLONGERUSED},
196 { 66, F_PDISPCORR },
197 { 76, F_ANHARM_POL },
198 { 79, F_DVDL_COUL },
199 { 79, F_DVDL_VDW, },
200 { 79, F_DVDL_BONDED, },
201 { 79, F_DVDL_RESTRAINT },
202 { 79, F_DVDL_TEMPERATURE },
204 #define NFTUPD asize(ftupd)
206 /* Needed for backward compatibility */
207 #define MAXNODES 256
209 /**************************************************************
211 * Now the higer level routines that do io of the structures and arrays
213 **************************************************************/
214 static void do_pullgrp_tpx_pre95(t_fileio *fio,
215 t_pull_group *pgrp,
216 t_pull_coord *pcrd,
217 gmx_bool bRead)
219 rvec tmp;
221 gmx_fio_do_int(fio, pgrp->nat);
222 if (bRead)
224 snew(pgrp->ind, pgrp->nat);
226 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
227 gmx_fio_do_int(fio, pgrp->nweight);
228 if (bRead)
230 snew(pgrp->weight, pgrp->nweight);
232 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
233 gmx_fio_do_int(fio, pgrp->pbcatom);
234 gmx_fio_do_rvec(fio, pcrd->vec);
235 clear_rvec(pcrd->origin);
236 gmx_fio_do_rvec(fio, tmp);
237 pcrd->init = tmp[0];
238 gmx_fio_do_real(fio, pcrd->rate);
239 gmx_fio_do_real(fio, pcrd->k);
240 gmx_fio_do_real(fio, pcrd->kB);
243 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
245 gmx_fio_do_int(fio, pgrp->nat);
246 if (bRead)
248 snew(pgrp->ind, pgrp->nat);
250 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
251 gmx_fio_do_int(fio, pgrp->nweight);
252 if (bRead)
254 snew(pgrp->weight, pgrp->nweight);
256 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
257 gmx_fio_do_int(fio, pgrp->pbcatom);
260 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd,
261 gmx_bool bRead, int file_version,
262 int ePullOld, int eGeomOld, ivec dimOld)
264 if (file_version >= tpxv_PullCoordNGroup)
266 gmx_fio_do_int(fio, pcrd->eType);
267 if (file_version >= tpxv_PullExternalPotential)
269 if (pcrd->eType == epullEXTERNAL)
271 if (bRead)
273 char buf[STRLEN];
275 gmx_fio_do_string(fio, buf);
276 pcrd->externalPotentialProvider = gmx_strdup(buf);
278 else
280 gmx_fio_do_string(fio, pcrd->externalPotentialProvider);
283 else
285 pcrd->externalPotentialProvider = nullptr;
288 else
290 if (bRead)
292 pcrd->externalPotentialProvider = nullptr;
295 /* Note that we try to support adding new geometries without
296 * changing the tpx version. This requires checks when printing the
297 * geometry string and a check and fatal_error in init_pull.
299 gmx_fio_do_int(fio, pcrd->eGeom);
300 gmx_fio_do_int(fio, pcrd->ngroup);
301 if (pcrd->ngroup <= c_pullCoordNgroupMax)
303 gmx_fio_ndo_int(fio, pcrd->group, pcrd->ngroup);
305 else
307 /* More groups in file than supported, this must be a new geometry
308 * that is not supported by our current code. Since we will not
309 * use the groups for this coord (checks in the pull and WHAM code
310 * ensure this), we can ignore the groups and set ngroup=0.
312 int *dum;
313 snew(dum, pcrd->ngroup);
314 gmx_fio_ndo_int(fio, dum, pcrd->ngroup);
315 sfree(dum);
317 pcrd->ngroup = 0;
319 gmx_fio_do_ivec(fio, pcrd->dim);
321 else
323 pcrd->ngroup = 2;
324 gmx_fio_do_int(fio, pcrd->group[0]);
325 gmx_fio_do_int(fio, pcrd->group[1]);
326 if (file_version >= tpxv_PullCoordTypeGeom)
328 pcrd->ngroup = (pcrd->eGeom == epullgDIRRELATIVE ? 4 : 2);
329 gmx_fio_do_int(fio, pcrd->eType);
330 gmx_fio_do_int(fio, pcrd->eGeom);
331 if (pcrd->ngroup == 4)
333 gmx_fio_do_int(fio, pcrd->group[2]);
334 gmx_fio_do_int(fio, pcrd->group[3]);
336 gmx_fio_do_ivec(fio, pcrd->dim);
338 else
340 pcrd->eType = ePullOld;
341 pcrd->eGeom = eGeomOld;
342 copy_ivec(dimOld, pcrd->dim);
345 gmx_fio_do_rvec(fio, pcrd->origin);
346 gmx_fio_do_rvec(fio, pcrd->vec);
347 if (file_version >= tpxv_PullCoordTypeGeom)
349 gmx_fio_do_gmx_bool(fio, pcrd->bStart);
351 else
353 /* This parameter is only printed, but not actually used by mdrun */
354 pcrd->bStart = FALSE;
356 gmx_fio_do_real(fio, pcrd->init);
357 gmx_fio_do_real(fio, pcrd->rate);
358 gmx_fio_do_real(fio, pcrd->k);
359 gmx_fio_do_real(fio, pcrd->kB);
362 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
364 int n_lambda = fepvals->n_lambda;
366 /* reset the lambda calculation window */
367 fepvals->lambda_start_n = 0;
368 fepvals->lambda_stop_n = n_lambda;
369 if (file_version >= 79)
371 if (n_lambda > 0)
373 if (bRead)
375 snew(expand->init_lambda_weights, n_lambda);
377 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
378 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
381 gmx_fio_do_int(fio, expand->nstexpanded);
382 gmx_fio_do_int(fio, expand->elmcmove);
383 gmx_fio_do_int(fio, expand->elamstats);
384 gmx_fio_do_int(fio, expand->lmc_repeats);
385 gmx_fio_do_int(fio, expand->gibbsdeltalam);
386 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
387 gmx_fio_do_int(fio, expand->lmc_seed);
388 gmx_fio_do_real(fio, expand->mc_temp);
389 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
390 gmx_fio_do_int(fio, expand->nstTij);
391 gmx_fio_do_int(fio, expand->minvarmin);
392 gmx_fio_do_int(fio, expand->c_range);
393 gmx_fio_do_real(fio, expand->wl_scale);
394 gmx_fio_do_real(fio, expand->wl_ratio);
395 gmx_fio_do_real(fio, expand->init_wl_delta);
396 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
397 gmx_fio_do_int(fio, expand->elmceq);
398 gmx_fio_do_int(fio, expand->equil_steps);
399 gmx_fio_do_int(fio, expand->equil_samples);
400 gmx_fio_do_int(fio, expand->equil_n_at_lam);
401 gmx_fio_do_real(fio, expand->equil_wl_delta);
402 gmx_fio_do_real(fio, expand->equil_ratio);
406 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
407 int file_version)
409 if (file_version >= 79)
411 gmx_fio_do_int(fio, simtemp->eSimTempScale);
412 gmx_fio_do_real(fio, simtemp->simtemp_high);
413 gmx_fio_do_real(fio, simtemp->simtemp_low);
414 if (n_lambda > 0)
416 if (bRead)
418 snew(simtemp->temperatures, n_lambda);
420 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
425 static void do_imd(t_fileio *fio, t_IMD *imd, gmx_bool bRead)
427 gmx_fio_do_int(fio, imd->nat);
428 if (bRead)
430 snew(imd->ind, imd->nat);
432 gmx_fio_ndo_int(fio, imd->ind, imd->nat);
435 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
437 /* i is defined in the ndo_double macro; use g to iterate. */
438 int g;
439 real rdum;
441 /* free energy values */
443 if (file_version >= 79)
445 gmx_fio_do_int(fio, fepvals->init_fep_state);
446 gmx_fio_do_double(fio, fepvals->init_lambda);
447 gmx_fio_do_double(fio, fepvals->delta_lambda);
449 else if (file_version >= 59)
451 gmx_fio_do_double(fio, fepvals->init_lambda);
452 gmx_fio_do_double(fio, fepvals->delta_lambda);
454 else
456 gmx_fio_do_real(fio, rdum);
457 fepvals->init_lambda = rdum;
458 gmx_fio_do_real(fio, rdum);
459 fepvals->delta_lambda = rdum;
461 if (file_version >= 79)
463 gmx_fio_do_int(fio, fepvals->n_lambda);
464 if (bRead)
466 snew(fepvals->all_lambda, efptNR);
468 for (g = 0; g < efptNR; g++)
470 if (fepvals->n_lambda > 0)
472 if (bRead)
474 snew(fepvals->all_lambda[g], fepvals->n_lambda);
476 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
477 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
479 else if (fepvals->init_lambda >= 0)
481 fepvals->separate_dvdl[efptFEP] = TRUE;
485 else if (file_version >= 64)
487 gmx_fio_do_int(fio, fepvals->n_lambda);
488 if (bRead)
490 int g;
492 snew(fepvals->all_lambda, efptNR);
493 /* still allocate the all_lambda array's contents. */
494 for (g = 0; g < efptNR; g++)
496 if (fepvals->n_lambda > 0)
498 snew(fepvals->all_lambda[g], fepvals->n_lambda);
502 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
503 fepvals->n_lambda);
504 if (fepvals->init_lambda >= 0)
506 int g, h;
508 fepvals->separate_dvdl[efptFEP] = TRUE;
510 if (bRead)
512 /* copy the contents of the efptFEP lambda component to all
513 the other components */
514 for (g = 0; g < efptNR; g++)
516 for (h = 0; h < fepvals->n_lambda; h++)
518 if (g != efptFEP)
520 fepvals->all_lambda[g][h] =
521 fepvals->all_lambda[efptFEP][h];
528 else
530 fepvals->n_lambda = 0;
531 fepvals->all_lambda = nullptr;
532 if (fepvals->init_lambda >= 0)
534 fepvals->separate_dvdl[efptFEP] = TRUE;
537 gmx_fio_do_real(fio, fepvals->sc_alpha);
538 gmx_fio_do_int(fio, fepvals->sc_power);
539 if (file_version >= 79)
541 gmx_fio_do_real(fio, fepvals->sc_r_power);
543 else
545 fepvals->sc_r_power = 6.0;
547 gmx_fio_do_real(fio, fepvals->sc_sigma);
548 if (bRead)
550 if (file_version >= 71)
552 fepvals->sc_sigma_min = fepvals->sc_sigma;
554 else
556 fepvals->sc_sigma_min = 0;
559 if (file_version >= 79)
561 gmx_fio_do_int(fio, fepvals->bScCoul);
563 else
565 fepvals->bScCoul = TRUE;
567 if (file_version >= 64)
569 gmx_fio_do_int(fio, fepvals->nstdhdl);
571 else
573 fepvals->nstdhdl = 1;
576 if (file_version >= 73)
578 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
579 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
581 else
583 fepvals->separate_dhdl_file = esepdhdlfileYES;
584 fepvals->dhdl_derivatives = edhdlderivativesYES;
586 if (file_version >= 71)
588 gmx_fio_do_int(fio, fepvals->dh_hist_size);
589 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
591 else
593 fepvals->dh_hist_size = 0;
594 fepvals->dh_hist_spacing = 0.1;
596 if (file_version >= 79)
598 gmx_fio_do_int(fio, fepvals->edHdLPrintEnergy);
600 else
602 fepvals->edHdLPrintEnergy = edHdLPrintEnergyNO;
605 /* handle lambda_neighbors */
606 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
608 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
609 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
610 (fepvals->init_lambda < 0) )
612 fepvals->lambda_start_n = (fepvals->init_fep_state -
613 fepvals->lambda_neighbors);
614 fepvals->lambda_stop_n = (fepvals->init_fep_state +
615 fepvals->lambda_neighbors + 1);
616 if (fepvals->lambda_start_n < 0)
618 fepvals->lambda_start_n = 0;;
620 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
622 fepvals->lambda_stop_n = fepvals->n_lambda;
625 else
627 fepvals->lambda_start_n = 0;
628 fepvals->lambda_stop_n = fepvals->n_lambda;
631 else
633 fepvals->lambda_start_n = 0;
634 fepvals->lambda_stop_n = fepvals->n_lambda;
638 static void do_awhBias(t_fileio *fio, gmx::AwhBiasParams *awhBiasParams, gmx_bool bRead,
639 int gmx_unused file_version)
641 gmx_fio_do_int(fio, awhBiasParams->eTarget);
642 gmx_fio_do_double(fio, awhBiasParams->targetBetaScaling);
643 gmx_fio_do_double(fio, awhBiasParams->targetCutoff);
644 gmx_fio_do_int(fio, awhBiasParams->eGrowth);
645 gmx_fio_do_int(fio, awhBiasParams->bUserData);
646 gmx_fio_do_double(fio, awhBiasParams->errorInitial);
647 gmx_fio_do_int(fio, awhBiasParams->ndim);
648 gmx_fio_do_int(fio, awhBiasParams->shareGroup);
649 gmx_fio_do_gmx_bool(fio, awhBiasParams->equilibrateHistogram);
651 if (bRead)
653 snew(awhBiasParams->dimParams, awhBiasParams->ndim);
656 for (int d = 0; d < awhBiasParams->ndim; d++)
658 gmx::AwhDimParams *dimParams = &awhBiasParams->dimParams[d];
660 gmx_fio_do_int(fio, dimParams->eCoordProvider);
661 gmx_fio_do_int(fio, dimParams->coordIndex);
662 gmx_fio_do_double(fio, dimParams->origin);
663 gmx_fio_do_double(fio, dimParams->end);
664 gmx_fio_do_double(fio, dimParams->period);
665 gmx_fio_do_double(fio, dimParams->forceConstant);
666 gmx_fio_do_double(fio, dimParams->diffusion);
667 gmx_fio_do_double(fio, dimParams->coordValueInit);
668 gmx_fio_do_double(fio, dimParams->coverDiameter);
672 static void do_awh(t_fileio *fio, gmx::AwhParams *awhParams, gmx_bool bRead,
673 int gmx_unused file_version)
675 gmx_fio_do_int(fio, awhParams->numBias);
676 gmx_fio_do_int(fio, awhParams->nstOut);
677 gmx_fio_do_int64(fio, awhParams->seed);
678 gmx_fio_do_int(fio, awhParams->nstSampleCoord);
679 gmx_fio_do_int(fio, awhParams->numSamplesUpdateFreeEnergy);
680 gmx_fio_do_int(fio, awhParams->ePotential);
681 gmx_fio_do_gmx_bool(fio, awhParams->shareBiasMultisim);
683 if (awhParams->numBias > 0)
685 if (bRead)
687 snew(awhParams->awhBiasParams, awhParams->numBias);
690 for (int k = 0; k < awhParams->numBias; k++)
692 do_awhBias(fio, &awhParams->awhBiasParams[k], bRead, file_version);
697 static void do_pull(t_fileio *fio, pull_params_t *pull, gmx_bool bRead,
698 int file_version, int ePullOld)
700 int eGeomOld = -1;
701 ivec dimOld;
702 int g;
704 if (file_version >= 95)
706 gmx_fio_do_int(fio, pull->ngroup);
708 gmx_fio_do_int(fio, pull->ncoord);
709 if (file_version < 95)
711 pull->ngroup = pull->ncoord + 1;
713 if (file_version < tpxv_PullCoordTypeGeom)
715 real dum;
717 gmx_fio_do_int(fio, eGeomOld);
718 gmx_fio_do_ivec(fio, dimOld);
719 /* The inner cylinder radius, now removed */
720 gmx_fio_do_real(fio, dum);
722 gmx_fio_do_real(fio, pull->cylinder_r);
723 gmx_fio_do_real(fio, pull->constr_tol);
724 if (file_version >= 95)
726 gmx_fio_do_int(fio, pull->bPrintCOM);
727 /* With file_version < 95 this value is set below */
729 if (file_version >= tpxv_ReplacePullPrintCOM12)
731 gmx_fio_do_int(fio, pull->bPrintRefValue);
732 gmx_fio_do_int(fio, pull->bPrintComp);
734 else if (file_version >= tpxv_PullCoordTypeGeom)
736 int idum;
737 gmx_fio_do_int(fio, idum); /* used to be bPrintCOM2 */
738 gmx_fio_do_int(fio, pull->bPrintRefValue);
739 gmx_fio_do_int(fio, pull->bPrintComp);
741 else
743 pull->bPrintRefValue = FALSE;
744 pull->bPrintComp = TRUE;
746 gmx_fio_do_int(fio, pull->nstxout);
747 gmx_fio_do_int(fio, pull->nstfout);
748 if (bRead)
750 snew(pull->group, pull->ngroup);
751 snew(pull->coord, pull->ncoord);
753 if (file_version < 95)
755 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
756 if (eGeomOld == epullgDIRPBC)
758 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
760 if (eGeomOld > epullgDIRPBC)
762 eGeomOld -= 1;
765 for (g = 0; g < pull->ngroup; g++)
767 /* We read and ignore a pull coordinate for group 0 */
768 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[std::max(g-1, 0)],
769 bRead);
770 if (g > 0)
772 pull->coord[g-1].group[0] = 0;
773 pull->coord[g-1].group[1] = g;
777 pull->bPrintCOM = (pull->group[0].nat > 0);
779 else
781 for (g = 0; g < pull->ngroup; g++)
783 do_pull_group(fio, &pull->group[g], bRead);
785 for (g = 0; g < pull->ncoord; g++)
787 do_pull_coord(fio, &pull->coord[g],
788 bRead, file_version, ePullOld, eGeomOld, dimOld);
794 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
796 gmx_fio_do_int(fio, rotg->eType);
797 gmx_fio_do_int(fio, rotg->bMassW);
798 gmx_fio_do_int(fio, rotg->nat);
799 if (bRead)
801 snew(rotg->ind, rotg->nat);
803 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
804 if (bRead)
806 snew(rotg->x_ref, rotg->nat);
808 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
809 gmx_fio_do_rvec(fio, rotg->vec);
810 gmx_fio_do_rvec(fio, rotg->pivot);
811 gmx_fio_do_real(fio, rotg->rate);
812 gmx_fio_do_real(fio, rotg->k);
813 gmx_fio_do_real(fio, rotg->slab_dist);
814 gmx_fio_do_real(fio, rotg->min_gaussian);
815 gmx_fio_do_real(fio, rotg->eps);
816 gmx_fio_do_int(fio, rotg->eFittype);
817 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
818 gmx_fio_do_real(fio, rotg->PotAngle_step);
821 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
823 int g;
825 gmx_fio_do_int(fio, rot->ngrp);
826 gmx_fio_do_int(fio, rot->nstrout);
827 gmx_fio_do_int(fio, rot->nstsout);
828 if (bRead)
830 snew(rot->grp, rot->ngrp);
832 for (g = 0; g < rot->ngrp; g++)
834 do_rotgrp(fio, &rot->grp[g], bRead);
839 static void do_swapgroup(t_fileio *fio, t_swapGroup *g, gmx_bool bRead)
842 /* Name of the group or molecule */
843 if (bRead)
845 char buf[STRLEN];
847 gmx_fio_do_string(fio, buf);
848 g->molname = gmx_strdup(buf);
850 else
852 gmx_fio_do_string(fio, g->molname);
855 /* Number of atoms in the group */
856 gmx_fio_do_int(fio, g->nat);
858 /* The group's atom indices */
859 if (bRead)
861 snew(g->ind, g->nat);
863 gmx_fio_ndo_int(fio, g->ind, g->nat);
865 /* Requested counts for compartments A and B */
866 gmx_fio_ndo_int(fio, g->nmolReq, eCompNR);
869 static void do_swapcoords_tpx(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead, int file_version)
871 /* Enums for better readability of the code */
872 enum {
873 eCompA = 0, eCompB
875 enum {
876 eChannel0 = 0, eChannel1
880 if (file_version >= tpxv_CompElPolyatomicIonsAndMultipleIonTypes)
882 /* The total number of swap groups is the sum of the fixed groups
883 * (split0, split1, solvent), and the user-defined groups (2+ types of ions)
885 gmx_fio_do_int(fio, swap->ngrp);
886 if (bRead)
888 snew(swap->grp, swap->ngrp);
890 for (int ig = 0; ig < swap->ngrp; ig++)
892 do_swapgroup(fio, &swap->grp[ig], bRead);
894 gmx_fio_do_int(fio, swap->massw_split[eChannel0]);
895 gmx_fio_do_int(fio, swap->massw_split[eChannel1]);
896 gmx_fio_do_int(fio, swap->nstswap);
897 gmx_fio_do_int(fio, swap->nAverage);
898 gmx_fio_do_real(fio, swap->threshold);
899 gmx_fio_do_real(fio, swap->cyl0r);
900 gmx_fio_do_real(fio, swap->cyl0u);
901 gmx_fio_do_real(fio, swap->cyl0l);
902 gmx_fio_do_real(fio, swap->cyl1r);
903 gmx_fio_do_real(fio, swap->cyl1u);
904 gmx_fio_do_real(fio, swap->cyl1l);
906 else
908 /*** Support reading older CompEl .tpr files ***/
910 /* In the original CompEl .tpr files, we always have 5 groups: */
911 swap->ngrp = 5;
912 snew(swap->grp, swap->ngrp);
914 swap->grp[eGrpSplit0 ].molname = gmx_strdup("split0" ); // group 0: split0
915 swap->grp[eGrpSplit1 ].molname = gmx_strdup("split1" ); // group 1: split1
916 swap->grp[eGrpSolvent].molname = gmx_strdup("solvent"); // group 2: solvent
917 swap->grp[3 ].molname = gmx_strdup("anions" ); // group 3: anions
918 swap->grp[4 ].molname = gmx_strdup("cations"); // group 4: cations
920 gmx_fio_do_int(fio, swap->grp[3].nat);
921 gmx_fio_do_int(fio, swap->grp[eGrpSolvent].nat);
922 gmx_fio_do_int(fio, swap->grp[eGrpSplit0].nat);
923 gmx_fio_do_int(fio, swap->massw_split[eChannel0]);
924 gmx_fio_do_int(fio, swap->grp[eGrpSplit1].nat);
925 gmx_fio_do_int(fio, swap->massw_split[eChannel1]);
926 gmx_fio_do_int(fio, swap->nstswap);
927 gmx_fio_do_int(fio, swap->nAverage);
928 gmx_fio_do_real(fio, swap->threshold);
929 gmx_fio_do_real(fio, swap->cyl0r);
930 gmx_fio_do_real(fio, swap->cyl0u);
931 gmx_fio_do_real(fio, swap->cyl0l);
932 gmx_fio_do_real(fio, swap->cyl1r);
933 gmx_fio_do_real(fio, swap->cyl1u);
934 gmx_fio_do_real(fio, swap->cyl1l);
936 // The order[] array keeps compatibility with older .tpr files
937 // by reading in the groups in the classic order
939 const int order[4] = {3, eGrpSolvent, eGrpSplit0, eGrpSplit1 };
941 for (int ig = 0; ig < 4; ig++)
943 int g = order[ig];
944 snew(swap->grp[g].ind, swap->grp[g].nat);
945 gmx_fio_ndo_int(fio, swap->grp[g].ind, swap->grp[g].nat);
949 for (int j = eCompA; j <= eCompB; j++)
951 gmx_fio_do_int(fio, swap->grp[3].nmolReq[j]); // group 3 = anions
952 gmx_fio_do_int(fio, swap->grp[4].nmolReq[j]); // group 4 = cations
954 } /* End support reading older CompEl .tpr files */
956 if (file_version >= tpxv_CompElWithSwapLayerOffset)
958 gmx_fio_do_real(fio, swap->bulkOffset[eCompA]);
959 gmx_fio_do_real(fio, swap->bulkOffset[eCompB]);
964 static void do_legacy_efield(t_fileio *fio, gmx::KeyValueTreeObjectBuilder *root)
966 const char *const dimName[] = { "x", "y", "z" };
968 auto appliedForcesObj = root->addObject("applied-forces");
969 auto efieldObj = appliedForcesObj.addObject("electric-field");
970 // The content of the tpr file for this feature has
971 // been the same since gromacs 4.0 that was used for
972 // developing.
973 for (int j = 0; j < DIM; ++j)
975 int n, nt;
976 gmx_fio_do_int(fio, n);
977 gmx_fio_do_int(fio, nt);
978 std::vector<real> aa(n+1), phi(nt+1), at(nt+1), phit(nt+1);
979 gmx_fio_ndo_real(fio, aa.data(), n);
980 gmx_fio_ndo_real(fio, phi.data(), n);
981 gmx_fio_ndo_real(fio, at.data(), nt);
982 gmx_fio_ndo_real(fio, phit.data(), nt);
983 if (n > 0)
985 if (n > 1 || nt > 1)
987 gmx_fatal(FARGS, "Can not handle tpr files with more than one electric field term per direction.");
989 auto dimObj = efieldObj.addObject(dimName[j]);
990 dimObj.addValue<real>("E0", aa[0]);
991 dimObj.addValue<real>("omega", at[0]);
992 dimObj.addValue<real>("t0", phi[0]);
993 dimObj.addValue<real>("sigma", phit[0]);
999 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
1000 int file_version)
1002 int i, j, k, idum = 0;
1003 real rdum;
1004 gmx_bool bdum = 0;
1006 if (file_version != tpx_version)
1008 /* Give a warning about features that are not accessible */
1009 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
1010 file_version, tpx_version);
1013 if (file_version == 0)
1015 return;
1018 gmx::KeyValueTreeBuilder paramsBuilder;
1019 gmx::KeyValueTreeObjectBuilder paramsObj = paramsBuilder.rootObject();
1021 /* Basic inputrec stuff */
1022 gmx_fio_do_int(fio, ir->eI);
1023 if (file_version >= 62)
1025 gmx_fio_do_int64(fio, ir->nsteps);
1027 else
1029 gmx_fio_do_int(fio, idum);
1030 ir->nsteps = idum;
1033 if (file_version >= 62)
1035 gmx_fio_do_int64(fio, ir->init_step);
1037 else
1039 gmx_fio_do_int(fio, idum);
1040 ir->init_step = idum;
1043 gmx_fio_do_int(fio, ir->simulation_part);
1045 if (file_version >= 67)
1047 gmx_fio_do_int(fio, ir->nstcalcenergy);
1049 else
1051 ir->nstcalcenergy = 1;
1053 if (file_version >= 81)
1055 gmx_fio_do_int(fio, ir->cutoff_scheme);
1056 if (file_version < 94)
1058 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
1061 else
1063 ir->cutoff_scheme = ecutsGROUP;
1065 gmx_fio_do_int(fio, ir->ns_type);
1066 gmx_fio_do_int(fio, ir->nstlist);
1067 gmx_fio_do_int(fio, idum); /* used to be ndelta; not used anymore */
1069 gmx_fio_do_real(fio, ir->rtpi);
1071 gmx_fio_do_int(fio, ir->nstcomm);
1072 gmx_fio_do_int(fio, ir->comm_mode);
1074 /* ignore nstcheckpoint */
1075 if (file_version < tpxv_RemoveObsoleteParameters1)
1077 gmx_fio_do_int(fio, idum);
1080 gmx_fio_do_int(fio, ir->nstcgsteep);
1082 gmx_fio_do_int(fio, ir->nbfgscorr);
1084 gmx_fio_do_int(fio, ir->nstlog);
1085 gmx_fio_do_int(fio, ir->nstxout);
1086 gmx_fio_do_int(fio, ir->nstvout);
1087 gmx_fio_do_int(fio, ir->nstfout);
1088 gmx_fio_do_int(fio, ir->nstenergy);
1089 gmx_fio_do_int(fio, ir->nstxout_compressed);
1090 if (file_version >= 59)
1092 gmx_fio_do_double(fio, ir->init_t);
1093 gmx_fio_do_double(fio, ir->delta_t);
1095 else
1097 gmx_fio_do_real(fio, rdum);
1098 ir->init_t = rdum;
1099 gmx_fio_do_real(fio, rdum);
1100 ir->delta_t = rdum;
1102 gmx_fio_do_real(fio, ir->x_compression_precision);
1103 if (file_version >= 81)
1105 gmx_fio_do_real(fio, ir->verletbuf_tol);
1107 else
1109 ir->verletbuf_tol = 0;
1111 gmx_fio_do_real(fio, ir->rlist);
1112 if (file_version >= 67 && file_version < tpxv_RemoveTwinRange)
1114 if (bRead)
1116 // Reading such a file version could invoke the twin-range
1117 // scheme, about which mdrun should give a fatal error.
1118 real dummy_rlistlong = -1;
1119 gmx_fio_do_real(fio, dummy_rlistlong);
1121 if (ir->rlist > 0 && (dummy_rlistlong == 0 || dummy_rlistlong > ir->rlist))
1123 // Get mdrun to issue an error (regardless of
1124 // ir->cutoff_scheme).
1125 ir->useTwinRange = true;
1127 else
1129 // grompp used to set rlistlong actively. Users were
1130 // probably also confused and set rlistlong == rlist.
1131 // However, in all remaining cases, it is safe to let
1132 // mdrun proceed normally.
1133 ir->useTwinRange = false;
1137 else
1139 // No need to read or write anything
1140 ir->useTwinRange = false;
1142 if (file_version >= 82 && file_version != 90)
1144 // Multiple time-stepping is no longer enabled, but the old
1145 // support required the twin-range scheme, for which mdrun
1146 // already emits a fatal error.
1147 int dummy_nstcalclr = -1;
1148 gmx_fio_do_int(fio, dummy_nstcalclr);
1150 gmx_fio_do_int(fio, ir->coulombtype);
1151 if (file_version >= 81)
1153 gmx_fio_do_int(fio, ir->coulomb_modifier);
1155 else
1157 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1159 gmx_fio_do_real(fio, ir->rcoulomb_switch);
1160 gmx_fio_do_real(fio, ir->rcoulomb);
1161 gmx_fio_do_int(fio, ir->vdwtype);
1162 if (file_version >= 81)
1164 gmx_fio_do_int(fio, ir->vdw_modifier);
1166 else
1168 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1170 gmx_fio_do_real(fio, ir->rvdw_switch);
1171 gmx_fio_do_real(fio, ir->rvdw);
1172 gmx_fio_do_int(fio, ir->eDispCorr);
1173 gmx_fio_do_real(fio, ir->epsilon_r);
1174 gmx_fio_do_real(fio, ir->epsilon_rf);
1175 gmx_fio_do_real(fio, ir->tabext);
1177 // This permits reading a .tpr file that used implicit solvent,
1178 // and later permitting mdrun to refuse to run it.
1179 if (bRead)
1181 if (file_version < tpxv_RemoveImplicitSolvation)
1183 gmx_fio_do_int(fio, idum);
1184 gmx_fio_do_int(fio, idum);
1185 gmx_fio_do_real(fio, rdum);
1186 gmx_fio_do_real(fio, rdum);
1187 gmx_fio_do_int(fio, idum);
1188 ir->implicit_solvent = (idum > 0);
1190 else
1192 ir->implicit_solvent = false;
1194 if (file_version < tpxv_RemoveImplicitSolvation)
1196 gmx_fio_do_real(fio, rdum);
1197 gmx_fio_do_real(fio, rdum);
1198 gmx_fio_do_real(fio, rdum);
1199 gmx_fio_do_real(fio, rdum);
1200 if (file_version >= 60)
1202 gmx_fio_do_real(fio, rdum);
1203 gmx_fio_do_int(fio, idum);
1205 gmx_fio_do_real(fio, rdum);
1209 if (file_version >= 81)
1211 gmx_fio_do_real(fio, ir->fourier_spacing);
1213 else
1215 ir->fourier_spacing = 0.0;
1217 gmx_fio_do_int(fio, ir->nkx);
1218 gmx_fio_do_int(fio, ir->nky);
1219 gmx_fio_do_int(fio, ir->nkz);
1220 gmx_fio_do_int(fio, ir->pme_order);
1221 gmx_fio_do_real(fio, ir->ewald_rtol);
1223 if (file_version >= 93)
1225 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1227 else
1229 ir->ewald_rtol_lj = ir->ewald_rtol;
1231 gmx_fio_do_int(fio, ir->ewald_geometry);
1232 gmx_fio_do_real(fio, ir->epsilon_surface);
1234 /* ignore bOptFFT */
1235 if (file_version < tpxv_RemoveObsoleteParameters1)
1237 gmx_fio_do_gmx_bool(fio, bdum);
1240 if (file_version >= 93)
1242 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1244 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1245 gmx_fio_do_int(fio, ir->etc);
1246 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1247 * but the values 0 and 1 still mean no and
1248 * berendsen temperature coupling, respectively.
1250 if (file_version >= 79)
1252 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1254 if (file_version >= 71)
1256 gmx_fio_do_int(fio, ir->nsttcouple);
1258 else
1260 ir->nsttcouple = ir->nstcalcenergy;
1262 gmx_fio_do_int(fio, ir->epc);
1263 gmx_fio_do_int(fio, ir->epct);
1264 if (file_version >= 71)
1266 gmx_fio_do_int(fio, ir->nstpcouple);
1268 else
1270 ir->nstpcouple = ir->nstcalcenergy;
1272 gmx_fio_do_real(fio, ir->tau_p);
1273 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1274 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1275 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1276 gmx_fio_do_rvec(fio, ir->compress[XX]);
1277 gmx_fio_do_rvec(fio, ir->compress[YY]);
1278 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1279 gmx_fio_do_int(fio, ir->refcoord_scaling);
1280 gmx_fio_do_rvec(fio, ir->posres_com);
1281 gmx_fio_do_rvec(fio, ir->posres_comB);
1283 if (file_version < 79)
1285 gmx_fio_do_int(fio, ir->andersen_seed);
1287 else
1289 ir->andersen_seed = 0;
1292 gmx_fio_do_real(fio, ir->shake_tol);
1294 gmx_fio_do_int(fio, ir->efep);
1295 do_fepvals(fio, ir->fepvals, bRead, file_version);
1297 if (file_version >= 79)
1299 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1300 if (ir->bSimTemp)
1302 ir->bSimTemp = TRUE;
1305 else
1307 ir->bSimTemp = FALSE;
1309 if (ir->bSimTemp)
1311 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1314 if (file_version >= 79)
1316 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1317 if (ir->bExpanded)
1319 ir->bExpanded = TRUE;
1321 else
1323 ir->bExpanded = FALSE;
1326 if (ir->bExpanded)
1328 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1331 gmx_fio_do_int(fio, ir->eDisre);
1332 gmx_fio_do_int(fio, ir->eDisreWeighting);
1333 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1334 gmx_fio_do_real(fio, ir->dr_fc);
1335 gmx_fio_do_real(fio, ir->dr_tau);
1336 gmx_fio_do_int(fio, ir->nstdisreout);
1337 gmx_fio_do_real(fio, ir->orires_fc);
1338 gmx_fio_do_real(fio, ir->orires_tau);
1339 gmx_fio_do_int(fio, ir->nstorireout);
1341 /* ignore dihre_fc */
1342 if (file_version < 79)
1344 gmx_fio_do_real(fio, rdum);
1347 gmx_fio_do_real(fio, ir->em_stepsize);
1348 gmx_fio_do_real(fio, ir->em_tol);
1349 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1350 gmx_fio_do_int(fio, ir->niter);
1351 gmx_fio_do_real(fio, ir->fc_stepsize);
1352 gmx_fio_do_int(fio, ir->eConstrAlg);
1353 gmx_fio_do_int(fio, ir->nProjOrder);
1354 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1355 gmx_fio_do_int(fio, ir->nLincsIter);
1356 gmx_fio_do_real(fio, ir->bd_fric);
1357 if (file_version >= tpxv_Use64BitRandomSeed)
1359 gmx_fio_do_int64(fio, ir->ld_seed);
1361 else
1363 gmx_fio_do_int(fio, idum);
1364 ir->ld_seed = idum;
1367 for (i = 0; i < DIM; i++)
1369 gmx_fio_do_rvec(fio, ir->deform[i]);
1371 gmx_fio_do_real(fio, ir->cos_accel);
1373 gmx_fio_do_int(fio, ir->userint1);
1374 gmx_fio_do_int(fio, ir->userint2);
1375 gmx_fio_do_int(fio, ir->userint3);
1376 gmx_fio_do_int(fio, ir->userint4);
1377 gmx_fio_do_real(fio, ir->userreal1);
1378 gmx_fio_do_real(fio, ir->userreal2);
1379 gmx_fio_do_real(fio, ir->userreal3);
1380 gmx_fio_do_real(fio, ir->userreal4);
1382 /* AdResS is removed, but we need to be able to read old files,
1383 and let mdrun refuse to run them */
1384 if (file_version >= 77 && file_version < tpxv_RemoveAdress)
1386 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1387 if (ir->bAdress)
1389 int idum, numThermoForceGroups, numEnergyGroups;
1390 real rdum;
1391 rvec rvecdum;
1392 gmx_fio_do_int(fio, idum);
1393 gmx_fio_do_real(fio, rdum);
1394 gmx_fio_do_real(fio, rdum);
1395 gmx_fio_do_real(fio, rdum);
1396 gmx_fio_do_int(fio, idum);
1397 gmx_fio_do_int(fio, idum);
1398 gmx_fio_do_rvec(fio, rvecdum);
1399 gmx_fio_do_int(fio, numThermoForceGroups);
1400 gmx_fio_do_real(fio, rdum);
1401 gmx_fio_do_int(fio, numEnergyGroups);
1402 gmx_fio_do_int(fio, idum);
1404 if (numThermoForceGroups > 0)
1406 std::vector<int> idumn(numThermoForceGroups);
1407 gmx_fio_ndo_int(fio, idumn.data(), idumn.size());
1409 if (numEnergyGroups > 0)
1411 std::vector<int> idumn(numEnergyGroups);
1412 gmx_fio_ndo_int(fio, idumn.data(), idumn.size());
1416 else
1418 ir->bAdress = FALSE;
1421 /* pull stuff */
1423 int ePullOld = 0;
1425 if (file_version >= tpxv_PullCoordTypeGeom)
1427 gmx_fio_do_gmx_bool(fio, ir->bPull);
1429 else
1431 gmx_fio_do_int(fio, ePullOld);
1432 ir->bPull = (ePullOld > 0);
1433 /* We removed the first ePull=ePullNo for the enum */
1434 ePullOld -= 1;
1436 if (ir->bPull)
1438 if (bRead)
1440 snew(ir->pull, 1);
1442 do_pull(fio, ir->pull, bRead, file_version, ePullOld);
1446 if (file_version >= tpxv_AcceleratedWeightHistogram)
1448 gmx_fio_do_gmx_bool(fio, ir->bDoAwh);
1450 if (ir->bDoAwh)
1452 if (bRead)
1454 snew(ir->awhParams, 1);
1456 do_awh(fio, ir->awhParams, bRead, file_version);
1459 else
1461 ir->bDoAwh = FALSE;
1464 /* Enforced rotation */
1465 if (file_version >= 74)
1467 gmx_fio_do_int(fio, ir->bRot);
1468 if (ir->bRot == TRUE)
1470 if (bRead)
1472 snew(ir->rot, 1);
1474 do_rot(fio, ir->rot, bRead);
1477 else
1479 ir->bRot = FALSE;
1482 /* Interactive molecular dynamics */
1483 if (file_version >= tpxv_InteractiveMolecularDynamics)
1485 gmx_fio_do_int(fio, ir->bIMD);
1486 if (TRUE == ir->bIMD)
1488 if (bRead)
1490 snew(ir->imd, 1);
1492 do_imd(fio, ir->imd, bRead);
1495 else
1497 /* We don't support IMD sessions for old .tpr files */
1498 ir->bIMD = FALSE;
1501 /* grpopts stuff */
1502 gmx_fio_do_int(fio, ir->opts.ngtc);
1503 if (file_version >= 69)
1505 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1507 else
1509 ir->opts.nhchainlength = 1;
1511 gmx_fio_do_int(fio, ir->opts.ngacc);
1512 gmx_fio_do_int(fio, ir->opts.ngfrz);
1513 gmx_fio_do_int(fio, ir->opts.ngener);
1515 if (bRead)
1517 snew(ir->opts.nrdf, ir->opts.ngtc);
1518 snew(ir->opts.ref_t, ir->opts.ngtc);
1519 snew(ir->opts.annealing, ir->opts.ngtc);
1520 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1521 snew(ir->opts.anneal_time, ir->opts.ngtc);
1522 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1523 snew(ir->opts.tau_t, ir->opts.ngtc);
1524 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1525 snew(ir->opts.acc, ir->opts.ngacc);
1526 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1528 if (ir->opts.ngtc > 0)
1530 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1531 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1532 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1534 if (ir->opts.ngfrz > 0)
1536 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1538 if (ir->opts.ngacc > 0)
1540 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1542 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1543 ir->opts.ngener*ir->opts.ngener);
1545 /* First read the lists with annealing and npoints for each group */
1546 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1547 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1548 for (j = 0; j < (ir->opts.ngtc); j++)
1550 k = ir->opts.anneal_npoints[j];
1551 if (bRead)
1553 snew(ir->opts.anneal_time[j], k);
1554 snew(ir->opts.anneal_temp[j], k);
1556 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1557 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1559 /* Walls */
1561 gmx_fio_do_int(fio, ir->nwall);
1562 gmx_fio_do_int(fio, ir->wall_type);
1563 gmx_fio_do_real(fio, ir->wall_r_linpot);
1564 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1565 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1566 gmx_fio_do_real(fio, ir->wall_density[0]);
1567 gmx_fio_do_real(fio, ir->wall_density[1]);
1568 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1571 /* Cosine stuff for electric fields */
1572 if (file_version < tpxv_GenericParamsForElectricField)
1574 do_legacy_efield(fio, &paramsObj);
1577 /* Swap ions */
1578 if (file_version >= tpxv_ComputationalElectrophysiology)
1580 gmx_fio_do_int(fio, ir->eSwapCoords);
1581 if (ir->eSwapCoords != eswapNO)
1583 if (bRead)
1585 snew(ir->swap, 1);
1587 do_swapcoords_tpx(fio, ir->swap, bRead, file_version);
1591 /* QMMM stuff */
1593 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1594 gmx_fio_do_int(fio, ir->QMMMscheme);
1595 gmx_fio_do_real(fio, ir->scalefactor);
1596 gmx_fio_do_int(fio, ir->opts.ngQM);
1597 if (bRead)
1599 snew(ir->opts.QMmethod, ir->opts.ngQM);
1600 snew(ir->opts.QMbasis, ir->opts.ngQM);
1601 snew(ir->opts.QMcharge, ir->opts.ngQM);
1602 snew(ir->opts.QMmult, ir->opts.ngQM);
1603 snew(ir->opts.bSH, ir->opts.ngQM);
1604 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1605 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1606 snew(ir->opts.SAon, ir->opts.ngQM);
1607 snew(ir->opts.SAoff, ir->opts.ngQM);
1608 snew(ir->opts.SAsteps, ir->opts.ngQM);
1610 if (ir->opts.ngQM > 0)
1612 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1613 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1614 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1615 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1616 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1617 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1618 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1619 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1620 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1621 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1622 /* We leave in dummy i/o for removed parameters to avoid
1623 * changing the tpr format for every QMMM change.
1625 std::vector<int> dummy(ir->opts.ngQM, 0);
1626 gmx_fio_ndo_gmx_bool(fio, dummy.data(), ir->opts.ngQM);
1627 gmx_fio_ndo_gmx_bool(fio, dummy.data(), ir->opts.ngQM);
1629 /* end of QMMM stuff */
1632 if (file_version >= tpxv_GenericParamsForElectricField)
1634 gmx::FileIOXdrSerializer serializer(fio);
1635 if (bRead)
1637 paramsObj.mergeObject(
1638 gmx::deserializeKeyValueTree(&serializer));
1640 else
1642 GMX_RELEASE_ASSERT(ir->params != nullptr,
1643 "Parameters should be present when writing inputrec");
1644 gmx::serializeKeyValueTree(*ir->params, &serializer);
1647 if (bRead)
1649 ir->params = new gmx::KeyValueTreeObject(paramsBuilder.build());
1654 static void do_harm(t_fileio *fio, t_iparams *iparams)
1656 gmx_fio_do_real(fio, iparams->harmonic.rA);
1657 gmx_fio_do_real(fio, iparams->harmonic.krA);
1658 gmx_fio_do_real(fio, iparams->harmonic.rB);
1659 gmx_fio_do_real(fio, iparams->harmonic.krB);
1662 static void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1663 gmx_bool bRead, int file_version)
1665 int idum;
1666 real rdum;
1668 switch (ftype)
1670 case F_ANGLES:
1671 case F_G96ANGLES:
1672 case F_BONDS:
1673 case F_G96BONDS:
1674 case F_HARMONIC:
1675 case F_IDIHS:
1676 do_harm(fio, iparams);
1677 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1679 /* Correct incorrect storage of parameters */
1680 iparams->pdihs.phiB = iparams->pdihs.phiA;
1681 iparams->pdihs.cpB = iparams->pdihs.cpA;
1683 break;
1684 case F_RESTRANGLES:
1685 gmx_fio_do_real(fio, iparams->harmonic.rA);
1686 gmx_fio_do_real(fio, iparams->harmonic.krA);
1687 break;
1688 case F_LINEAR_ANGLES:
1689 gmx_fio_do_real(fio, iparams->linangle.klinA);
1690 gmx_fio_do_real(fio, iparams->linangle.aA);
1691 gmx_fio_do_real(fio, iparams->linangle.klinB);
1692 gmx_fio_do_real(fio, iparams->linangle.aB);
1693 break;
1694 case F_FENEBONDS:
1695 gmx_fio_do_real(fio, iparams->fene.bm);
1696 gmx_fio_do_real(fio, iparams->fene.kb);
1697 break;
1699 case F_RESTRBONDS:
1700 gmx_fio_do_real(fio, iparams->restraint.lowA);
1701 gmx_fio_do_real(fio, iparams->restraint.up1A);
1702 gmx_fio_do_real(fio, iparams->restraint.up2A);
1703 gmx_fio_do_real(fio, iparams->restraint.kA);
1704 gmx_fio_do_real(fio, iparams->restraint.lowB);
1705 gmx_fio_do_real(fio, iparams->restraint.up1B);
1706 gmx_fio_do_real(fio, iparams->restraint.up2B);
1707 gmx_fio_do_real(fio, iparams->restraint.kB);
1708 break;
1709 case F_TABBONDS:
1710 case F_TABBONDSNC:
1711 case F_TABANGLES:
1712 case F_TABDIHS:
1713 gmx_fio_do_real(fio, iparams->tab.kA);
1714 gmx_fio_do_int(fio, iparams->tab.table);
1715 gmx_fio_do_real(fio, iparams->tab.kB);
1716 break;
1717 case F_CROSS_BOND_BONDS:
1718 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1719 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1720 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1721 break;
1722 case F_CROSS_BOND_ANGLES:
1723 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1724 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1725 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1726 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1727 break;
1728 case F_UREY_BRADLEY:
1729 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1730 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1731 gmx_fio_do_real(fio, iparams->u_b.r13A);
1732 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1733 if (file_version >= 79)
1735 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1736 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1737 gmx_fio_do_real(fio, iparams->u_b.r13B);
1738 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1740 else
1742 iparams->u_b.thetaB = iparams->u_b.thetaA;
1743 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1744 iparams->u_b.r13B = iparams->u_b.r13A;
1745 iparams->u_b.kUBB = iparams->u_b.kUBA;
1747 break;
1748 case F_QUARTIC_ANGLES:
1749 gmx_fio_do_real(fio, iparams->qangle.theta);
1750 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1751 break;
1752 case F_BHAM:
1753 gmx_fio_do_real(fio, iparams->bham.a);
1754 gmx_fio_do_real(fio, iparams->bham.b);
1755 gmx_fio_do_real(fio, iparams->bham.c);
1756 break;
1757 case F_MORSE:
1758 gmx_fio_do_real(fio, iparams->morse.b0A);
1759 gmx_fio_do_real(fio, iparams->morse.cbA);
1760 gmx_fio_do_real(fio, iparams->morse.betaA);
1761 if (file_version >= 79)
1763 gmx_fio_do_real(fio, iparams->morse.b0B);
1764 gmx_fio_do_real(fio, iparams->morse.cbB);
1765 gmx_fio_do_real(fio, iparams->morse.betaB);
1767 else
1769 iparams->morse.b0B = iparams->morse.b0A;
1770 iparams->morse.cbB = iparams->morse.cbA;
1771 iparams->morse.betaB = iparams->morse.betaA;
1773 break;
1774 case F_CUBICBONDS:
1775 gmx_fio_do_real(fio, iparams->cubic.b0);
1776 gmx_fio_do_real(fio, iparams->cubic.kb);
1777 gmx_fio_do_real(fio, iparams->cubic.kcub);
1778 break;
1779 case F_CONNBONDS:
1780 break;
1781 case F_POLARIZATION:
1782 gmx_fio_do_real(fio, iparams->polarize.alpha);
1783 break;
1784 case F_ANHARM_POL:
1785 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1786 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1787 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1788 break;
1789 case F_WATER_POL:
1790 gmx_fio_do_real(fio, iparams->wpol.al_x);
1791 gmx_fio_do_real(fio, iparams->wpol.al_y);
1792 gmx_fio_do_real(fio, iparams->wpol.al_z);
1793 gmx_fio_do_real(fio, iparams->wpol.rOH);
1794 gmx_fio_do_real(fio, iparams->wpol.rHH);
1795 gmx_fio_do_real(fio, iparams->wpol.rOD);
1796 break;
1797 case F_THOLE_POL:
1798 gmx_fio_do_real(fio, iparams->thole.a);
1799 gmx_fio_do_real(fio, iparams->thole.alpha1);
1800 gmx_fio_do_real(fio, iparams->thole.alpha2);
1801 gmx_fio_do_real(fio, iparams->thole.rfac);
1802 break;
1803 case F_LJ:
1804 gmx_fio_do_real(fio, iparams->lj.c6);
1805 gmx_fio_do_real(fio, iparams->lj.c12);
1806 break;
1807 case F_LJ14:
1808 gmx_fio_do_real(fio, iparams->lj14.c6A);
1809 gmx_fio_do_real(fio, iparams->lj14.c12A);
1810 gmx_fio_do_real(fio, iparams->lj14.c6B);
1811 gmx_fio_do_real(fio, iparams->lj14.c12B);
1812 break;
1813 case F_LJC14_Q:
1814 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1815 gmx_fio_do_real(fio, iparams->ljc14.qi);
1816 gmx_fio_do_real(fio, iparams->ljc14.qj);
1817 gmx_fio_do_real(fio, iparams->ljc14.c6);
1818 gmx_fio_do_real(fio, iparams->ljc14.c12);
1819 break;
1820 case F_LJC_PAIRS_NB:
1821 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1822 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1823 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1824 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1825 break;
1826 case F_PDIHS:
1827 case F_PIDIHS:
1828 case F_ANGRES:
1829 case F_ANGRESZ:
1830 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1831 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1832 gmx_fio_do_real(fio, iparams->pdihs.phiB);
1833 gmx_fio_do_real(fio, iparams->pdihs.cpB);
1834 gmx_fio_do_int(fio, iparams->pdihs.mult);
1835 break;
1836 case F_RESTRDIHS:
1837 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1838 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1839 break;
1840 case F_DISRES:
1841 gmx_fio_do_int(fio, iparams->disres.label);
1842 gmx_fio_do_int(fio, iparams->disres.type);
1843 gmx_fio_do_real(fio, iparams->disres.low);
1844 gmx_fio_do_real(fio, iparams->disres.up1);
1845 gmx_fio_do_real(fio, iparams->disres.up2);
1846 gmx_fio_do_real(fio, iparams->disres.kfac);
1847 break;
1848 case F_ORIRES:
1849 gmx_fio_do_int(fio, iparams->orires.ex);
1850 gmx_fio_do_int(fio, iparams->orires.label);
1851 gmx_fio_do_int(fio, iparams->orires.power);
1852 gmx_fio_do_real(fio, iparams->orires.c);
1853 gmx_fio_do_real(fio, iparams->orires.obs);
1854 gmx_fio_do_real(fio, iparams->orires.kfac);
1855 break;
1856 case F_DIHRES:
1857 if (file_version < 82)
1859 gmx_fio_do_int(fio, idum);
1860 gmx_fio_do_int(fio, idum);
1862 gmx_fio_do_real(fio, iparams->dihres.phiA);
1863 gmx_fio_do_real(fio, iparams->dihres.dphiA);
1864 gmx_fio_do_real(fio, iparams->dihres.kfacA);
1865 if (file_version >= 82)
1867 gmx_fio_do_real(fio, iparams->dihres.phiB);
1868 gmx_fio_do_real(fio, iparams->dihres.dphiB);
1869 gmx_fio_do_real(fio, iparams->dihres.kfacB);
1871 else
1873 iparams->dihres.phiB = iparams->dihres.phiA;
1874 iparams->dihres.dphiB = iparams->dihres.dphiA;
1875 iparams->dihres.kfacB = iparams->dihres.kfacA;
1877 break;
1878 case F_POSRES:
1879 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
1880 gmx_fio_do_rvec(fio, iparams->posres.fcA);
1881 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
1882 gmx_fio_do_rvec(fio, iparams->posres.fcB);
1883 break;
1884 case F_FBPOSRES:
1885 gmx_fio_do_int(fio, iparams->fbposres.geom);
1886 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
1887 gmx_fio_do_real(fio, iparams->fbposres.r);
1888 gmx_fio_do_real(fio, iparams->fbposres.k);
1889 break;
1890 case F_CBTDIHS:
1891 gmx_fio_ndo_real(fio, iparams->cbtdihs.cbtcA, NR_CBTDIHS);
1892 break;
1893 case F_RBDIHS:
1894 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1895 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1896 break;
1897 case F_FOURDIHS:
1898 /* Fourier dihedrals are internally represented
1899 * as Ryckaert-Bellemans since those are faster to compute.
1901 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1902 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1903 break;
1904 case F_CONSTR:
1905 case F_CONSTRNC:
1906 gmx_fio_do_real(fio, iparams->constr.dA);
1907 gmx_fio_do_real(fio, iparams->constr.dB);
1908 break;
1909 case F_SETTLE:
1910 gmx_fio_do_real(fio, iparams->settle.doh);
1911 gmx_fio_do_real(fio, iparams->settle.dhh);
1912 break;
1913 case F_VSITE2:
1914 gmx_fio_do_real(fio, iparams->vsite.a);
1915 break;
1916 case F_VSITE3:
1917 case F_VSITE3FD:
1918 case F_VSITE3FAD:
1919 gmx_fio_do_real(fio, iparams->vsite.a);
1920 gmx_fio_do_real(fio, iparams->vsite.b);
1921 break;
1922 case F_VSITE3OUT:
1923 case F_VSITE4FD:
1924 case F_VSITE4FDN:
1925 gmx_fio_do_real(fio, iparams->vsite.a);
1926 gmx_fio_do_real(fio, iparams->vsite.b);
1927 gmx_fio_do_real(fio, iparams->vsite.c);
1928 break;
1929 case F_VSITEN:
1930 gmx_fio_do_int(fio, iparams->vsiten.n);
1931 gmx_fio_do_real(fio, iparams->vsiten.a);
1932 break;
1933 case F_GB12_NOLONGERUSED:
1934 case F_GB13_NOLONGERUSED:
1935 case F_GB14_NOLONGERUSED:
1936 // Implicit solvent parameters can still be read, but never used
1937 if (bRead)
1939 if (file_version < 68)
1941 gmx_fio_do_real(fio, rdum);
1942 gmx_fio_do_real(fio, rdum);
1943 gmx_fio_do_real(fio, rdum);
1944 gmx_fio_do_real(fio, rdum);
1946 if (file_version < tpxv_RemoveImplicitSolvation)
1948 gmx_fio_do_real(fio, rdum);
1949 gmx_fio_do_real(fio, rdum);
1950 gmx_fio_do_real(fio, rdum);
1951 gmx_fio_do_real(fio, rdum);
1952 gmx_fio_do_real(fio, rdum);
1955 break;
1956 case F_CMAP:
1957 gmx_fio_do_int(fio, iparams->cmap.cmapA);
1958 gmx_fio_do_int(fio, iparams->cmap.cmapB);
1959 break;
1960 default:
1961 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1962 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1966 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead)
1968 gmx_fio_do_int(fio, ilist->nr);
1969 if (bRead)
1971 snew(ilist->iatoms, ilist->nr);
1973 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
1976 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1977 gmx_bool bRead, int file_version)
1979 gmx_fio_do_int(fio, ffparams->atnr);
1980 gmx_fio_do_int(fio, ffparams->ntypes);
1981 if (bRead)
1983 snew(ffparams->functype, ffparams->ntypes);
1984 snew(ffparams->iparams, ffparams->ntypes);
1986 /* Read/write all the function types */
1987 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
1989 if (file_version >= 66)
1991 gmx_fio_do_double(fio, ffparams->reppow);
1993 else
1995 ffparams->reppow = 12.0;
1998 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2000 /* Check whether all these function types are supported by the code.
2001 * In practice the code is backwards compatible, which means that the
2002 * numbering may have to be altered from old numbering to new numbering
2004 for (int i = 0; i < ffparams->ntypes; i++)
2006 if (bRead)
2008 /* Loop over file versions */
2009 for (int k = 0; k < NFTUPD; k++)
2011 /* Compare the read file_version to the update table */
2012 if ((file_version < ftupd[k].fvnr) &&
2013 (ffparams->functype[i] >= ftupd[k].ftype))
2015 ffparams->functype[i] += 1;
2020 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2021 file_version);
2025 static void add_settle_atoms(t_ilist *ilist)
2027 int i;
2029 /* Settle used to only store the first atom: add the other two */
2030 srenew(ilist->iatoms, 2*ilist->nr);
2031 for (i = ilist->nr/2-1; i >= 0; i--)
2033 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2034 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2035 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2036 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2038 ilist->nr = 2*ilist->nr;
2041 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2042 int file_version)
2044 int j;
2045 gmx_bool bClear;
2046 unsigned int k;
2048 for (j = 0; (j < F_NRE); j++)
2050 bClear = FALSE;
2051 if (bRead)
2053 for (k = 0; k < NFTUPD; k++)
2055 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2057 bClear = TRUE;
2061 if (bClear)
2063 ilist[j].nr = 0;
2064 ilist[j].iatoms = nullptr;
2066 else
2068 do_ilist(fio, &ilist[j], bRead);
2069 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2071 add_settle_atoms(&ilist[j]);
2077 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead)
2079 gmx_fio_do_int(fio, block->nr);
2080 if (bRead)
2082 if ((block->nalloc_index > 0) && (nullptr != block->index))
2084 sfree(block->index);
2086 block->nalloc_index = block->nr+1;
2087 snew(block->index, block->nalloc_index);
2089 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2092 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead)
2094 gmx_fio_do_int(fio, block->nr);
2095 gmx_fio_do_int(fio, block->nra);
2096 if (bRead)
2098 block->nalloc_index = block->nr+1;
2099 snew(block->index, block->nalloc_index);
2100 block->nalloc_a = block->nra;
2101 snew(block->a, block->nalloc_a);
2103 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2104 gmx_fio_ndo_int(fio, block->a, block->nra);
2107 /* This is a primitive routine to make it possible to translate atomic numbers
2108 * to element names when reading TPR files, without making the Gromacs library
2109 * directory a dependency on mdrun (which is the case if we need elements.dat).
2111 static const char *
2112 atomicnumber_to_element(int atomicnumber)
2114 const char * p;
2116 /* This does not have to be complete, so we only include elements likely
2117 * to occur in PDB files.
2119 switch (atomicnumber)
2121 case 1: p = "H"; break;
2122 case 5: p = "B"; break;
2123 case 6: p = "C"; break;
2124 case 7: p = "N"; break;
2125 case 8: p = "O"; break;
2126 case 9: p = "F"; break;
2127 case 11: p = "Na"; break;
2128 case 12: p = "Mg"; break;
2129 case 15: p = "P"; break;
2130 case 16: p = "S"; break;
2131 case 17: p = "Cl"; break;
2132 case 18: p = "Ar"; break;
2133 case 19: p = "K"; break;
2134 case 20: p = "Ca"; break;
2135 case 25: p = "Mn"; break;
2136 case 26: p = "Fe"; break;
2137 case 28: p = "Ni"; break;
2138 case 29: p = "Cu"; break;
2139 case 30: p = "Zn"; break;
2140 case 35: p = "Br"; break;
2141 case 47: p = "Ag"; break;
2142 default: p = ""; break;
2144 return p;
2148 static void do_atom(t_fileio *fio, t_atom *atom, gmx_bool bRead)
2150 gmx_fio_do_real(fio, atom->m);
2151 gmx_fio_do_real(fio, atom->q);
2152 gmx_fio_do_real(fio, atom->mB);
2153 gmx_fio_do_real(fio, atom->qB);
2154 gmx_fio_do_ushort(fio, atom->type);
2155 gmx_fio_do_ushort(fio, atom->typeB);
2156 gmx_fio_do_int(fio, atom->ptype);
2157 gmx_fio_do_int(fio, atom->resind);
2158 gmx_fio_do_int(fio, atom->atomnumber);
2159 if (bRead)
2161 /* Set element string from atomic number if present.
2162 * This routine returns an empty string if the name is not found.
2164 std::strncpy(atom->elem, atomicnumber_to_element(atom->atomnumber), 4);
2165 /* avoid warnings about potentially unterminated string */
2166 atom->elem[3] = '\0';
2170 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead)
2172 for (int j = 0; j < ngrp; j++)
2174 gmx_fio_do_int(fio, grps[j].nr);
2175 if (bRead)
2177 snew(grps[j].nm_ind, grps[j].nr);
2179 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2183 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2185 int ls;
2187 if (bRead)
2189 gmx_fio_do_int(fio, ls);
2190 *nm = get_symtab_handle(symtab, ls);
2192 else
2194 ls = lookup_symtab(symtab, *nm);
2195 gmx_fio_do_int(fio, ls);
2199 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2200 t_symtab *symtab)
2202 int j;
2204 for (j = 0; (j < nstr); j++)
2206 do_symstr(fio, &(nm[j]), bRead, symtab);
2210 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2211 t_symtab *symtab, int file_version)
2213 int j;
2215 for (j = 0; (j < n); j++)
2217 do_symstr(fio, &(ri[j].name), bRead, symtab);
2218 if (file_version >= 63)
2220 gmx_fio_do_int(fio, ri[j].nr);
2221 gmx_fio_do_uchar(fio, ri[j].ic);
2223 else
2225 ri[j].nr = j + 1;
2226 ri[j].ic = ' ';
2231 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2232 int file_version)
2234 int i;
2236 gmx_fio_do_int(fio, atoms->nr);
2237 gmx_fio_do_int(fio, atoms->nres);
2238 if (bRead)
2240 /* Since we have always written all t_atom properties in the tpr file
2241 * (at least for all backward compatible versions), we don't store
2242 * but simple set the booleans here.
2244 atoms->haveMass = TRUE;
2245 atoms->haveCharge = TRUE;
2246 atoms->haveType = TRUE;
2247 atoms->haveBState = TRUE;
2248 atoms->havePdbInfo = FALSE;
2250 snew(atoms->atom, atoms->nr);
2251 snew(atoms->atomname, atoms->nr);
2252 snew(atoms->atomtype, atoms->nr);
2253 snew(atoms->atomtypeB, atoms->nr);
2254 snew(atoms->resinfo, atoms->nres);
2255 atoms->pdbinfo = nullptr;
2257 else
2259 GMX_RELEASE_ASSERT(atoms->haveMass && atoms->haveCharge && atoms->haveType && atoms->haveBState, "Mass, charge, atomtype and B-state parameters should be present in t_atoms when writing a tpr file");
2261 for (i = 0; (i < atoms->nr); i++)
2263 do_atom(fio, &atoms->atom[i], bRead);
2265 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2266 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2267 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2269 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2272 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2273 gmx_bool bRead, t_symtab *symtab)
2275 int g;
2277 do_grps(fio, egcNR, groups->grps, bRead);
2278 gmx_fio_do_int(fio, groups->ngrpname);
2279 if (bRead)
2281 snew(groups->grpname, groups->ngrpname);
2283 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2284 for (g = 0; g < egcNR; g++)
2286 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2287 if (groups->ngrpnr[g] == 0)
2289 if (bRead)
2291 groups->grpnr[g] = nullptr;
2294 else
2296 if (bRead)
2298 snew(groups->grpnr[g], groups->ngrpnr[g]);
2300 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2305 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2306 int file_version)
2308 int j;
2310 gmx_fio_do_int(fio, atomtypes->nr);
2311 j = atomtypes->nr;
2312 if (bRead)
2314 snew(atomtypes->atomnumber, j);
2316 if (bRead && file_version < tpxv_RemoveImplicitSolvation)
2318 std::vector<real> dummy(atomtypes->nr, 0);
2319 gmx_fio_ndo_real(fio, dummy.data(), dummy.size());
2320 gmx_fio_ndo_real(fio, dummy.data(), dummy.size());
2321 gmx_fio_ndo_real(fio, dummy.data(), dummy.size());
2323 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2325 if (bRead && file_version >= 60 && file_version < tpxv_RemoveImplicitSolvation)
2327 std::vector<real> dummy(atomtypes->nr, 0);
2328 gmx_fio_ndo_real(fio, dummy.data(), dummy.size());
2329 gmx_fio_ndo_real(fio, dummy.data(), dummy.size());
2333 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2335 int i, nr;
2336 t_symbuf *symbuf;
2337 char buf[STRLEN];
2339 gmx_fio_do_int(fio, symtab->nr);
2340 nr = symtab->nr;
2341 if (bRead)
2343 snew(symtab->symbuf, 1);
2344 symbuf = symtab->symbuf;
2345 symbuf->bufsize = nr;
2346 snew(symbuf->buf, nr);
2347 for (i = 0; (i < nr); i++)
2349 gmx_fio_do_string(fio, buf);
2350 symbuf->buf[i] = gmx_strdup(buf);
2353 else
2355 symbuf = symtab->symbuf;
2356 while (symbuf != nullptr)
2358 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2360 gmx_fio_do_string(fio, symbuf->buf[i]);
2362 nr -= i;
2363 symbuf = symbuf->next;
2365 if (nr != 0)
2367 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2372 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2374 int i, j, ngrid, gs, nelem;
2376 gmx_fio_do_int(fio, cmap_grid->ngrid);
2377 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2379 ngrid = cmap_grid->ngrid;
2380 gs = cmap_grid->grid_spacing;
2381 nelem = gs * gs;
2383 if (bRead)
2385 snew(cmap_grid->cmapdata, ngrid);
2387 for (i = 0; i < cmap_grid->ngrid; i++)
2389 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2393 for (i = 0; i < cmap_grid->ngrid; i++)
2395 for (j = 0; j < nelem; j++)
2397 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2398 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2399 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2400 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2406 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2407 t_symtab *symtab, int file_version)
2409 do_symstr(fio, &(molt->name), bRead, symtab);
2411 do_atoms(fio, &molt->atoms, bRead, symtab, file_version);
2413 do_ilists(fio, molt->ilist, bRead, file_version);
2415 do_block(fio, &molt->cgs, bRead);
2417 /* This used to be in the atoms struct */
2418 do_blocka(fio, &molt->excls, bRead);
2421 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb,
2422 int numAtomsPerMolecule,
2423 gmx_bool bRead)
2425 gmx_fio_do_int(fio, molb->type);
2426 gmx_fio_do_int(fio, molb->nmol);
2427 /* To maintain forward topology reading compatibility, we store #atoms.
2428 * TODO: Change this to conditional reading of a dummy int when we
2429 * increase tpx_generation.
2431 gmx_fio_do_int(fio, numAtomsPerMolecule);
2432 /* Position restraint coordinates */
2433 int numPosres_xA = molb->posres_xA.size();
2434 gmx_fio_do_int(fio, numPosres_xA);
2435 if (numPosres_xA > 0)
2437 if (bRead)
2439 molb->posres_xA.resize(numPosres_xA);
2441 gmx_fio_ndo_rvec(fio, as_rvec_array(molb->posres_xA.data()), numPosres_xA);
2443 int numPosres_xB = molb->posres_xB.size();
2444 gmx_fio_do_int(fio, numPosres_xB);
2445 if (numPosres_xB > 0)
2447 if (bRead)
2449 molb->posres_xB.resize(numPosres_xB);
2451 gmx_fio_ndo_rvec(fio, as_rvec_array(molb->posres_xB.data()), numPosres_xB);
2456 static void set_disres_npair(gmx_mtop_t *mtop)
2458 t_iparams *ip;
2459 gmx_mtop_ilistloop_t iloop;
2460 const t_ilist *ilist, *il;
2461 int nmol, i, npair;
2462 t_iatom *a;
2464 ip = mtop->ffparams.iparams;
2466 iloop = gmx_mtop_ilistloop_init(mtop);
2467 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
2469 il = &ilist[F_DISRES];
2471 if (il->nr > 0)
2473 a = il->iatoms;
2474 npair = 0;
2475 for (i = 0; i < il->nr; i += 3)
2477 npair++;
2478 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2480 ip[a[i]].disres.npair = npair;
2481 npair = 0;
2488 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2489 int file_version)
2491 if (bRead)
2493 init_mtop(mtop);
2495 do_symtab(fio, &(mtop->symtab), bRead);
2497 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2499 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
2501 int nmoltype = mtop->moltype.size();
2502 gmx_fio_do_int(fio, nmoltype);
2503 if (bRead)
2505 mtop->moltype.resize(nmoltype);
2507 for (gmx_moltype_t &moltype : mtop->moltype)
2509 do_moltype(fio, &moltype, bRead, &mtop->symtab, file_version);
2512 int nmolblock = mtop->molblock.size();
2513 gmx_fio_do_int(fio, nmolblock);
2514 if (bRead)
2516 mtop->molblock.resize(nmolblock);
2518 for (gmx_molblock_t &molblock : mtop->molblock)
2520 int numAtomsPerMolecule = (bRead ? 0 : mtop->moltype[molblock.type].atoms.nr);
2521 do_molblock(fio, &molblock, numAtomsPerMolecule, bRead);
2523 gmx_fio_do_int(fio, mtop->natoms);
2525 if (file_version >= tpxv_IntermolecularBondeds)
2527 gmx_fio_do_gmx_bool(fio, mtop->bIntermolecularInteractions);
2528 if (mtop->bIntermolecularInteractions)
2530 if (bRead)
2532 snew(mtop->intermolecular_ilist, F_NRE);
2534 do_ilists(fio, mtop->intermolecular_ilist, bRead, file_version);
2537 else
2539 mtop->bIntermolecularInteractions = FALSE;
2542 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
2544 if (file_version >= 65)
2546 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
2548 else
2550 mtop->ffparams.cmap_grid.ngrid = 0;
2551 mtop->ffparams.cmap_grid.grid_spacing = 0;
2552 mtop->ffparams.cmap_grid.cmapdata = nullptr;
2555 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab));
2557 mtop->haveMoleculeIndices = true;
2559 if (bRead)
2561 close_symtab(&(mtop->symtab));
2565 /* If TopOnlyOK is TRUE then we can read even future versions
2566 * of tpx files, provided the fileGeneration hasn't changed.
2567 * If it is FALSE, we need the inputrecord too, and bail out
2568 * if the file is newer than the program.
2570 * The version and generation of the topology (see top of this file)
2571 * are returned in the two last arguments, if those arguments are non-NULL.
2573 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2575 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
2576 gmx_bool TopOnlyOK, int *fileVersionPointer, int *fileGenerationPointer)
2578 char buf[STRLEN];
2579 char file_tag[STRLEN];
2580 gmx_bool bDouble;
2581 int precision;
2582 int idum = 0;
2583 real rdum = 0;
2584 int fileVersion; /* Version number of the code that wrote the file */
2585 int fileGeneration; /* Generation version number of the code that wrote the file */
2587 /* XDR binary topology file */
2588 precision = sizeof(real);
2589 if (bRead)
2591 gmx_fio_do_string(fio, buf);
2592 if (std::strncmp(buf, "VERSION", 7))
2594 gmx_fatal(FARGS, "Can not read file %s,\n"
2595 " this file is from a GROMACS version which is older than 2.0\n"
2596 " Make a new one with grompp or use a gro or pdb file, if possible",
2597 gmx_fio_getname(fio));
2599 gmx_fio_do_int(fio, precision);
2600 bDouble = (precision == sizeof(double));
2601 if ((precision != sizeof(float)) && !bDouble)
2603 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
2604 "instead of %d or %d",
2605 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
2607 gmx_fio_setprecision(fio, bDouble);
2608 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
2609 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
2611 else
2613 snprintf(buf, STRLEN, "VERSION %s", gmx_version());
2614 gmx_fio_write_string(fio, buf);
2615 bDouble = (precision == sizeof(double));
2616 gmx_fio_setprecision(fio, bDouble);
2617 gmx_fio_do_int(fio, precision);
2618 fileVersion = tpx_version;
2619 sprintf(file_tag, "%s", tpx_tag);
2620 fileGeneration = tpx_generation;
2623 /* Check versions! */
2624 gmx_fio_do_int(fio, fileVersion);
2626 /* This is for backward compatibility with development versions 77-79
2627 * where the tag was, mistakenly, placed before the generation,
2628 * which would cause a segv instead of a proper error message
2629 * when reading the topology only from tpx with <77 code.
2631 if (fileVersion >= 77 && fileVersion <= 79)
2633 gmx_fio_do_string(fio, file_tag);
2636 gmx_fio_do_int(fio, fileGeneration);
2638 if (fileVersion >= 81)
2640 gmx_fio_do_string(fio, file_tag);
2642 if (bRead)
2644 if (fileVersion < 77)
2646 /* Versions before 77 don't have the tag, set it to release */
2647 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
2650 if (std::strcmp(file_tag, tpx_tag) != 0)
2652 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
2653 file_tag, tpx_tag);
2655 /* We only support reading tpx files with the same tag as the code
2656 * or tpx files with the release tag and with lower version number.
2658 if (std::strcmp(file_tag, TPX_TAG_RELEASE) != 0 && fileVersion < tpx_version)
2660 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2661 gmx_fio_getname(fio), fileVersion, file_tag,
2662 tpx_version, tpx_tag);
2667 if (fileVersionPointer)
2669 *fileVersionPointer = fileVersion;
2671 if (fileGenerationPointer)
2673 *fileGenerationPointer = fileGeneration;
2676 if ((fileVersion <= tpx_incompatible_version) ||
2677 ((fileVersion > tpx_version) && !TopOnlyOK) ||
2678 (fileGeneration > tpx_generation) ||
2679 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
2681 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
2682 gmx_fio_getname(fio), fileVersion, tpx_version);
2685 gmx_fio_do_int(fio, tpx->natoms);
2686 gmx_fio_do_int(fio, tpx->ngtc);
2688 if (fileVersion < 62)
2690 gmx_fio_do_int(fio, idum);
2691 gmx_fio_do_real(fio, rdum);
2693 if (fileVersion >= 79)
2695 gmx_fio_do_int(fio, tpx->fep_state);
2697 gmx_fio_do_real(fio, tpx->lambda);
2698 gmx_fio_do_int(fio, tpx->bIr);
2699 gmx_fio_do_int(fio, tpx->bTop);
2700 gmx_fio_do_int(fio, tpx->bX);
2701 gmx_fio_do_int(fio, tpx->bV);
2702 gmx_fio_do_int(fio, tpx->bF);
2703 gmx_fio_do_int(fio, tpx->bBox);
2705 if ((fileGeneration > tpx_generation))
2707 /* This can only happen if TopOnlyOK=TRUE */
2708 tpx->bIr = FALSE;
2712 static int do_tpx(t_fileio *fio, gmx_bool bRead,
2713 t_inputrec *ir, t_state *state, rvec *x, rvec *v,
2714 gmx_mtop_t *mtop)
2716 t_tpxheader tpx;
2717 gmx_bool TopOnlyOK;
2718 int ePBC;
2719 gmx_bool bPeriodicMols;
2721 if (!bRead)
2723 GMX_RELEASE_ASSERT(state != nullptr, "Cannot write a null state");
2724 GMX_RELEASE_ASSERT(x == nullptr && v == nullptr, "Passing separate x and v pointers to do_tpx() is not supported when writing");
2726 tpx.natoms = state->natoms;
2727 tpx.ngtc = state->ngtc;
2728 tpx.fep_state = state->fep_state;
2729 tpx.lambda = state->lambda[efptFEP];
2730 tpx.bIr = (ir != nullptr);
2731 tpx.bTop = (mtop != nullptr);
2732 tpx.bX = (state->flags & (1 << estX));
2733 tpx.bV = (state->flags & (1 << estV));
2734 tpx.bF = FALSE;
2735 tpx.bBox = TRUE;
2737 else
2739 GMX_RELEASE_ASSERT(!(x == nullptr && v != nullptr), "Passing x==NULL and v!=NULL is not supported");
2742 TopOnlyOK = (ir == nullptr);
2744 int fileVersion; /* Version number of the code that wrote the file */
2745 int fileGeneration; /* Generation version number of the code that wrote the file */
2746 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &fileVersion, &fileGeneration);
2748 if (bRead)
2750 state->flags = 0;
2751 init_gtc_state(state, tpx.ngtc, 0, 0);
2752 if (x == nullptr)
2754 // v is also nullptr by the above assertion, so we may
2755 // need to make memory in state for storing the contents
2756 // of the tpx file.
2757 if (tpx.bX)
2759 state->flags |= (1 << estX);
2761 if (tpx.bV)
2763 state->flags |= (1 << estV);
2765 state_change_natoms(state, tpx.natoms);
2769 if (x == nullptr)
2771 x = as_rvec_array(state->x.data());
2772 v = as_rvec_array(state->v.data());
2775 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
2777 do_test(fio, tpx.bBox, state->box);
2778 if (tpx.bBox)
2780 gmx_fio_ndo_rvec(fio, state->box, DIM);
2781 if (fileVersion >= 51)
2783 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
2785 else
2787 /* We initialize box_rel after reading the inputrec */
2788 clear_mat(state->box_rel);
2790 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
2791 if (fileVersion < 56)
2793 matrix mdum;
2794 gmx_fio_ndo_rvec(fio, mdum, DIM);
2798 if (state->ngtc > 0)
2800 real *dumv;
2801 snew(dumv, state->ngtc);
2802 if (fileVersion < 69)
2804 gmx_fio_ndo_real(fio, dumv, state->ngtc);
2806 /* These used to be the Berendsen tcoupl_lambda's */
2807 gmx_fio_ndo_real(fio, dumv, state->ngtc);
2808 sfree(dumv);
2811 do_test(fio, tpx.bTop, mtop);
2812 if (tpx.bTop)
2814 if (mtop)
2816 do_mtop(fio, mtop, bRead, fileVersion);
2818 else
2820 gmx_mtop_t dum_top;
2821 do_mtop(fio, &dum_top, bRead, fileVersion);
2824 do_test(fio, tpx.bX, x);
2825 if (tpx.bX)
2827 if (bRead)
2829 state->flags |= (1<<estX);
2831 gmx_fio_ndo_rvec(fio, x, tpx.natoms);
2834 do_test(fio, tpx.bV, v);
2835 if (tpx.bV)
2837 if (bRead)
2839 state->flags |= (1<<estV);
2841 gmx_fio_ndo_rvec(fio, v, tpx.natoms);
2844 // No need to run do_test when the last argument is NULL
2845 if (tpx.bF)
2847 rvec *dummyForces;
2848 snew(dummyForces, state->natoms);
2849 gmx_fio_ndo_rvec(fio, dummyForces, tpx.natoms);
2850 sfree(dummyForces);
2853 /* Starting with tpx version 26, we have the inputrec
2854 * at the end of the file, so we can ignore it
2855 * if the file is never than the software (but still the
2856 * same generation - see comments at the top of this file.
2860 ePBC = -1;
2861 bPeriodicMols = FALSE;
2863 do_test(fio, tpx.bIr, ir);
2864 if (tpx.bIr)
2866 if (fileVersion >= 53)
2868 /* Removed the pbc info from do_inputrec, since we always want it */
2869 if (!bRead)
2871 ePBC = ir->ePBC;
2872 bPeriodicMols = ir->bPeriodicMols;
2874 gmx_fio_do_int(fio, ePBC);
2875 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
2877 if (fileGeneration <= tpx_generation && ir)
2879 do_inputrec(fio, ir, bRead, fileVersion);
2880 if (fileVersion < 51)
2882 set_box_rel(ir, state);
2884 if (fileVersion < 53)
2886 ePBC = ir->ePBC;
2887 bPeriodicMols = ir->bPeriodicMols;
2890 if (bRead && ir && fileVersion >= 53)
2892 /* We need to do this after do_inputrec, since that initializes ir */
2893 ir->ePBC = ePBC;
2894 ir->bPeriodicMols = bPeriodicMols;
2898 if (bRead)
2900 if (tpx.bIr && ir)
2902 if (state->ngtc == 0)
2904 /* Reading old version without tcoupl state data: set it */
2905 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
2907 if (tpx.bTop && mtop)
2909 if (fileVersion < 57)
2911 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
2913 ir->eDisre = edrSimple;
2915 else
2917 ir->eDisre = edrNone;
2920 set_disres_npair(mtop);
2924 if (tpx.bTop && mtop)
2926 gmx_mtop_finalize(mtop);
2930 return ePBC;
2933 static t_fileio *open_tpx(const char *fn, const char *mode)
2935 return gmx_fio_open(fn, mode);
2938 static void close_tpx(t_fileio *fio)
2940 gmx_fio_close(fio);
2943 /************************************************************
2945 * The following routines are the exported ones
2947 ************************************************************/
2949 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK)
2951 t_fileio *fio;
2953 fio = open_tpx(fn, "r");
2954 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, nullptr, nullptr);
2955 close_tpx(fio);
2958 void write_tpx_state(const char *fn,
2959 const t_inputrec *ir, const t_state *state, const gmx_mtop_t *mtop)
2961 t_fileio *fio;
2963 fio = open_tpx(fn, "w");
2964 do_tpx(fio, FALSE,
2965 const_cast<t_inputrec *>(ir),
2966 const_cast<t_state *>(state), nullptr, nullptr,
2967 const_cast<gmx_mtop_t *>(mtop));
2968 close_tpx(fio);
2971 void read_tpx_state(const char *fn,
2972 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
2974 t_fileio *fio;
2976 fio = open_tpx(fn, "r");
2977 do_tpx(fio, TRUE, ir, state, nullptr, nullptr, mtop);
2978 close_tpx(fio);
2981 int read_tpx(const char *fn,
2982 t_inputrec *ir, matrix box, int *natoms,
2983 rvec *x, rvec *v, gmx_mtop_t *mtop)
2985 t_fileio *fio;
2986 t_state state;
2987 int ePBC;
2989 fio = open_tpx(fn, "r");
2990 ePBC = do_tpx(fio, TRUE, ir, &state, x, v, mtop);
2991 close_tpx(fio);
2992 *natoms = mtop->natoms;
2993 if (box)
2995 copy_mat(state.box, box);
2998 return ePBC;
3001 int read_tpx_top(const char *fn,
3002 t_inputrec *ir, matrix box, int *natoms,
3003 rvec *x, rvec *v, t_topology *top)
3005 gmx_mtop_t mtop;
3006 int ePBC;
3008 ePBC = read_tpx(fn, ir, box, natoms, x, v, &mtop);
3010 *top = gmx_mtop_t_to_t_topology(&mtop, true);
3012 return ePBC;
3015 gmx_bool fn2bTPX(const char *file)
3017 return (efTPR == fn2ftp(file));
3020 void pr_tpxheader(FILE *fp, int indent, const char *title, const t_tpxheader *sh)
3022 if (available(fp, sh, indent, title))
3024 indent = pr_title(fp, indent, title);
3025 pr_indent(fp, indent);
3026 fprintf(fp, "bIr = %spresent\n", sh->bIr ? "" : "not ");
3027 pr_indent(fp, indent);
3028 fprintf(fp, "bBox = %spresent\n", sh->bBox ? "" : "not ");
3029 pr_indent(fp, indent);
3030 fprintf(fp, "bTop = %spresent\n", sh->bTop ? "" : "not ");
3031 pr_indent(fp, indent);
3032 fprintf(fp, "bX = %spresent\n", sh->bX ? "" : "not ");
3033 pr_indent(fp, indent);
3034 fprintf(fp, "bV = %spresent\n", sh->bV ? "" : "not ");
3035 pr_indent(fp, indent);
3036 fprintf(fp, "bF = %spresent\n", sh->bF ? "" : "not ");
3038 pr_indent(fp, indent);
3039 fprintf(fp, "natoms = %d\n", sh->natoms);
3040 pr_indent(fp, indent);
3041 fprintf(fp, "lambda = %e\n", sh->lambda);