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.
39 /* This file is completely threadsafe - keep it that way! */
50 #include "gromacs/compat/make_unique.h"
51 #include "gromacs/fileio/filetypes.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/gmxfio-xdr.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdtypes/awh-history.h"
57 #include "gromacs/mdtypes/awh-params.h"
58 #include "gromacs/mdtypes/inputrec.h"
59 #include "gromacs/mdtypes/md_enums.h"
60 #include "gromacs/mdtypes/pull-params.h"
61 #include "gromacs/mdtypes/state.h"
62 #include "gromacs/pbcutil/boxutilities.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/topology/block.h"
65 #include "gromacs/topology/ifunc.h"
66 #include "gromacs/topology/mtop_util.h"
67 #include "gromacs/topology/symtab.h"
68 #include "gromacs/topology/topology.h"
69 #include "gromacs/utility/arraysize.h"
70 #include "gromacs/utility/baseversion.h"
71 #include "gromacs/utility/cstringutil.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/futil.h"
74 #include "gromacs/utility/gmxassert.h"
75 #include "gromacs/utility/keyvaluetreebuilder.h"
76 #include "gromacs/utility/keyvaluetreeserializer.h"
77 #include "gromacs/utility/smalloc.h"
78 #include "gromacs/utility/snprintf.h"
79 #include "gromacs/utility/txtdump.h"
81 #define TPX_TAG_RELEASE "release"
83 /*! \brief Tag string for the file format written to run input files
84 * written by this version of the code.
86 * Change this if you want to change the run input format in a feature
87 * branch. This ensures that there will not be different run input
88 * formats around which cannot be distinguished, while not causing
89 * problems rebasing the feature branch onto upstream changes. When
90 * merging with mainstream GROMACS, set this tag string back to
91 * TPX_TAG_RELEASE, and instead add an element to tpxv.
93 static const char *tpx_tag = TPX_TAG_RELEASE;
95 /*! \brief Enum of values that describe the contents of a tpr file
96 * whose format matches a version number
98 * The enum helps the code be more self-documenting and ensure merges
99 * do not silently resolve when two patches make the same bump. When
100 * adding new functionality, add a new element just above tpxv_Count
101 * in this enumeration, and write code below that does the right thing
102 * according to the value of file_version.
105 tpxv_ComputationalElectrophysiology = 96, /**< support for ion/water position swaps (computational electrophysiology) */
106 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to int64_t */
107 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
108 tpxv_InteractiveMolecularDynamics, /**< interactive molecular dynamics (IMD) */
109 tpxv_RemoveObsoleteParameters1, /**< remove optimize_fft, dihre_fc, nstcheckpoint */
110 tpxv_PullCoordTypeGeom, /**< add pull type and geometry per group and flat-bottom */
111 tpxv_PullGeomDirRel, /**< add pull geometry direction-relative */
112 tpxv_IntermolecularBondeds, /**< permit inter-molecular bonded interactions in the topology */
113 tpxv_CompElWithSwapLayerOffset, /**< added parameters for improved CompEl setups */
114 tpxv_CompElPolyatomicIonsAndMultipleIonTypes, /**< CompEl now can handle polyatomic ions and more than two types of ions */
115 tpxv_RemoveAdress, /**< removed support for AdResS */
116 tpxv_PullCoordNGroup, /**< add ngroup to pull coord */
117 tpxv_RemoveTwinRange, /**< removed support for twin-range interactions */
118 tpxv_ReplacePullPrintCOM12, /**< Replaced print-com-1, 2 with pull-print-com */
119 tpxv_PullExternalPotential, /**< Added pull type external potential */
120 tpxv_GenericParamsForElectricField, /**< Introduced KeyValueTree and moved electric field parameters */
121 tpxv_AcceleratedWeightHistogram, /**< sampling with accelerated weight histogram method (AWH) */
122 tpxv_RemoveImplicitSolvation, /**< removed support for implicit solvation */
123 tpxv_Count /**< the total number of tpxv versions */
126 /*! \brief Version number of the file format written to run input
127 * files by this version of the code.
129 * The tpx_version increases whenever the file format in the main
130 * development branch changes, due to an extension of the tpxv enum above.
131 * Backward compatibility for reading old run input files is maintained
132 * by checking this version number against that of the file and then using
133 * the correct code path.
135 * When developing a feature branch that needs to change the run input
136 * file format, change tpx_tag instead. */
137 static const int tpx_version = tpxv_Count - 1;
140 /* This number should only be increased when you edit the TOPOLOGY section
141 * or the HEADER of the tpx format.
142 * This way we can maintain forward compatibility too for all analysis tools
143 * and/or external programs that only need to know the atom/residue names,
144 * charges, and bond connectivity.
146 * It first appeared in tpx version 26, when I also moved the inputrecord
147 * to the end of the tpx file, so we can just skip it if we only
150 * In particular, it must be increased when adding new elements to
151 * ftupd, so that old code can read new .tpr files.
153 static const int tpx_generation = 26;
155 /* This number should be the most recent backwards incompatible version
156 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
158 static const int tpx_incompatible_version = 57; // GMX4.0 has version 58
162 /* Struct used to maintain tpx compatibility when function types are added */
164 int fvnr; /* file version number in which the function type first appeared */
165 int ftype; /* function type */
169 * TODO The following three lines make little sense, please clarify if
170 * you've had to work out how ftupd works.
172 * The entries should be ordered in:
173 * 1. ascending function type number
174 * 2. ascending file version number
176 * Because we support reading of old .tpr file versions (even when
177 * mdrun can no longer run the simulation), we need to be able to read
178 * obsolete t_interaction_function types. Any data read from such
179 * fields is discarded. Their names have _NOLONGERUSED appended to
180 * them to make things clear.
182 static const t_ftupd ftupd[] = {
183 { 70, F_RESTRBONDS },
184 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
185 { 76, F_LINEAR_ANGLES },
186 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
187 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
189 { 60, F_GB12_NOLONGERUSED },
190 { 61, F_GB13_NOLONGERUSED },
191 { 61, F_GB14_NOLONGERUSED },
192 { 72, F_GBPOL_NOLONGERUSED },
193 { 72, F_NPSOLVATION_NOLONGERUSED },
196 { 69, F_VTEMP_NOLONGERUSED},
198 { 76, F_ANHARM_POL },
201 { 79, F_DVDL_BONDED, },
202 { 79, F_DVDL_RESTRAINT },
203 { 79, F_DVDL_TEMPERATURE },
205 #define NFTUPD asize(ftupd)
207 /* Needed for backward compatibility */
210 /**************************************************************
212 * Now the higer level routines that do io of the structures and arrays
214 **************************************************************/
215 static void do_pullgrp_tpx_pre95(t_fileio *fio,
222 gmx_fio_do_int(fio, pgrp->nat);
225 snew(pgrp->ind, pgrp->nat);
227 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
228 gmx_fio_do_int(fio, pgrp->nweight);
231 snew(pgrp->weight, pgrp->nweight);
233 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
234 gmx_fio_do_int(fio, pgrp->pbcatom);
235 gmx_fio_do_rvec(fio, pcrd->vec);
236 clear_rvec(pcrd->origin);
237 gmx_fio_do_rvec(fio, tmp);
239 gmx_fio_do_real(fio, pcrd->rate);
240 gmx_fio_do_real(fio, pcrd->k);
241 gmx_fio_do_real(fio, pcrd->kB);
244 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
246 gmx_fio_do_int(fio, pgrp->nat);
249 snew(pgrp->ind, pgrp->nat);
251 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
252 gmx_fio_do_int(fio, pgrp->nweight);
255 snew(pgrp->weight, pgrp->nweight);
257 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
258 gmx_fio_do_int(fio, pgrp->pbcatom);
261 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd,
262 gmx_bool bRead, int file_version,
263 int ePullOld, int eGeomOld, ivec dimOld)
265 if (file_version >= tpxv_PullCoordNGroup)
267 gmx_fio_do_int(fio, pcrd->eType);
268 if (file_version >= tpxv_PullExternalPotential)
270 if (pcrd->eType == epullEXTERNAL)
276 gmx_fio_do_string(fio, buf);
277 pcrd->externalPotentialProvider = gmx_strdup(buf);
281 gmx_fio_do_string(fio, pcrd->externalPotentialProvider);
286 pcrd->externalPotentialProvider = nullptr;
293 pcrd->externalPotentialProvider = nullptr;
296 /* Note that we try to support adding new geometries without
297 * changing the tpx version. This requires checks when printing the
298 * geometry string and a check and fatal_error in init_pull.
300 gmx_fio_do_int(fio, pcrd->eGeom);
301 gmx_fio_do_int(fio, pcrd->ngroup);
302 if (pcrd->ngroup <= c_pullCoordNgroupMax)
304 gmx_fio_ndo_int(fio, pcrd->group, pcrd->ngroup);
308 /* More groups in file than supported, this must be a new geometry
309 * that is not supported by our current code. Since we will not
310 * use the groups for this coord (checks in the pull and WHAM code
311 * ensure this), we can ignore the groups and set ngroup=0.
314 snew(dum, pcrd->ngroup);
315 gmx_fio_ndo_int(fio, dum, pcrd->ngroup);
320 gmx_fio_do_ivec(fio, pcrd->dim);
325 gmx_fio_do_int(fio, pcrd->group[0]);
326 gmx_fio_do_int(fio, pcrd->group[1]);
327 if (file_version >= tpxv_PullCoordTypeGeom)
329 pcrd->ngroup = (pcrd->eGeom == epullgDIRRELATIVE ? 4 : 2);
330 gmx_fio_do_int(fio, pcrd->eType);
331 gmx_fio_do_int(fio, pcrd->eGeom);
332 if (pcrd->ngroup == 4)
334 gmx_fio_do_int(fio, pcrd->group[2]);
335 gmx_fio_do_int(fio, pcrd->group[3]);
337 gmx_fio_do_ivec(fio, pcrd->dim);
341 pcrd->eType = ePullOld;
342 pcrd->eGeom = eGeomOld;
343 copy_ivec(dimOld, pcrd->dim);
346 gmx_fio_do_rvec(fio, pcrd->origin);
347 gmx_fio_do_rvec(fio, pcrd->vec);
348 if (file_version >= tpxv_PullCoordTypeGeom)
350 gmx_fio_do_gmx_bool(fio, pcrd->bStart);
354 /* This parameter is only printed, but not actually used by mdrun */
355 pcrd->bStart = FALSE;
357 gmx_fio_do_real(fio, pcrd->init);
358 gmx_fio_do_real(fio, pcrd->rate);
359 gmx_fio_do_real(fio, pcrd->k);
360 gmx_fio_do_real(fio, pcrd->kB);
363 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
365 int n_lambda = fepvals->n_lambda;
367 /* reset the lambda calculation window */
368 fepvals->lambda_start_n = 0;
369 fepvals->lambda_stop_n = n_lambda;
370 if (file_version >= 79)
376 snew(expand->init_lambda_weights, n_lambda);
378 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
379 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
382 gmx_fio_do_int(fio, expand->nstexpanded);
383 gmx_fio_do_int(fio, expand->elmcmove);
384 gmx_fio_do_int(fio, expand->elamstats);
385 gmx_fio_do_int(fio, expand->lmc_repeats);
386 gmx_fio_do_int(fio, expand->gibbsdeltalam);
387 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
388 gmx_fio_do_int(fio, expand->lmc_seed);
389 gmx_fio_do_real(fio, expand->mc_temp);
390 gmx_fio_do_gmx_bool(fio, expand->bSymmetrizedTMatrix);
391 gmx_fio_do_int(fio, expand->nstTij);
392 gmx_fio_do_int(fio, expand->minvarmin);
393 gmx_fio_do_int(fio, expand->c_range);
394 gmx_fio_do_real(fio, expand->wl_scale);
395 gmx_fio_do_real(fio, expand->wl_ratio);
396 gmx_fio_do_real(fio, expand->init_wl_delta);
397 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
398 gmx_fio_do_int(fio, expand->elmceq);
399 gmx_fio_do_int(fio, expand->equil_steps);
400 gmx_fio_do_int(fio, expand->equil_samples);
401 gmx_fio_do_int(fio, expand->equil_n_at_lam);
402 gmx_fio_do_real(fio, expand->equil_wl_delta);
403 gmx_fio_do_real(fio, expand->equil_ratio);
407 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
410 if (file_version >= 79)
412 gmx_fio_do_int(fio, simtemp->eSimTempScale);
413 gmx_fio_do_real(fio, simtemp->simtemp_high);
414 gmx_fio_do_real(fio, simtemp->simtemp_low);
419 snew(simtemp->temperatures, n_lambda);
421 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
426 static void do_imd(t_fileio *fio, t_IMD *imd, gmx_bool bRead)
428 gmx_fio_do_int(fio, imd->nat);
431 snew(imd->ind, imd->nat);
433 gmx_fio_ndo_int(fio, imd->ind, imd->nat);
436 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
438 /* i is defined in the ndo_double macro; use g to iterate. */
442 /* free energy values */
444 if (file_version >= 79)
446 gmx_fio_do_int(fio, fepvals->init_fep_state);
447 gmx_fio_do_double(fio, fepvals->init_lambda);
448 gmx_fio_do_double(fio, fepvals->delta_lambda);
450 else if (file_version >= 59)
452 gmx_fio_do_double(fio, fepvals->init_lambda);
453 gmx_fio_do_double(fio, fepvals->delta_lambda);
457 gmx_fio_do_real(fio, rdum);
458 fepvals->init_lambda = rdum;
459 gmx_fio_do_real(fio, rdum);
460 fepvals->delta_lambda = rdum;
462 if (file_version >= 79)
464 gmx_fio_do_int(fio, fepvals->n_lambda);
467 snew(fepvals->all_lambda, efptNR);
469 for (g = 0; g < efptNR; g++)
471 if (fepvals->n_lambda > 0)
475 snew(fepvals->all_lambda[g], fepvals->n_lambda);
477 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
478 gmx_fio_ndo_gmx_bool(fio, fepvals->separate_dvdl, efptNR);
480 else if (fepvals->init_lambda >= 0)
482 fepvals->separate_dvdl[efptFEP] = TRUE;
486 else if (file_version >= 64)
488 gmx_fio_do_int(fio, fepvals->n_lambda);
493 snew(fepvals->all_lambda, efptNR);
494 /* still allocate the all_lambda array's contents. */
495 for (g = 0; g < efptNR; g++)
497 if (fepvals->n_lambda > 0)
499 snew(fepvals->all_lambda[g], fepvals->n_lambda);
503 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
505 if (fepvals->init_lambda >= 0)
509 fepvals->separate_dvdl[efptFEP] = TRUE;
513 /* copy the contents of the efptFEP lambda component to all
514 the other components */
515 for (g = 0; g < efptNR; g++)
517 for (h = 0; h < fepvals->n_lambda; h++)
521 fepvals->all_lambda[g][h] =
522 fepvals->all_lambda[efptFEP][h];
531 fepvals->n_lambda = 0;
532 fepvals->all_lambda = nullptr;
533 if (fepvals->init_lambda >= 0)
535 fepvals->separate_dvdl[efptFEP] = TRUE;
538 gmx_fio_do_real(fio, fepvals->sc_alpha);
539 gmx_fio_do_int(fio, fepvals->sc_power);
540 if (file_version >= 79)
542 gmx_fio_do_real(fio, fepvals->sc_r_power);
546 fepvals->sc_r_power = 6.0;
548 gmx_fio_do_real(fio, fepvals->sc_sigma);
551 if (file_version >= 71)
553 fepvals->sc_sigma_min = fepvals->sc_sigma;
557 fepvals->sc_sigma_min = 0;
560 if (file_version >= 79)
562 gmx_fio_do_gmx_bool(fio, fepvals->bScCoul);
566 fepvals->bScCoul = TRUE;
568 if (file_version >= 64)
570 gmx_fio_do_int(fio, fepvals->nstdhdl);
574 fepvals->nstdhdl = 1;
577 if (file_version >= 73)
579 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
580 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
584 fepvals->separate_dhdl_file = esepdhdlfileYES;
585 fepvals->dhdl_derivatives = edhdlderivativesYES;
587 if (file_version >= 71)
589 gmx_fio_do_int(fio, fepvals->dh_hist_size);
590 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
594 fepvals->dh_hist_size = 0;
595 fepvals->dh_hist_spacing = 0.1;
597 if (file_version >= 79)
599 gmx_fio_do_int(fio, fepvals->edHdLPrintEnergy);
603 fepvals->edHdLPrintEnergy = edHdLPrintEnergyNO;
606 /* handle lambda_neighbors */
607 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
609 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
610 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
611 (fepvals->init_lambda < 0) )
613 fepvals->lambda_start_n = (fepvals->init_fep_state -
614 fepvals->lambda_neighbors);
615 fepvals->lambda_stop_n = (fepvals->init_fep_state +
616 fepvals->lambda_neighbors + 1);
617 if (fepvals->lambda_start_n < 0)
619 fepvals->lambda_start_n = 0;;
621 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
623 fepvals->lambda_stop_n = fepvals->n_lambda;
628 fepvals->lambda_start_n = 0;
629 fepvals->lambda_stop_n = fepvals->n_lambda;
634 fepvals->lambda_start_n = 0;
635 fepvals->lambda_stop_n = fepvals->n_lambda;
639 static void do_awhBias(t_fileio *fio, gmx::AwhBiasParams *awhBiasParams, gmx_bool bRead,
640 int gmx_unused file_version)
642 gmx_fio_do_int(fio, awhBiasParams->eTarget);
643 gmx_fio_do_double(fio, awhBiasParams->targetBetaScaling);
644 gmx_fio_do_double(fio, awhBiasParams->targetCutoff);
645 gmx_fio_do_int(fio, awhBiasParams->eGrowth);
646 gmx_fio_do_int(fio, awhBiasParams->bUserData);
647 gmx_fio_do_double(fio, awhBiasParams->errorInitial);
648 gmx_fio_do_int(fio, awhBiasParams->ndim);
649 gmx_fio_do_int(fio, awhBiasParams->shareGroup);
650 gmx_fio_do_gmx_bool(fio, awhBiasParams->equilibrateHistogram);
654 snew(awhBiasParams->dimParams, awhBiasParams->ndim);
657 for (int d = 0; d < awhBiasParams->ndim; d++)
659 gmx::AwhDimParams *dimParams = &awhBiasParams->dimParams[d];
661 gmx_fio_do_int(fio, dimParams->eCoordProvider);
662 gmx_fio_do_int(fio, dimParams->coordIndex);
663 gmx_fio_do_double(fio, dimParams->origin);
664 gmx_fio_do_double(fio, dimParams->end);
665 gmx_fio_do_double(fio, dimParams->period);
666 gmx_fio_do_double(fio, dimParams->forceConstant);
667 gmx_fio_do_double(fio, dimParams->diffusion);
668 gmx_fio_do_double(fio, dimParams->coordValueInit);
669 gmx_fio_do_double(fio, dimParams->coverDiameter);
673 static void do_awh(t_fileio *fio, gmx::AwhParams *awhParams, gmx_bool bRead,
674 int gmx_unused file_version)
676 gmx_fio_do_int(fio, awhParams->numBias);
677 gmx_fio_do_int(fio, awhParams->nstOut);
678 gmx_fio_do_int64(fio, awhParams->seed);
679 gmx_fio_do_int(fio, awhParams->nstSampleCoord);
680 gmx_fio_do_int(fio, awhParams->numSamplesUpdateFreeEnergy);
681 gmx_fio_do_int(fio, awhParams->ePotential);
682 gmx_fio_do_gmx_bool(fio, awhParams->shareBiasMultisim);
684 if (awhParams->numBias > 0)
688 snew(awhParams->awhBiasParams, awhParams->numBias);
691 for (int k = 0; k < awhParams->numBias; k++)
693 do_awhBias(fio, &awhParams->awhBiasParams[k], bRead, file_version);
698 static void do_pull(t_fileio *fio, pull_params_t *pull, gmx_bool bRead,
699 int file_version, int ePullOld)
705 if (file_version >= 95)
707 gmx_fio_do_int(fio, pull->ngroup);
709 gmx_fio_do_int(fio, pull->ncoord);
710 if (file_version < 95)
712 pull->ngroup = pull->ncoord + 1;
714 if (file_version < tpxv_PullCoordTypeGeom)
718 gmx_fio_do_int(fio, eGeomOld);
719 gmx_fio_do_ivec(fio, dimOld);
720 /* The inner cylinder radius, now removed */
721 gmx_fio_do_real(fio, dum);
723 gmx_fio_do_real(fio, pull->cylinder_r);
724 gmx_fio_do_real(fio, pull->constr_tol);
725 if (file_version >= 95)
727 gmx_fio_do_gmx_bool(fio, pull->bPrintCOM);
728 /* With file_version < 95 this value is set below */
730 if (file_version >= tpxv_ReplacePullPrintCOM12)
732 gmx_fio_do_gmx_bool(fio, pull->bPrintRefValue);
733 gmx_fio_do_gmx_bool(fio, pull->bPrintComp);
735 else if (file_version >= tpxv_PullCoordTypeGeom)
738 gmx_fio_do_int(fio, idum); /* used to be bPrintCOM2 */
739 gmx_fio_do_gmx_bool(fio, pull->bPrintRefValue);
740 gmx_fio_do_gmx_bool(fio, pull->bPrintComp);
744 pull->bPrintRefValue = FALSE;
745 pull->bPrintComp = TRUE;
747 gmx_fio_do_int(fio, pull->nstxout);
748 gmx_fio_do_int(fio, pull->nstfout);
751 snew(pull->group, pull->ngroup);
752 snew(pull->coord, pull->ncoord);
754 if (file_version < 95)
756 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
757 if (eGeomOld == epullgDIRPBC)
759 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
761 if (eGeomOld > epullgDIRPBC)
766 for (g = 0; g < pull->ngroup; g++)
768 /* We read and ignore a pull coordinate for group 0 */
769 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[std::max(g-1, 0)],
773 pull->coord[g-1].group[0] = 0;
774 pull->coord[g-1].group[1] = g;
778 pull->bPrintCOM = (pull->group[0].nat > 0);
782 for (g = 0; g < pull->ngroup; g++)
784 do_pull_group(fio, &pull->group[g], bRead);
786 for (g = 0; g < pull->ncoord; g++)
788 do_pull_coord(fio, &pull->coord[g],
789 bRead, file_version, ePullOld, eGeomOld, dimOld);
795 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
797 gmx_fio_do_int(fio, rotg->eType);
798 gmx_fio_do_int(fio, rotg->bMassW);
799 gmx_fio_do_int(fio, rotg->nat);
802 snew(rotg->ind, rotg->nat);
804 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
807 snew(rotg->x_ref, rotg->nat);
809 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
810 gmx_fio_do_rvec(fio, rotg->inputVec);
811 gmx_fio_do_rvec(fio, rotg->pivot);
812 gmx_fio_do_real(fio, rotg->rate);
813 gmx_fio_do_real(fio, rotg->k);
814 gmx_fio_do_real(fio, rotg->slab_dist);
815 gmx_fio_do_real(fio, rotg->min_gaussian);
816 gmx_fio_do_real(fio, rotg->eps);
817 gmx_fio_do_int(fio, rotg->eFittype);
818 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
819 gmx_fio_do_real(fio, rotg->PotAngle_step);
822 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
826 gmx_fio_do_int(fio, rot->ngrp);
827 gmx_fio_do_int(fio, rot->nstrout);
828 gmx_fio_do_int(fio, rot->nstsout);
831 snew(rot->grp, rot->ngrp);
833 for (g = 0; g < rot->ngrp; g++)
835 do_rotgrp(fio, &rot->grp[g], bRead);
840 static void do_swapgroup(t_fileio *fio, t_swapGroup *g, gmx_bool bRead)
843 /* Name of the group or molecule */
848 gmx_fio_do_string(fio, buf);
849 g->molname = gmx_strdup(buf);
853 gmx_fio_do_string(fio, g->molname);
856 /* Number of atoms in the group */
857 gmx_fio_do_int(fio, g->nat);
859 /* The group's atom indices */
862 snew(g->ind, g->nat);
864 gmx_fio_ndo_int(fio, g->ind, g->nat);
866 /* Requested counts for compartments A and B */
867 gmx_fio_ndo_int(fio, g->nmolReq, eCompNR);
870 static void do_swapcoords_tpx(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead, int file_version)
872 /* Enums for better readability of the code */
877 eChannel0 = 0, eChannel1
881 if (file_version >= tpxv_CompElPolyatomicIonsAndMultipleIonTypes)
883 /* The total number of swap groups is the sum of the fixed groups
884 * (split0, split1, solvent), and the user-defined groups (2+ types of ions)
886 gmx_fio_do_int(fio, swap->ngrp);
889 snew(swap->grp, swap->ngrp);
891 for (int ig = 0; ig < swap->ngrp; ig++)
893 do_swapgroup(fio, &swap->grp[ig], bRead);
895 gmx_fio_do_gmx_bool(fio, swap->massw_split[eChannel0]);
896 gmx_fio_do_gmx_bool(fio, swap->massw_split[eChannel1]);
897 gmx_fio_do_int(fio, swap->nstswap);
898 gmx_fio_do_int(fio, swap->nAverage);
899 gmx_fio_do_real(fio, swap->threshold);
900 gmx_fio_do_real(fio, swap->cyl0r);
901 gmx_fio_do_real(fio, swap->cyl0u);
902 gmx_fio_do_real(fio, swap->cyl0l);
903 gmx_fio_do_real(fio, swap->cyl1r);
904 gmx_fio_do_real(fio, swap->cyl1u);
905 gmx_fio_do_real(fio, swap->cyl1l);
909 /*** Support reading older CompEl .tpr files ***/
911 /* In the original CompEl .tpr files, we always have 5 groups: */
913 snew(swap->grp, swap->ngrp);
915 swap->grp[eGrpSplit0 ].molname = gmx_strdup("split0" ); // group 0: split0
916 swap->grp[eGrpSplit1 ].molname = gmx_strdup("split1" ); // group 1: split1
917 swap->grp[eGrpSolvent].molname = gmx_strdup("solvent"); // group 2: solvent
918 swap->grp[3 ].molname = gmx_strdup("anions" ); // group 3: anions
919 swap->grp[4 ].molname = gmx_strdup("cations"); // group 4: cations
921 gmx_fio_do_int(fio, swap->grp[3].nat);
922 gmx_fio_do_int(fio, swap->grp[eGrpSolvent].nat);
923 gmx_fio_do_int(fio, swap->grp[eGrpSplit0].nat);
924 gmx_fio_do_gmx_bool(fio, swap->massw_split[eChannel0]);
925 gmx_fio_do_int(fio, swap->grp[eGrpSplit1].nat);
926 gmx_fio_do_gmx_bool(fio, swap->massw_split[eChannel1]);
927 gmx_fio_do_int(fio, swap->nstswap);
928 gmx_fio_do_int(fio, swap->nAverage);
929 gmx_fio_do_real(fio, swap->threshold);
930 gmx_fio_do_real(fio, swap->cyl0r);
931 gmx_fio_do_real(fio, swap->cyl0u);
932 gmx_fio_do_real(fio, swap->cyl0l);
933 gmx_fio_do_real(fio, swap->cyl1r);
934 gmx_fio_do_real(fio, swap->cyl1u);
935 gmx_fio_do_real(fio, swap->cyl1l);
937 // The order[] array keeps compatibility with older .tpr files
938 // by reading in the groups in the classic order
940 const int order[4] = {3, eGrpSolvent, eGrpSplit0, eGrpSplit1 };
942 for (int ig = 0; ig < 4; ig++)
945 snew(swap->grp[g].ind, swap->grp[g].nat);
946 gmx_fio_ndo_int(fio, swap->grp[g].ind, swap->grp[g].nat);
950 for (int j = eCompA; j <= eCompB; j++)
952 gmx_fio_do_int(fio, swap->grp[3].nmolReq[j]); // group 3 = anions
953 gmx_fio_do_int(fio, swap->grp[4].nmolReq[j]); // group 4 = cations
955 } /* End support reading older CompEl .tpr files */
957 if (file_version >= tpxv_CompElWithSwapLayerOffset)
959 gmx_fio_do_real(fio, swap->bulkOffset[eCompA]);
960 gmx_fio_do_real(fio, swap->bulkOffset[eCompB]);
965 static void do_legacy_efield(t_fileio *fio, gmx::KeyValueTreeObjectBuilder *root)
967 const char *const dimName[] = { "x", "y", "z" };
969 auto appliedForcesObj = root->addObject("applied-forces");
970 auto efieldObj = appliedForcesObj.addObject("electric-field");
971 // The content of the tpr file for this feature has
972 // been the same since gromacs 4.0 that was used for
974 for (int j = 0; j < DIM; ++j)
977 gmx_fio_do_int(fio, n);
978 gmx_fio_do_int(fio, nt);
979 std::vector<real> aa(n+1), phi(nt+1), at(nt+1), phit(nt+1);
980 gmx_fio_ndo_real(fio, aa.data(), n);
981 gmx_fio_ndo_real(fio, phi.data(), n);
982 gmx_fio_ndo_real(fio, at.data(), nt);
983 gmx_fio_ndo_real(fio, phit.data(), nt);
988 gmx_fatal(FARGS, "Can not handle tpr files with more than one electric field term per direction.");
990 auto dimObj = efieldObj.addObject(dimName[j]);
991 dimObj.addValue<real>("E0", aa[0]);
992 dimObj.addValue<real>("omega", at[0]);
993 dimObj.addValue<real>("t0", phi[0]);
994 dimObj.addValue<real>("sigma", phit[0]);
1000 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
1003 int i, j, k, idum = 0;
1005 gmx_bool bdum = false;
1007 if (file_version != tpx_version)
1009 /* Give a warning about features that are not accessible */
1010 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
1011 file_version, tpx_version);
1014 if (file_version == 0)
1019 gmx::KeyValueTreeBuilder paramsBuilder;
1020 gmx::KeyValueTreeObjectBuilder paramsObj = paramsBuilder.rootObject();
1022 /* Basic inputrec stuff */
1023 gmx_fio_do_int(fio, ir->eI);
1024 if (file_version >= 62)
1026 gmx_fio_do_int64(fio, ir->nsteps);
1030 gmx_fio_do_int(fio, idum);
1034 if (file_version >= 62)
1036 gmx_fio_do_int64(fio, ir->init_step);
1040 gmx_fio_do_int(fio, idum);
1041 ir->init_step = idum;
1044 gmx_fio_do_int(fio, ir->simulation_part);
1046 if (file_version >= 67)
1048 gmx_fio_do_int(fio, ir->nstcalcenergy);
1052 ir->nstcalcenergy = 1;
1054 if (file_version >= 81)
1056 gmx_fio_do_int(fio, ir->cutoff_scheme);
1057 if (file_version < 94)
1059 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
1064 ir->cutoff_scheme = ecutsGROUP;
1066 gmx_fio_do_int(fio, ir->ns_type);
1067 gmx_fio_do_int(fio, ir->nstlist);
1068 gmx_fio_do_int(fio, idum); /* used to be ndelta; not used anymore */
1070 gmx_fio_do_real(fio, ir->rtpi);
1072 gmx_fio_do_int(fio, ir->nstcomm);
1073 gmx_fio_do_int(fio, ir->comm_mode);
1075 /* ignore nstcheckpoint */
1076 if (file_version < tpxv_RemoveObsoleteParameters1)
1078 gmx_fio_do_int(fio, idum);
1081 gmx_fio_do_int(fio, ir->nstcgsteep);
1083 gmx_fio_do_int(fio, ir->nbfgscorr);
1085 gmx_fio_do_int(fio, ir->nstlog);
1086 gmx_fio_do_int(fio, ir->nstxout);
1087 gmx_fio_do_int(fio, ir->nstvout);
1088 gmx_fio_do_int(fio, ir->nstfout);
1089 gmx_fio_do_int(fio, ir->nstenergy);
1090 gmx_fio_do_int(fio, ir->nstxout_compressed);
1091 if (file_version >= 59)
1093 gmx_fio_do_double(fio, ir->init_t);
1094 gmx_fio_do_double(fio, ir->delta_t);
1098 gmx_fio_do_real(fio, rdum);
1100 gmx_fio_do_real(fio, rdum);
1103 gmx_fio_do_real(fio, ir->x_compression_precision);
1104 if (file_version >= 81)
1106 gmx_fio_do_real(fio, ir->verletbuf_tol);
1110 ir->verletbuf_tol = 0;
1112 gmx_fio_do_real(fio, ir->rlist);
1113 if (file_version >= 67 && file_version < tpxv_RemoveTwinRange)
1117 // Reading such a file version could invoke the twin-range
1118 // scheme, about which mdrun should give a fatal error.
1119 real dummy_rlistlong = -1;
1120 gmx_fio_do_real(fio, dummy_rlistlong);
1122 if (ir->rlist > 0 && (dummy_rlistlong == 0 || dummy_rlistlong > ir->rlist))
1124 // Get mdrun to issue an error (regardless of
1125 // ir->cutoff_scheme).
1126 ir->useTwinRange = true;
1130 // grompp used to set rlistlong actively. Users were
1131 // probably also confused and set rlistlong == rlist.
1132 // However, in all remaining cases, it is safe to let
1133 // mdrun proceed normally.
1134 ir->useTwinRange = false;
1140 // No need to read or write anything
1141 ir->useTwinRange = false;
1143 if (file_version >= 82 && file_version != 90)
1145 // Multiple time-stepping is no longer enabled, but the old
1146 // support required the twin-range scheme, for which mdrun
1147 // already emits a fatal error.
1148 int dummy_nstcalclr = -1;
1149 gmx_fio_do_int(fio, dummy_nstcalclr);
1151 gmx_fio_do_int(fio, ir->coulombtype);
1152 if (file_version >= 81)
1154 gmx_fio_do_int(fio, ir->coulomb_modifier);
1158 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1160 gmx_fio_do_real(fio, ir->rcoulomb_switch);
1161 gmx_fio_do_real(fio, ir->rcoulomb);
1162 gmx_fio_do_int(fio, ir->vdwtype);
1163 if (file_version >= 81)
1165 gmx_fio_do_int(fio, ir->vdw_modifier);
1169 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1171 gmx_fio_do_real(fio, ir->rvdw_switch);
1172 gmx_fio_do_real(fio, ir->rvdw);
1173 gmx_fio_do_int(fio, ir->eDispCorr);
1174 gmx_fio_do_real(fio, ir->epsilon_r);
1175 gmx_fio_do_real(fio, ir->epsilon_rf);
1176 gmx_fio_do_real(fio, ir->tabext);
1178 // This permits reading a .tpr file that used implicit solvent,
1179 // and later permitting mdrun to refuse to run it.
1182 if (file_version < tpxv_RemoveImplicitSolvation)
1184 gmx_fio_do_int(fio, idum);
1185 gmx_fio_do_int(fio, idum);
1186 gmx_fio_do_real(fio, rdum);
1187 gmx_fio_do_real(fio, rdum);
1188 gmx_fio_do_int(fio, idum);
1189 ir->implicit_solvent = (idum > 0);
1193 ir->implicit_solvent = false;
1195 if (file_version < tpxv_RemoveImplicitSolvation)
1197 gmx_fio_do_real(fio, rdum);
1198 gmx_fio_do_real(fio, rdum);
1199 gmx_fio_do_real(fio, rdum);
1200 gmx_fio_do_real(fio, rdum);
1201 if (file_version >= 60)
1203 gmx_fio_do_real(fio, rdum);
1204 gmx_fio_do_int(fio, idum);
1206 gmx_fio_do_real(fio, rdum);
1210 if (file_version >= 81)
1212 gmx_fio_do_real(fio, ir->fourier_spacing);
1216 ir->fourier_spacing = 0.0;
1218 gmx_fio_do_int(fio, ir->nkx);
1219 gmx_fio_do_int(fio, ir->nky);
1220 gmx_fio_do_int(fio, ir->nkz);
1221 gmx_fio_do_int(fio, ir->pme_order);
1222 gmx_fio_do_real(fio, ir->ewald_rtol);
1224 if (file_version >= 93)
1226 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1230 ir->ewald_rtol_lj = ir->ewald_rtol;
1232 gmx_fio_do_int(fio, ir->ewald_geometry);
1233 gmx_fio_do_real(fio, ir->epsilon_surface);
1235 /* ignore bOptFFT */
1236 if (file_version < tpxv_RemoveObsoleteParameters1)
1238 gmx_fio_do_gmx_bool(fio, bdum);
1241 if (file_version >= 93)
1243 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1245 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1246 gmx_fio_do_int(fio, ir->etc);
1247 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1248 * but the values 0 and 1 still mean no and
1249 * berendsen temperature coupling, respectively.
1251 if (file_version >= 79)
1253 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1255 if (file_version >= 71)
1257 gmx_fio_do_int(fio, ir->nsttcouple);
1261 ir->nsttcouple = ir->nstcalcenergy;
1263 gmx_fio_do_int(fio, ir->epc);
1264 gmx_fio_do_int(fio, ir->epct);
1265 if (file_version >= 71)
1267 gmx_fio_do_int(fio, ir->nstpcouple);
1271 ir->nstpcouple = ir->nstcalcenergy;
1273 gmx_fio_do_real(fio, ir->tau_p);
1274 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1275 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1276 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1277 gmx_fio_do_rvec(fio, ir->compress[XX]);
1278 gmx_fio_do_rvec(fio, ir->compress[YY]);
1279 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1280 gmx_fio_do_int(fio, ir->refcoord_scaling);
1281 gmx_fio_do_rvec(fio, ir->posres_com);
1282 gmx_fio_do_rvec(fio, ir->posres_comB);
1284 if (file_version < 79)
1286 gmx_fio_do_int(fio, ir->andersen_seed);
1290 ir->andersen_seed = 0;
1293 gmx_fio_do_real(fio, ir->shake_tol);
1295 gmx_fio_do_int(fio, ir->efep);
1296 do_fepvals(fio, ir->fepvals, bRead, file_version);
1298 if (file_version >= 79)
1300 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1303 ir->bSimTemp = TRUE;
1308 ir->bSimTemp = FALSE;
1312 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1315 if (file_version >= 79)
1317 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1320 ir->bExpanded = TRUE;
1324 ir->bExpanded = FALSE;
1329 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1332 gmx_fio_do_int(fio, ir->eDisre);
1333 gmx_fio_do_int(fio, ir->eDisreWeighting);
1334 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1335 gmx_fio_do_real(fio, ir->dr_fc);
1336 gmx_fio_do_real(fio, ir->dr_tau);
1337 gmx_fio_do_int(fio, ir->nstdisreout);
1338 gmx_fio_do_real(fio, ir->orires_fc);
1339 gmx_fio_do_real(fio, ir->orires_tau);
1340 gmx_fio_do_int(fio, ir->nstorireout);
1342 /* ignore dihre_fc */
1343 if (file_version < 79)
1345 gmx_fio_do_real(fio, rdum);
1348 gmx_fio_do_real(fio, ir->em_stepsize);
1349 gmx_fio_do_real(fio, ir->em_tol);
1350 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1351 gmx_fio_do_int(fio, ir->niter);
1352 gmx_fio_do_real(fio, ir->fc_stepsize);
1353 gmx_fio_do_int(fio, ir->eConstrAlg);
1354 gmx_fio_do_int(fio, ir->nProjOrder);
1355 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1356 gmx_fio_do_int(fio, ir->nLincsIter);
1357 gmx_fio_do_real(fio, ir->bd_fric);
1358 if (file_version >= tpxv_Use64BitRandomSeed)
1360 gmx_fio_do_int64(fio, ir->ld_seed);
1364 gmx_fio_do_int(fio, idum);
1368 for (i = 0; i < DIM; i++)
1370 gmx_fio_do_rvec(fio, ir->deform[i]);
1372 gmx_fio_do_real(fio, ir->cos_accel);
1374 gmx_fio_do_int(fio, ir->userint1);
1375 gmx_fio_do_int(fio, ir->userint2);
1376 gmx_fio_do_int(fio, ir->userint3);
1377 gmx_fio_do_int(fio, ir->userint4);
1378 gmx_fio_do_real(fio, ir->userreal1);
1379 gmx_fio_do_real(fio, ir->userreal2);
1380 gmx_fio_do_real(fio, ir->userreal3);
1381 gmx_fio_do_real(fio, ir->userreal4);
1383 /* AdResS is removed, but we need to be able to read old files,
1384 and let mdrun refuse to run them */
1385 if (file_version >= 77 && file_version < tpxv_RemoveAdress)
1387 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1390 int idum, numThermoForceGroups, numEnergyGroups;
1393 gmx_fio_do_int(fio, idum);
1394 gmx_fio_do_real(fio, rdum);
1395 gmx_fio_do_real(fio, rdum);
1396 gmx_fio_do_real(fio, rdum);
1397 gmx_fio_do_int(fio, idum);
1398 gmx_fio_do_int(fio, idum);
1399 gmx_fio_do_rvec(fio, rvecdum);
1400 gmx_fio_do_int(fio, numThermoForceGroups);
1401 gmx_fio_do_real(fio, rdum);
1402 gmx_fio_do_int(fio, numEnergyGroups);
1403 gmx_fio_do_int(fio, idum);
1405 if (numThermoForceGroups > 0)
1407 std::vector<int> idumn(numThermoForceGroups);
1408 gmx_fio_ndo_int(fio, idumn.data(), idumn.size());
1410 if (numEnergyGroups > 0)
1412 std::vector<int> idumn(numEnergyGroups);
1413 gmx_fio_ndo_int(fio, idumn.data(), idumn.size());
1419 ir->bAdress = FALSE;
1426 if (file_version >= tpxv_PullCoordTypeGeom)
1428 gmx_fio_do_gmx_bool(fio, ir->bPull);
1432 gmx_fio_do_int(fio, ePullOld);
1433 ir->bPull = (ePullOld > 0);
1434 /* We removed the first ePull=ePullNo for the enum */
1443 do_pull(fio, ir->pull, bRead, file_version, ePullOld);
1447 if (file_version >= tpxv_AcceleratedWeightHistogram)
1449 gmx_fio_do_gmx_bool(fio, ir->bDoAwh);
1455 snew(ir->awhParams, 1);
1457 do_awh(fio, ir->awhParams, bRead, file_version);
1465 /* Enforced rotation */
1466 if (file_version >= 74)
1468 gmx_fio_do_gmx_bool(fio, ir->bRot);
1475 do_rot(fio, ir->rot, bRead);
1483 /* Interactive molecular dynamics */
1484 if (file_version >= tpxv_InteractiveMolecularDynamics)
1486 gmx_fio_do_gmx_bool(fio, ir->bIMD);
1493 do_imd(fio, ir->imd, bRead);
1498 /* We don't support IMD sessions for old .tpr files */
1503 gmx_fio_do_int(fio, ir->opts.ngtc);
1504 if (file_version >= 69)
1506 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1510 ir->opts.nhchainlength = 1;
1512 gmx_fio_do_int(fio, ir->opts.ngacc);
1513 gmx_fio_do_int(fio, ir->opts.ngfrz);
1514 gmx_fio_do_int(fio, ir->opts.ngener);
1518 snew(ir->opts.nrdf, ir->opts.ngtc);
1519 snew(ir->opts.ref_t, ir->opts.ngtc);
1520 snew(ir->opts.annealing, ir->opts.ngtc);
1521 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1522 snew(ir->opts.anneal_time, ir->opts.ngtc);
1523 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1524 snew(ir->opts.tau_t, ir->opts.ngtc);
1525 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1526 snew(ir->opts.acc, ir->opts.ngacc);
1527 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1529 if (ir->opts.ngtc > 0)
1531 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1532 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1533 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1535 if (ir->opts.ngfrz > 0)
1537 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1539 if (ir->opts.ngacc > 0)
1541 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1543 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1544 ir->opts.ngener*ir->opts.ngener);
1546 /* First read the lists with annealing and npoints for each group */
1547 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1548 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1549 for (j = 0; j < (ir->opts.ngtc); j++)
1551 k = ir->opts.anneal_npoints[j];
1554 snew(ir->opts.anneal_time[j], k);
1555 snew(ir->opts.anneal_temp[j], k);
1557 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1558 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1562 gmx_fio_do_int(fio, ir->nwall);
1563 gmx_fio_do_int(fio, ir->wall_type);
1564 gmx_fio_do_real(fio, ir->wall_r_linpot);
1565 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1566 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1567 gmx_fio_do_real(fio, ir->wall_density[0]);
1568 gmx_fio_do_real(fio, ir->wall_density[1]);
1569 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1572 /* Cosine stuff for electric fields */
1573 if (file_version < tpxv_GenericParamsForElectricField)
1575 do_legacy_efield(fio, ¶msObj);
1579 if (file_version >= tpxv_ComputationalElectrophysiology)
1581 gmx_fio_do_int(fio, ir->eSwapCoords);
1582 if (ir->eSwapCoords != eswapNO)
1588 do_swapcoords_tpx(fio, ir->swap, bRead, file_version);
1594 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1595 gmx_fio_do_int(fio, ir->QMMMscheme);
1596 gmx_fio_do_real(fio, ir->scalefactor);
1597 gmx_fio_do_int(fio, ir->opts.ngQM);
1600 snew(ir->opts.QMmethod, ir->opts.ngQM);
1601 snew(ir->opts.QMbasis, ir->opts.ngQM);
1602 snew(ir->opts.QMcharge, ir->opts.ngQM);
1603 snew(ir->opts.QMmult, ir->opts.ngQM);
1604 snew(ir->opts.bSH, ir->opts.ngQM);
1605 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1606 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1607 snew(ir->opts.SAon, ir->opts.ngQM);
1608 snew(ir->opts.SAoff, ir->opts.ngQM);
1609 snew(ir->opts.SAsteps, ir->opts.ngQM);
1611 if (ir->opts.ngQM > 0)
1613 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1614 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1615 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1616 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1617 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1618 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1619 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1620 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1621 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1622 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1623 /* We leave in dummy i/o for removed parameters to avoid
1624 * changing the tpr format for every QMMM change.
1626 std::vector<int> dummy(ir->opts.ngQM, 0);
1627 gmx_fio_ndo_int(fio, dummy.data(), ir->opts.ngQM);
1628 gmx_fio_ndo_int(fio, dummy.data(), ir->opts.ngQM);
1630 /* end of QMMM stuff */
1633 if (file_version >= tpxv_GenericParamsForElectricField)
1635 gmx::FileIOXdrSerializer serializer(fio);
1638 paramsObj.mergeObject(
1639 gmx::deserializeKeyValueTree(&serializer));
1643 GMX_RELEASE_ASSERT(ir->params != nullptr,
1644 "Parameters should be present when writing inputrec");
1645 gmx::serializeKeyValueTree(*ir->params, &serializer);
1650 ir->params = new gmx::KeyValueTreeObject(paramsBuilder.build());
1655 static void do_harm(t_fileio *fio, t_iparams *iparams)
1657 gmx_fio_do_real(fio, iparams->harmonic.rA);
1658 gmx_fio_do_real(fio, iparams->harmonic.krA);
1659 gmx_fio_do_real(fio, iparams->harmonic.rB);
1660 gmx_fio_do_real(fio, iparams->harmonic.krB);
1663 static void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1664 gmx_bool bRead, int file_version)
1677 do_harm(fio, iparams);
1678 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1680 /* Correct incorrect storage of parameters */
1681 iparams->pdihs.phiB = iparams->pdihs.phiA;
1682 iparams->pdihs.cpB = iparams->pdihs.cpA;
1686 gmx_fio_do_real(fio, iparams->harmonic.rA);
1687 gmx_fio_do_real(fio, iparams->harmonic.krA);
1689 case F_LINEAR_ANGLES:
1690 gmx_fio_do_real(fio, iparams->linangle.klinA);
1691 gmx_fio_do_real(fio, iparams->linangle.aA);
1692 gmx_fio_do_real(fio, iparams->linangle.klinB);
1693 gmx_fio_do_real(fio, iparams->linangle.aB);
1696 gmx_fio_do_real(fio, iparams->fene.bm);
1697 gmx_fio_do_real(fio, iparams->fene.kb);
1701 gmx_fio_do_real(fio, iparams->restraint.lowA);
1702 gmx_fio_do_real(fio, iparams->restraint.up1A);
1703 gmx_fio_do_real(fio, iparams->restraint.up2A);
1704 gmx_fio_do_real(fio, iparams->restraint.kA);
1705 gmx_fio_do_real(fio, iparams->restraint.lowB);
1706 gmx_fio_do_real(fio, iparams->restraint.up1B);
1707 gmx_fio_do_real(fio, iparams->restraint.up2B);
1708 gmx_fio_do_real(fio, iparams->restraint.kB);
1714 gmx_fio_do_real(fio, iparams->tab.kA);
1715 gmx_fio_do_int(fio, iparams->tab.table);
1716 gmx_fio_do_real(fio, iparams->tab.kB);
1718 case F_CROSS_BOND_BONDS:
1719 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1720 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1721 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1723 case F_CROSS_BOND_ANGLES:
1724 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1725 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1726 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1727 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1729 case F_UREY_BRADLEY:
1730 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1731 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1732 gmx_fio_do_real(fio, iparams->u_b.r13A);
1733 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1734 if (file_version >= 79)
1736 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1737 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1738 gmx_fio_do_real(fio, iparams->u_b.r13B);
1739 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1743 iparams->u_b.thetaB = iparams->u_b.thetaA;
1744 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1745 iparams->u_b.r13B = iparams->u_b.r13A;
1746 iparams->u_b.kUBB = iparams->u_b.kUBA;
1749 case F_QUARTIC_ANGLES:
1750 gmx_fio_do_real(fio, iparams->qangle.theta);
1751 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1754 gmx_fio_do_real(fio, iparams->bham.a);
1755 gmx_fio_do_real(fio, iparams->bham.b);
1756 gmx_fio_do_real(fio, iparams->bham.c);
1759 gmx_fio_do_real(fio, iparams->morse.b0A);
1760 gmx_fio_do_real(fio, iparams->morse.cbA);
1761 gmx_fio_do_real(fio, iparams->morse.betaA);
1762 if (file_version >= 79)
1764 gmx_fio_do_real(fio, iparams->morse.b0B);
1765 gmx_fio_do_real(fio, iparams->morse.cbB);
1766 gmx_fio_do_real(fio, iparams->morse.betaB);
1770 iparams->morse.b0B = iparams->morse.b0A;
1771 iparams->morse.cbB = iparams->morse.cbA;
1772 iparams->morse.betaB = iparams->morse.betaA;
1776 gmx_fio_do_real(fio, iparams->cubic.b0);
1777 gmx_fio_do_real(fio, iparams->cubic.kb);
1778 gmx_fio_do_real(fio, iparams->cubic.kcub);
1782 case F_POLARIZATION:
1783 gmx_fio_do_real(fio, iparams->polarize.alpha);
1786 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1787 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1788 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1791 gmx_fio_do_real(fio, iparams->wpol.al_x);
1792 gmx_fio_do_real(fio, iparams->wpol.al_y);
1793 gmx_fio_do_real(fio, iparams->wpol.al_z);
1794 gmx_fio_do_real(fio, iparams->wpol.rOH);
1795 gmx_fio_do_real(fio, iparams->wpol.rHH);
1796 gmx_fio_do_real(fio, iparams->wpol.rOD);
1799 gmx_fio_do_real(fio, iparams->thole.a);
1800 gmx_fio_do_real(fio, iparams->thole.alpha1);
1801 gmx_fio_do_real(fio, iparams->thole.alpha2);
1802 gmx_fio_do_real(fio, iparams->thole.rfac);
1805 gmx_fio_do_real(fio, iparams->lj.c6);
1806 gmx_fio_do_real(fio, iparams->lj.c12);
1809 gmx_fio_do_real(fio, iparams->lj14.c6A);
1810 gmx_fio_do_real(fio, iparams->lj14.c12A);
1811 gmx_fio_do_real(fio, iparams->lj14.c6B);
1812 gmx_fio_do_real(fio, iparams->lj14.c12B);
1815 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1816 gmx_fio_do_real(fio, iparams->ljc14.qi);
1817 gmx_fio_do_real(fio, iparams->ljc14.qj);
1818 gmx_fio_do_real(fio, iparams->ljc14.c6);
1819 gmx_fio_do_real(fio, iparams->ljc14.c12);
1821 case F_LJC_PAIRS_NB:
1822 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1823 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1824 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1825 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1831 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1832 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1833 gmx_fio_do_real(fio, iparams->pdihs.phiB);
1834 gmx_fio_do_real(fio, iparams->pdihs.cpB);
1835 gmx_fio_do_int(fio, iparams->pdihs.mult);
1838 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1839 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1842 gmx_fio_do_int(fio, iparams->disres.label);
1843 gmx_fio_do_int(fio, iparams->disres.type);
1844 gmx_fio_do_real(fio, iparams->disres.low);
1845 gmx_fio_do_real(fio, iparams->disres.up1);
1846 gmx_fio_do_real(fio, iparams->disres.up2);
1847 gmx_fio_do_real(fio, iparams->disres.kfac);
1850 gmx_fio_do_int(fio, iparams->orires.ex);
1851 gmx_fio_do_int(fio, iparams->orires.label);
1852 gmx_fio_do_int(fio, iparams->orires.power);
1853 gmx_fio_do_real(fio, iparams->orires.c);
1854 gmx_fio_do_real(fio, iparams->orires.obs);
1855 gmx_fio_do_real(fio, iparams->orires.kfac);
1858 if (file_version < 82)
1860 gmx_fio_do_int(fio, idum);
1861 gmx_fio_do_int(fio, idum);
1863 gmx_fio_do_real(fio, iparams->dihres.phiA);
1864 gmx_fio_do_real(fio, iparams->dihres.dphiA);
1865 gmx_fio_do_real(fio, iparams->dihres.kfacA);
1866 if (file_version >= 82)
1868 gmx_fio_do_real(fio, iparams->dihres.phiB);
1869 gmx_fio_do_real(fio, iparams->dihres.dphiB);
1870 gmx_fio_do_real(fio, iparams->dihres.kfacB);
1874 iparams->dihres.phiB = iparams->dihres.phiA;
1875 iparams->dihres.dphiB = iparams->dihres.dphiA;
1876 iparams->dihres.kfacB = iparams->dihres.kfacA;
1880 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
1881 gmx_fio_do_rvec(fio, iparams->posres.fcA);
1882 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
1883 gmx_fio_do_rvec(fio, iparams->posres.fcB);
1886 gmx_fio_do_int(fio, iparams->fbposres.geom);
1887 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
1888 gmx_fio_do_real(fio, iparams->fbposres.r);
1889 gmx_fio_do_real(fio, iparams->fbposres.k);
1892 gmx_fio_ndo_real(fio, iparams->cbtdihs.cbtcA, NR_CBTDIHS);
1895 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1896 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1899 /* Fourier dihedrals are internally represented
1900 * as Ryckaert-Bellemans since those are faster to compute.
1902 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1903 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1907 gmx_fio_do_real(fio, iparams->constr.dA);
1908 gmx_fio_do_real(fio, iparams->constr.dB);
1911 gmx_fio_do_real(fio, iparams->settle.doh);
1912 gmx_fio_do_real(fio, iparams->settle.dhh);
1915 gmx_fio_do_real(fio, iparams->vsite.a);
1920 gmx_fio_do_real(fio, iparams->vsite.a);
1921 gmx_fio_do_real(fio, iparams->vsite.b);
1926 gmx_fio_do_real(fio, iparams->vsite.a);
1927 gmx_fio_do_real(fio, iparams->vsite.b);
1928 gmx_fio_do_real(fio, iparams->vsite.c);
1931 gmx_fio_do_int(fio, iparams->vsiten.n);
1932 gmx_fio_do_real(fio, iparams->vsiten.a);
1934 case F_GB12_NOLONGERUSED:
1935 case F_GB13_NOLONGERUSED:
1936 case F_GB14_NOLONGERUSED:
1937 // Implicit solvent parameters can still be read, but never used
1940 if (file_version < 68)
1942 gmx_fio_do_real(fio, rdum);
1943 gmx_fio_do_real(fio, rdum);
1944 gmx_fio_do_real(fio, rdum);
1945 gmx_fio_do_real(fio, rdum);
1947 if (file_version < tpxv_RemoveImplicitSolvation)
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);
1953 gmx_fio_do_real(fio, rdum);
1958 gmx_fio_do_int(fio, iparams->cmap.cmapA);
1959 gmx_fio_do_int(fio, iparams->cmap.cmapB);
1962 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1963 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1967 static void do_ilist(t_fileio *fio, InteractionList *ilist, gmx_bool bRead)
1969 int nr = ilist->size();
1970 gmx_fio_do_int(fio, nr);
1973 ilist->iatoms.resize(nr);
1975 gmx_fio_ndo_int(fio, ilist->iatoms.data(), ilist->size());
1978 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1979 gmx_bool bRead, int file_version)
1981 gmx_fio_do_int(fio, ffparams->atnr);
1982 gmx_fio_do_int(fio, ffparams->ntypes);
1985 snew(ffparams->functype, ffparams->ntypes);
1986 snew(ffparams->iparams, ffparams->ntypes);
1988 /* Read/write all the function types */
1989 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
1991 if (file_version >= 66)
1993 gmx_fio_do_double(fio, ffparams->reppow);
1997 ffparams->reppow = 12.0;
2000 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2002 /* Check whether all these function types are supported by the code.
2003 * In practice the code is backwards compatible, which means that the
2004 * numbering may have to be altered from old numbering to new numbering
2006 for (int i = 0; i < ffparams->ntypes; i++)
2010 /* Loop over file versions */
2011 for (int k = 0; k < NFTUPD; k++)
2013 /* Compare the read file_version to the update table */
2014 if ((file_version < ftupd[k].fvnr) &&
2015 (ffparams->functype[i] >= ftupd[k].ftype))
2017 ffparams->functype[i] += 1;
2022 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2027 static void add_settle_atoms(InteractionList *ilist)
2031 /* Settle used to only store the first atom: add the other two */
2032 ilist->iatoms.resize(2*ilist->size());
2033 for (i = ilist->size()/4 - 1; i >= 0; i--)
2035 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2036 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2037 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2038 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2042 static void do_ilists(t_fileio *fio, InteractionLists *ilists, gmx_bool bRead,
2045 GMX_RELEASE_ASSERT(ilists, "Need a valid ilists object");
2046 GMX_RELEASE_ASSERT(ilists->size() == F_NRE, "The code needs to be in sync with InteractionLists");
2048 for (int j = 0; j < F_NRE; j++)
2050 InteractionList &ilist = (*ilists)[j];
2051 gmx_bool bClear = FALSE;
2054 for (int k = 0; k < NFTUPD; k++)
2056 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2064 ilist.iatoms.clear();
2068 do_ilist(fio, &ilist, bRead);
2069 if (file_version < 78 && j == F_SETTLE && ilist.size() > 0)
2071 add_settle_atoms(&ilist);
2077 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead)
2079 gmx_fio_do_int(fio, block->nr);
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);
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).
2112 atomicnumber_to_element(int atomicnumber)
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;
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);
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);
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)
2189 gmx_fio_do_int(fio, ls);
2190 *nm = get_symtab_handle(symtab, ls);
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,
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)
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);
2231 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2236 gmx_fio_do_int(fio, atoms->nr);
2237 gmx_fio_do_int(fio, atoms->nres);
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;
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)
2277 do_grps(fio, egcNR, groups->grps, bRead);
2278 gmx_fio_do_int(fio, groups->ngrpname);
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)
2291 groups->grpnr[g] = nullptr;
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,
2310 gmx_fio_do_int(fio, atomtypes->nr);
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)
2339 gmx_fio_do_int(fio, symtab->nr);
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);
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]);
2363 symbuf = symbuf->next;
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;
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,
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)
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)
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)
2459 gmx_mtop_ilistloop_t iloop;
2462 ip = mtop->ffparams.iparams;
2464 iloop = gmx_mtop_ilistloop_init(mtop);
2465 while (const InteractionLists *ilist = gmx_mtop_ilistloop_next(iloop, &nmol))
2467 const InteractionList &il = (*ilist)[F_DISRES];
2471 gmx::ArrayRef<const int> a = il.iatoms;
2473 for (int i = 0; i < il.size(); i += 3)
2476 if (i+3 == il.size() || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2478 ip[a[i]].disres.npair = npair;
2486 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2493 do_symtab(fio, &(mtop->symtab), bRead);
2495 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2497 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
2499 int nmoltype = mtop->moltype.size();
2500 gmx_fio_do_int(fio, nmoltype);
2503 mtop->moltype.resize(nmoltype);
2505 for (gmx_moltype_t &moltype : mtop->moltype)
2507 do_moltype(fio, &moltype, bRead, &mtop->symtab, file_version);
2510 int nmolblock = mtop->molblock.size();
2511 gmx_fio_do_int(fio, nmolblock);
2514 mtop->molblock.resize(nmolblock);
2516 for (gmx_molblock_t &molblock : mtop->molblock)
2518 int numAtomsPerMolecule = (bRead ? 0 : mtop->moltype[molblock.type].atoms.nr);
2519 do_molblock(fio, &molblock, numAtomsPerMolecule, bRead);
2521 gmx_fio_do_int(fio, mtop->natoms);
2523 if (file_version >= tpxv_IntermolecularBondeds)
2525 gmx_fio_do_gmx_bool(fio, mtop->bIntermolecularInteractions);
2526 if (mtop->bIntermolecularInteractions)
2530 mtop->intermolecular_ilist = gmx::compat::make_unique<InteractionLists>();
2532 do_ilists(fio, mtop->intermolecular_ilist.get(), bRead, file_version);
2537 mtop->bIntermolecularInteractions = FALSE;
2540 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
2542 if (file_version >= 65)
2544 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
2548 mtop->ffparams.cmap_grid.ngrid = 0;
2549 mtop->ffparams.cmap_grid.grid_spacing = 0;
2550 mtop->ffparams.cmap_grid.cmapdata = nullptr;
2553 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab));
2555 mtop->haveMoleculeIndices = true;
2559 close_symtab(&(mtop->symtab));
2563 /* If TopOnlyOK is TRUE then we can read even future versions
2564 * of tpx files, provided the fileGeneration hasn't changed.
2565 * If it is FALSE, we need the inputrecord too, and bail out
2566 * if the file is newer than the program.
2568 * The version and generation of the topology (see top of this file)
2569 * are returned in the two last arguments, if those arguments are non-NULL.
2571 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2573 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
2574 gmx_bool TopOnlyOK, int *fileVersionPointer, int *fileGenerationPointer)
2577 char file_tag[STRLEN];
2582 int fileVersion; /* Version number of the code that wrote the file */
2583 int fileGeneration; /* Generation version number of the code that wrote the file */
2585 /* XDR binary topology file */
2586 precision = sizeof(real);
2589 gmx_fio_do_string(fio, buf);
2590 if (std::strncmp(buf, "VERSION", 7) != 0)
2592 gmx_fatal(FARGS, "Can not read file %s,\n"
2593 " this file is from a GROMACS version which is older than 2.0\n"
2594 " Make a new one with grompp or use a gro or pdb file, if possible",
2595 gmx_fio_getname(fio));
2597 gmx_fio_do_int(fio, precision);
2598 bDouble = (precision == sizeof(double));
2599 if ((precision != sizeof(float)) && !bDouble)
2601 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
2602 "instead of %zu or %zu",
2603 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
2605 gmx_fio_setprecision(fio, bDouble);
2606 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
2607 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
2611 snprintf(buf, STRLEN, "VERSION %s", gmx_version());
2612 gmx_fio_write_string(fio, buf);
2613 bDouble = (precision == sizeof(double));
2614 gmx_fio_setprecision(fio, bDouble);
2615 gmx_fio_do_int(fio, precision);
2616 fileVersion = tpx_version;
2617 sprintf(file_tag, "%s", tpx_tag);
2618 fileGeneration = tpx_generation;
2621 /* Check versions! */
2622 gmx_fio_do_int(fio, fileVersion);
2624 /* This is for backward compatibility with development versions 77-79
2625 * where the tag was, mistakenly, placed before the generation,
2626 * which would cause a segv instead of a proper error message
2627 * when reading the topology only from tpx with <77 code.
2629 if (fileVersion >= 77 && fileVersion <= 79)
2631 gmx_fio_do_string(fio, file_tag);
2634 gmx_fio_do_int(fio, fileGeneration);
2636 if (fileVersion >= 81)
2638 gmx_fio_do_string(fio, file_tag);
2642 if (fileVersion < 77)
2644 /* Versions before 77 don't have the tag, set it to release */
2645 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
2648 if (std::strcmp(file_tag, tpx_tag) != 0)
2650 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
2653 /* We only support reading tpx files with the same tag as the code
2654 * or tpx files with the release tag and with lower version number.
2656 if (std::strcmp(file_tag, TPX_TAG_RELEASE) != 0 && fileVersion < tpx_version)
2658 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2659 gmx_fio_getname(fio), fileVersion, file_tag,
2660 tpx_version, tpx_tag);
2665 if (fileVersionPointer)
2667 *fileVersionPointer = fileVersion;
2669 if (fileGenerationPointer)
2671 *fileGenerationPointer = fileGeneration;
2674 if ((fileVersion <= tpx_incompatible_version) ||
2675 ((fileVersion > tpx_version) && !TopOnlyOK) ||
2676 (fileGeneration > tpx_generation) ||
2677 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
2679 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
2680 gmx_fio_getname(fio), fileVersion, tpx_version);
2683 gmx_fio_do_int(fio, tpx->natoms);
2684 gmx_fio_do_int(fio, tpx->ngtc);
2686 if (fileVersion < 62)
2688 gmx_fio_do_int(fio, idum);
2689 gmx_fio_do_real(fio, rdum);
2691 if (fileVersion >= 79)
2693 gmx_fio_do_int(fio, tpx->fep_state);
2695 gmx_fio_do_real(fio, tpx->lambda);
2696 gmx_fio_do_gmx_bool(fio, tpx->bIr);
2697 gmx_fio_do_gmx_bool(fio, tpx->bTop);
2698 gmx_fio_do_gmx_bool(fio, tpx->bX);
2699 gmx_fio_do_gmx_bool(fio, tpx->bV);
2700 gmx_fio_do_gmx_bool(fio, tpx->bF);
2701 gmx_fio_do_gmx_bool(fio, tpx->bBox);
2703 if ((fileGeneration > tpx_generation))
2705 /* This can only happen if TopOnlyOK=TRUE */
2710 static int do_tpx(t_fileio *fio, gmx_bool bRead,
2711 t_inputrec *ir, t_state *state, rvec *x, rvec *v,
2717 gmx_bool bPeriodicMols;
2721 GMX_RELEASE_ASSERT(state != nullptr, "Cannot write a null state");
2722 GMX_RELEASE_ASSERT(x == nullptr && v == nullptr, "Passing separate x and v pointers to do_tpx() is not supported when writing");
2724 tpx.natoms = state->natoms;
2725 tpx.ngtc = state->ngtc;
2726 tpx.fep_state = state->fep_state;
2727 tpx.lambda = state->lambda[efptFEP];
2728 tpx.bIr = ir != nullptr;
2729 tpx.bTop = mtop != nullptr;
2730 tpx.bX = (state->flags & (1 << estX)) != 0;
2731 tpx.bV = (state->flags & (1 << estV)) != 0;
2737 GMX_RELEASE_ASSERT(!(x == nullptr && v != nullptr), "Passing x==NULL and v!=NULL is not supported");
2740 TopOnlyOK = (ir == nullptr);
2742 int fileVersion; /* Version number of the code that wrote the file */
2743 int fileGeneration; /* Generation version number of the code that wrote the file */
2744 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &fileVersion, &fileGeneration);
2749 init_gtc_state(state, tpx.ngtc, 0, 0);
2752 // v is also nullptr by the above assertion, so we may
2753 // need to make memory in state for storing the contents
2757 state->flags |= (1 << estX);
2761 state->flags |= (1 << estV);
2763 state_change_natoms(state, tpx.natoms);
2769 x = as_rvec_array(state->x.data());
2770 v = as_rvec_array(state->v.data());
2773 #define do_test(fio, b, p) if (bRead && ((p) != NULL) && !(b)) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
2775 do_test(fio, tpx.bBox, state->box);
2778 gmx_fio_ndo_rvec(fio, state->box, DIM);
2779 if (fileVersion >= 51)
2781 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
2785 /* We initialize box_rel after reading the inputrec */
2786 clear_mat(state->box_rel);
2788 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
2789 if (fileVersion < 56)
2792 gmx_fio_ndo_rvec(fio, mdum, DIM);
2796 if (state->ngtc > 0)
2799 snew(dumv, state->ngtc);
2800 if (fileVersion < 69)
2802 gmx_fio_ndo_real(fio, dumv, state->ngtc);
2804 /* These used to be the Berendsen tcoupl_lambda's */
2805 gmx_fio_ndo_real(fio, dumv, state->ngtc);
2809 do_test(fio, tpx.bTop, mtop);
2814 do_mtop(fio, mtop, bRead, fileVersion);
2819 do_mtop(fio, &dum_top, bRead, fileVersion);
2822 do_test(fio, tpx.bX, x);
2827 state->flags |= (1<<estX);
2829 gmx_fio_ndo_rvec(fio, x, tpx.natoms);
2832 do_test(fio, tpx.bV, v);
2837 state->flags |= (1<<estV);
2839 gmx_fio_ndo_rvec(fio, v, tpx.natoms);
2842 // No need to run do_test when the last argument is NULL
2846 snew(dummyForces, state->natoms);
2847 gmx_fio_ndo_rvec(fio, dummyForces, tpx.natoms);
2851 /* Starting with tpx version 26, we have the inputrec
2852 * at the end of the file, so we can ignore it
2853 * if the file is never than the software (but still the
2854 * same generation - see comments at the top of this file.
2859 bPeriodicMols = FALSE;
2861 do_test(fio, tpx.bIr, ir);
2864 if (fileVersion >= 53)
2866 /* Removed the pbc info from do_inputrec, since we always want it */
2870 bPeriodicMols = ir->bPeriodicMols;
2872 gmx_fio_do_int(fio, ePBC);
2873 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
2875 if (fileGeneration <= tpx_generation && ir)
2877 do_inputrec(fio, ir, bRead, fileVersion);
2878 if (fileVersion < 51)
2880 set_box_rel(ir, state);
2882 if (fileVersion < 53)
2885 bPeriodicMols = ir->bPeriodicMols;
2888 if (bRead && ir && fileVersion >= 53)
2890 /* We need to do this after do_inputrec, since that initializes ir */
2892 ir->bPeriodicMols = bPeriodicMols;
2900 if (state->ngtc == 0)
2902 /* Reading old version without tcoupl state data: set it */
2903 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
2905 if (tpx.bTop && mtop)
2907 if (fileVersion < 57)
2909 if (mtop->moltype[0].ilist[F_DISRES].size() > 0)
2911 ir->eDisre = edrSimple;
2915 ir->eDisre = edrNone;
2918 set_disres_npair(mtop);
2922 if (tpx.bTop && mtop)
2924 gmx_mtop_finalize(mtop);
2931 static t_fileio *open_tpx(const char *fn, const char *mode)
2933 return gmx_fio_open(fn, mode);
2936 static void close_tpx(t_fileio *fio)
2941 /************************************************************
2943 * The following routines are the exported ones
2945 ************************************************************/
2947 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK)
2951 fio = open_tpx(fn, "r");
2952 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, nullptr, nullptr);
2956 void write_tpx_state(const char *fn,
2957 const t_inputrec *ir, const t_state *state, const gmx_mtop_t *mtop)
2961 fio = open_tpx(fn, "w");
2963 const_cast<t_inputrec *>(ir),
2964 const_cast<t_state *>(state), nullptr, nullptr,
2965 const_cast<gmx_mtop_t *>(mtop));
2969 void read_tpx_state(const char *fn,
2970 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
2974 fio = open_tpx(fn, "r");
2975 do_tpx(fio, TRUE, ir, state, nullptr, nullptr, mtop);
2979 int read_tpx(const char *fn,
2980 t_inputrec *ir, matrix box, int *natoms,
2981 rvec *x, rvec *v, gmx_mtop_t *mtop)
2987 fio = open_tpx(fn, "r");
2988 ePBC = do_tpx(fio, TRUE, ir, &state, x, v, mtop);
2990 *natoms = mtop->natoms;
2993 copy_mat(state.box, box);
2999 int read_tpx_top(const char *fn,
3000 t_inputrec *ir, matrix box, int *natoms,
3001 rvec *x, rvec *v, t_topology *top)
3006 ePBC = read_tpx(fn, ir, box, natoms, x, v, &mtop);
3008 *top = gmx_mtop_t_to_t_topology(&mtop, true);
3013 gmx_bool fn2bTPX(const char *file)
3015 return (efTPR == fn2ftp(file));
3018 void pr_tpxheader(FILE *fp, int indent, const char *title, const t_tpxheader *sh)
3020 if (available(fp, sh, indent, title))
3022 indent = pr_title(fp, indent, title);
3023 pr_indent(fp, indent);
3024 fprintf(fp, "bIr = %spresent\n", sh->bIr ? "" : "not ");
3025 pr_indent(fp, indent);
3026 fprintf(fp, "bBox = %spresent\n", sh->bBox ? "" : "not ");
3027 pr_indent(fp, indent);
3028 fprintf(fp, "bTop = %spresent\n", sh->bTop ? "" : "not ");
3029 pr_indent(fp, indent);
3030 fprintf(fp, "bX = %spresent\n", sh->bX ? "" : "not ");
3031 pr_indent(fp, indent);
3032 fprintf(fp, "bV = %spresent\n", sh->bV ? "" : "not ");
3033 pr_indent(fp, indent);
3034 fprintf(fp, "bF = %spresent\n", sh->bF ? "" : "not ");
3036 pr_indent(fp, indent);
3037 fprintf(fp, "natoms = %d\n", sh->natoms);
3038 pr_indent(fp, indent);
3039 fprintf(fp, "lambda = %e\n", sh->lambda);