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, 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/fileio/filetypes.h"
51 #include "gromacs/fileio/gmxfio.h"
52 #include "gromacs/fileio/gmxfio-xdr.h"
53 #include "gromacs/fileio/txtdump.h"
54 #include "gromacs/gmxlib/ifunc.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdtypes/inputrec.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/mdtypes/pull-params.h"
59 #include "gromacs/mdtypes/state.h"
60 #include "gromacs/pbcutil/boxutilities.h"
61 #include "gromacs/pbcutil/pbc.h"
62 #include "gromacs/topology/block.h"
63 #include "gromacs/topology/mtop_util.h"
64 #include "gromacs/topology/symtab.h"
65 #include "gromacs/topology/topology.h"
66 #include "gromacs/utility/arraysize.h"
67 #include "gromacs/utility/baseversion.h"
68 #include "gromacs/utility/cstringutil.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/futil.h"
71 #include "gromacs/utility/smalloc.h"
73 #define TPX_TAG_RELEASE "release"
75 /*! \brief Tag string for the file format written to run input files
76 * written by this version of the code.
78 * Change this if you want to change the run input format in a feature
79 * branch. This ensures that there will not be different run input
80 * formats around which cannot be distinguished, while not causing
81 * problems rebasing the feature branch onto upstream changes. When
82 * merging with mainstream GROMACS, set this tag string back to
83 * TPX_TAG_RELEASE, and instead add an element to tpxv.
85 static const char *tpx_tag
= TPX_TAG_RELEASE
;
87 /*! \brief Enum of values that describe the contents of a tpr file
88 * whose format matches a version number
90 * The enum helps the code be more self-documenting and ensure merges
91 * do not silently resolve when two patches make the same bump. When
92 * adding new functionality, add a new element just above tpxv_Count
93 * in this enumeration, and write code below that does the right thing
94 * according to the value of file_version.
97 tpxv_ComputationalElectrophysiology
= 96, /**< support for ion/water position swaps (computational electrophysiology) */
98 tpxv_Use64BitRandomSeed
, /**< change ld_seed from int to gmx_int64_t */
99 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials
, /**< potentials for supporting coarse-grained force fields */
100 tpxv_InteractiveMolecularDynamics
, /**< interactive molecular dynamics (IMD) */
101 tpxv_RemoveObsoleteParameters1
, /**< remove optimize_fft, dihre_fc, nstcheckpoint */
102 tpxv_PullCoordTypeGeom
, /**< add pull type and geometry per group and flat-bottom */
103 tpxv_PullGeomDirRel
, /**< add pull geometry direction-relative */
104 tpxv_IntermolecularBondeds
, /**< permit inter-molecular bonded interactions in the topology */
105 tpxv_CompElWithSwapLayerOffset
, /**< added parameters for improved CompEl setups */
106 tpxv_CompElPolyatomicIonsAndMultipleIonTypes
, /**< CompEl now can handle polyatomic ions and more than two types of ions */
107 tpxv_RemoveAdress
, /**< removed support for AdResS */
108 tpxv_PullCoordNGroup
, /**< add ngroup to pull coord */
109 tpxv_RemoveTwinRange
, /**< removed support for twin-range interactions */
110 tpxv_Count
/**< the total number of tpxv versions */
113 /*! \brief Version number of the file format written to run input
114 * files by this version of the code.
116 * The tpx_version increases whenever the file format in the main
117 * development branch changes, due to an extension of the tpxv enum above.
118 * Backward compatibility for reading old run input files is maintained
119 * by checking this version number against that of the file and then using
120 * the correct code path.
122 * When developing a feature branch that needs to change the run input
123 * file format, change tpx_tag instead. */
124 static const int tpx_version
= tpxv_Count
- 1;
127 /* This number should only be increased when you edit the TOPOLOGY section
128 * or the HEADER of the tpx format.
129 * This way we can maintain forward compatibility too for all analysis tools
130 * and/or external programs that only need to know the atom/residue names,
131 * charges, and bond connectivity.
133 * It first appeared in tpx version 26, when I also moved the inputrecord
134 * to the end of the tpx file, so we can just skip it if we only
137 * In particular, it must be increased when adding new elements to
138 * ftupd, so that old code can read new .tpr files.
140 static const int tpx_generation
= 26;
142 /* This number should be the most recent backwards incompatible version
143 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
145 static const int tpx_incompatible_version
= 9;
149 /* Struct used to maintain tpx compatibility when function types are added */
151 int fvnr
; /* file version number in which the function type first appeared */
152 int ftype
; /* function type */
156 * The entries should be ordered in:
157 * 1. ascending function type number
158 * 2. ascending file version number
160 static const t_ftupd ftupd
[] = {
161 { 20, F_CUBICBONDS
},
166 { 43, F_TABBONDSNC
},
167 { 70, F_RESTRBONDS
},
168 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials
, F_RESTRANGLES
},
169 { 76, F_LINEAR_ANGLES
},
170 { 30, F_CROSS_BOND_BONDS
},
171 { 30, F_CROSS_BOND_ANGLES
},
172 { 30, F_UREY_BRADLEY
},
173 { 34, F_QUARTIC_ANGLES
},
175 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials
, F_RESTRDIHS
},
176 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials
, F_CBTDIHS
},
185 { 72, F_NPSOLVATION
},
187 { 41, F_LJC_PAIRS_NB
},
188 { 32, F_BHAM_LR_NOLONGERUSED
},
190 { 32, F_COUL_RECIP
},
193 { 30, F_POLARIZATION
},
196 { 22, F_DISRESVIOL
},
200 { 26, F_DIHRESVIOL
},
205 { 46, F_ECONSERVED
},
206 { 69, F_VTEMP_NOLONGERUSED
},
208 { 54, F_DVDL_CONSTR
},
209 { 76, F_ANHARM_POL
},
212 { 79, F_DVDL_BONDED
, },
213 { 79, F_DVDL_RESTRAINT
},
214 { 79, F_DVDL_TEMPERATURE
},
216 #define NFTUPD asize(ftupd)
218 /* Needed for backward compatibility */
221 /**************************************************************
223 * Now the higer level routines that do io of the structures and arrays
225 **************************************************************/
226 static void do_pullgrp_tpx_pre95(t_fileio
*fio
,
234 gmx_fio_do_int(fio
, pgrp
->nat
);
237 snew(pgrp
->ind
, pgrp
->nat
);
239 gmx_fio_ndo_int(fio
, pgrp
->ind
, pgrp
->nat
);
240 gmx_fio_do_int(fio
, pgrp
->nweight
);
243 snew(pgrp
->weight
, pgrp
->nweight
);
245 gmx_fio_ndo_real(fio
, pgrp
->weight
, pgrp
->nweight
);
246 gmx_fio_do_int(fio
, pgrp
->pbcatom
);
247 gmx_fio_do_rvec(fio
, pcrd
->vec
);
248 clear_rvec(pcrd
->origin
);
249 gmx_fio_do_rvec(fio
, tmp
);
251 gmx_fio_do_real(fio
, pcrd
->rate
);
252 gmx_fio_do_real(fio
, pcrd
->k
);
253 if (file_version
>= 56)
255 gmx_fio_do_real(fio
, pcrd
->kB
);
263 static void do_pull_group(t_fileio
*fio
, t_pull_group
*pgrp
, gmx_bool bRead
)
265 gmx_fio_do_int(fio
, pgrp
->nat
);
268 snew(pgrp
->ind
, pgrp
->nat
);
270 gmx_fio_ndo_int(fio
, pgrp
->ind
, pgrp
->nat
);
271 gmx_fio_do_int(fio
, pgrp
->nweight
);
274 snew(pgrp
->weight
, pgrp
->nweight
);
276 gmx_fio_ndo_real(fio
, pgrp
->weight
, pgrp
->nweight
);
277 gmx_fio_do_int(fio
, pgrp
->pbcatom
);
280 static void do_pull_coord(t_fileio
*fio
, t_pull_coord
*pcrd
, int file_version
,
281 int ePullOld
, int eGeomOld
, ivec dimOld
)
283 if (file_version
>= tpxv_PullCoordNGroup
)
285 gmx_fio_do_int(fio
, pcrd
->eType
);
286 /* Note that we try to support adding new geometries without
287 * changing the tpx version. This requires checks when printing the
288 * geometry string and a check and fatal_error in init_pull.
290 gmx_fio_do_int(fio
, pcrd
->eGeom
);
291 gmx_fio_do_int(fio
, pcrd
->ngroup
);
292 if (pcrd
->ngroup
<= PULL_COORD_NGROUP_MAX
)
294 gmx_fio_ndo_int(fio
, pcrd
->group
, pcrd
->ngroup
);
298 /* More groups in file than supported, this must be a new geometry
299 * that is not supported by our current code. Since we will not
300 * use the groups for this coord (checks in the pull and WHAM code
301 * ensure this), we can ignore the groups and set ngroup=0.
304 snew(dum
, pcrd
->ngroup
);
305 gmx_fio_ndo_int(fio
, dum
, pcrd
->ngroup
);
310 gmx_fio_do_ivec(fio
, pcrd
->dim
);
315 gmx_fio_do_int(fio
, pcrd
->group
[0]);
316 gmx_fio_do_int(fio
, pcrd
->group
[1]);
317 if (file_version
>= tpxv_PullCoordTypeGeom
)
319 pcrd
->ngroup
= (pcrd
->eGeom
== epullgDIRRELATIVE
? 4 : 2);
320 gmx_fio_do_int(fio
, pcrd
->eType
);
321 gmx_fio_do_int(fio
, pcrd
->eGeom
);
322 if (pcrd
->ngroup
== 4)
324 gmx_fio_do_int(fio
, pcrd
->group
[2]);
325 gmx_fio_do_int(fio
, pcrd
->group
[3]);
327 gmx_fio_do_ivec(fio
, pcrd
->dim
);
331 pcrd
->eType
= ePullOld
;
332 pcrd
->eGeom
= eGeomOld
;
333 copy_ivec(dimOld
, pcrd
->dim
);
336 gmx_fio_do_rvec(fio
, pcrd
->origin
);
337 gmx_fio_do_rvec(fio
, pcrd
->vec
);
338 if (file_version
>= tpxv_PullCoordTypeGeom
)
340 gmx_fio_do_gmx_bool(fio
, pcrd
->bStart
);
344 /* This parameter is only printed, but not actually used by mdrun */
345 pcrd
->bStart
= FALSE
;
347 gmx_fio_do_real(fio
, pcrd
->init
);
348 gmx_fio_do_real(fio
, pcrd
->rate
);
349 gmx_fio_do_real(fio
, pcrd
->k
);
350 gmx_fio_do_real(fio
, pcrd
->kB
);
353 static void do_expandedvals(t_fileio
*fio
, t_expanded
*expand
, t_lambda
*fepvals
, gmx_bool bRead
, int file_version
)
355 int n_lambda
= fepvals
->n_lambda
;
357 /* reset the lambda calculation window */
358 fepvals
->lambda_start_n
= 0;
359 fepvals
->lambda_stop_n
= n_lambda
;
360 if (file_version
>= 79)
366 snew(expand
->init_lambda_weights
, n_lambda
);
368 gmx_fio_ndo_real(fio
, expand
->init_lambda_weights
, n_lambda
);
369 gmx_fio_do_gmx_bool(fio
, expand
->bInit_weights
);
372 gmx_fio_do_int(fio
, expand
->nstexpanded
);
373 gmx_fio_do_int(fio
, expand
->elmcmove
);
374 gmx_fio_do_int(fio
, expand
->elamstats
);
375 gmx_fio_do_int(fio
, expand
->lmc_repeats
);
376 gmx_fio_do_int(fio
, expand
->gibbsdeltalam
);
377 gmx_fio_do_int(fio
, expand
->lmc_forced_nstart
);
378 gmx_fio_do_int(fio
, expand
->lmc_seed
);
379 gmx_fio_do_real(fio
, expand
->mc_temp
);
380 gmx_fio_do_int(fio
, expand
->bSymmetrizedTMatrix
);
381 gmx_fio_do_int(fio
, expand
->nstTij
);
382 gmx_fio_do_int(fio
, expand
->minvarmin
);
383 gmx_fio_do_int(fio
, expand
->c_range
);
384 gmx_fio_do_real(fio
, expand
->wl_scale
);
385 gmx_fio_do_real(fio
, expand
->wl_ratio
);
386 gmx_fio_do_real(fio
, expand
->init_wl_delta
);
387 gmx_fio_do_gmx_bool(fio
, expand
->bWLoneovert
);
388 gmx_fio_do_int(fio
, expand
->elmceq
);
389 gmx_fio_do_int(fio
, expand
->equil_steps
);
390 gmx_fio_do_int(fio
, expand
->equil_samples
);
391 gmx_fio_do_int(fio
, expand
->equil_n_at_lam
);
392 gmx_fio_do_real(fio
, expand
->equil_wl_delta
);
393 gmx_fio_do_real(fio
, expand
->equil_ratio
);
397 static void do_simtempvals(t_fileio
*fio
, t_simtemp
*simtemp
, int n_lambda
, gmx_bool bRead
,
400 if (file_version
>= 79)
402 gmx_fio_do_int(fio
, simtemp
->eSimTempScale
);
403 gmx_fio_do_real(fio
, simtemp
->simtemp_high
);
404 gmx_fio_do_real(fio
, simtemp
->simtemp_low
);
409 snew(simtemp
->temperatures
, n_lambda
);
411 gmx_fio_ndo_real(fio
, simtemp
->temperatures
, n_lambda
);
416 static void do_imd(t_fileio
*fio
, t_IMD
*imd
, gmx_bool bRead
)
418 gmx_fio_do_int(fio
, imd
->nat
);
421 snew(imd
->ind
, imd
->nat
);
423 gmx_fio_ndo_int(fio
, imd
->ind
, imd
->nat
);
426 static void do_fepvals(t_fileio
*fio
, t_lambda
*fepvals
, gmx_bool bRead
, int file_version
)
428 /* i is defined in the ndo_double macro; use g to iterate. */
432 /* free energy values */
434 if (file_version
>= 79)
436 gmx_fio_do_int(fio
, fepvals
->init_fep_state
);
437 gmx_fio_do_double(fio
, fepvals
->init_lambda
);
438 gmx_fio_do_double(fio
, fepvals
->delta_lambda
);
440 else if (file_version
>= 59)
442 gmx_fio_do_double(fio
, fepvals
->init_lambda
);
443 gmx_fio_do_double(fio
, fepvals
->delta_lambda
);
447 gmx_fio_do_real(fio
, rdum
);
448 fepvals
->init_lambda
= rdum
;
449 gmx_fio_do_real(fio
, rdum
);
450 fepvals
->delta_lambda
= rdum
;
452 if (file_version
>= 79)
454 gmx_fio_do_int(fio
, fepvals
->n_lambda
);
457 snew(fepvals
->all_lambda
, efptNR
);
459 for (g
= 0; g
< efptNR
; g
++)
461 if (fepvals
->n_lambda
> 0)
465 snew(fepvals
->all_lambda
[g
], fepvals
->n_lambda
);
467 gmx_fio_ndo_double(fio
, fepvals
->all_lambda
[g
], fepvals
->n_lambda
);
468 gmx_fio_ndo_int(fio
, fepvals
->separate_dvdl
, efptNR
);
470 else if (fepvals
->init_lambda
>= 0)
472 fepvals
->separate_dvdl
[efptFEP
] = TRUE
;
476 else if (file_version
>= 64)
478 gmx_fio_do_int(fio
, fepvals
->n_lambda
);
483 snew(fepvals
->all_lambda
, efptNR
);
484 /* still allocate the all_lambda array's contents. */
485 for (g
= 0; g
< efptNR
; g
++)
487 if (fepvals
->n_lambda
> 0)
489 snew(fepvals
->all_lambda
[g
], fepvals
->n_lambda
);
493 gmx_fio_ndo_double(fio
, fepvals
->all_lambda
[efptFEP
],
495 if (fepvals
->init_lambda
>= 0)
499 fepvals
->separate_dvdl
[efptFEP
] = TRUE
;
503 /* copy the contents of the efptFEP lambda component to all
504 the other components */
505 for (g
= 0; g
< efptNR
; g
++)
507 for (h
= 0; h
< fepvals
->n_lambda
; h
++)
511 fepvals
->all_lambda
[g
][h
] =
512 fepvals
->all_lambda
[efptFEP
][h
];
521 fepvals
->n_lambda
= 0;
522 fepvals
->all_lambda
= NULL
;
523 if (fepvals
->init_lambda
>= 0)
525 fepvals
->separate_dvdl
[efptFEP
] = TRUE
;
528 if (file_version
>= 13)
530 gmx_fio_do_real(fio
, fepvals
->sc_alpha
);
534 fepvals
->sc_alpha
= 0;
536 if (file_version
>= 38)
538 gmx_fio_do_int(fio
, fepvals
->sc_power
);
542 fepvals
->sc_power
= 2;
544 if (file_version
>= 79)
546 gmx_fio_do_real(fio
, fepvals
->sc_r_power
);
550 fepvals
->sc_r_power
= 6.0;
552 if (file_version
>= 15)
554 gmx_fio_do_real(fio
, fepvals
->sc_sigma
);
558 fepvals
->sc_sigma
= 0.3;
562 if (file_version
>= 71)
564 fepvals
->sc_sigma_min
= fepvals
->sc_sigma
;
568 fepvals
->sc_sigma_min
= 0;
571 if (file_version
>= 79)
573 gmx_fio_do_int(fio
, fepvals
->bScCoul
);
577 fepvals
->bScCoul
= TRUE
;
579 if (file_version
>= 64)
581 gmx_fio_do_int(fio
, fepvals
->nstdhdl
);
585 fepvals
->nstdhdl
= 1;
588 if (file_version
>= 73)
590 gmx_fio_do_int(fio
, fepvals
->separate_dhdl_file
);
591 gmx_fio_do_int(fio
, fepvals
->dhdl_derivatives
);
595 fepvals
->separate_dhdl_file
= esepdhdlfileYES
;
596 fepvals
->dhdl_derivatives
= edhdlderivativesYES
;
598 if (file_version
>= 71)
600 gmx_fio_do_int(fio
, fepvals
->dh_hist_size
);
601 gmx_fio_do_double(fio
, fepvals
->dh_hist_spacing
);
605 fepvals
->dh_hist_size
= 0;
606 fepvals
->dh_hist_spacing
= 0.1;
608 if (file_version
>= 79)
610 gmx_fio_do_int(fio
, fepvals
->edHdLPrintEnergy
);
614 fepvals
->edHdLPrintEnergy
= edHdLPrintEnergyNO
;
617 /* handle lambda_neighbors */
618 if ((file_version
>= 83 && file_version
< 90) || file_version
>= 92)
620 gmx_fio_do_int(fio
, fepvals
->lambda_neighbors
);
621 if ( (fepvals
->lambda_neighbors
>= 0) && (fepvals
->init_fep_state
>= 0) &&
622 (fepvals
->init_lambda
< 0) )
624 fepvals
->lambda_start_n
= (fepvals
->init_fep_state
-
625 fepvals
->lambda_neighbors
);
626 fepvals
->lambda_stop_n
= (fepvals
->init_fep_state
+
627 fepvals
->lambda_neighbors
+ 1);
628 if (fepvals
->lambda_start_n
< 0)
630 fepvals
->lambda_start_n
= 0;;
632 if (fepvals
->lambda_stop_n
>= fepvals
->n_lambda
)
634 fepvals
->lambda_stop_n
= fepvals
->n_lambda
;
639 fepvals
->lambda_start_n
= 0;
640 fepvals
->lambda_stop_n
= fepvals
->n_lambda
;
645 fepvals
->lambda_start_n
= 0;
646 fepvals
->lambda_stop_n
= fepvals
->n_lambda
;
650 static void do_pull(t_fileio
*fio
, pull_params_t
*pull
, gmx_bool bRead
,
651 int file_version
, int ePullOld
)
657 if (file_version
>= 95)
659 gmx_fio_do_int(fio
, pull
->ngroup
);
661 gmx_fio_do_int(fio
, pull
->ncoord
);
662 if (file_version
< 95)
664 pull
->ngroup
= pull
->ncoord
+ 1;
666 if (file_version
< tpxv_PullCoordTypeGeom
)
670 gmx_fio_do_int(fio
, eGeomOld
);
671 gmx_fio_do_ivec(fio
, dimOld
);
672 /* The inner cylinder radius, now removed */
673 gmx_fio_do_real(fio
, dum
);
675 gmx_fio_do_real(fio
, pull
->cylinder_r
);
676 gmx_fio_do_real(fio
, pull
->constr_tol
);
677 if (file_version
>= 95)
679 gmx_fio_do_int(fio
, pull
->bPrintCOM1
);
680 /* With file_version < 95 this value is set below */
682 if (file_version
>= tpxv_PullCoordTypeGeom
)
684 gmx_fio_do_int(fio
, pull
->bPrintCOM2
);
685 gmx_fio_do_int(fio
, pull
->bPrintRefValue
);
686 gmx_fio_do_int(fio
, pull
->bPrintComp
);
690 pull
->bPrintCOM2
= FALSE
;
691 pull
->bPrintRefValue
= FALSE
;
692 pull
->bPrintComp
= TRUE
;
694 gmx_fio_do_int(fio
, pull
->nstxout
);
695 gmx_fio_do_int(fio
, pull
->nstfout
);
698 snew(pull
->group
, pull
->ngroup
);
699 snew(pull
->coord
, pull
->ncoord
);
701 if (file_version
< 95)
703 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
704 if (eGeomOld
== epullgDIRPBC
)
706 gmx_fatal(FARGS
, "pull-geometry=position is no longer supported");
708 if (eGeomOld
> epullgDIRPBC
)
713 for (g
= 0; g
< pull
->ngroup
; g
++)
715 /* We read and ignore a pull coordinate for group 0 */
716 do_pullgrp_tpx_pre95(fio
, &pull
->group
[g
], &pull
->coord
[std::max(g
-1, 0)],
717 bRead
, file_version
);
720 pull
->coord
[g
-1].group
[0] = 0;
721 pull
->coord
[g
-1].group
[1] = g
;
725 pull
->bPrintCOM1
= (pull
->group
[0].nat
> 0);
729 for (g
= 0; g
< pull
->ngroup
; g
++)
731 do_pull_group(fio
, &pull
->group
[g
], bRead
);
733 for (g
= 0; g
< pull
->ncoord
; g
++)
735 do_pull_coord(fio
, &pull
->coord
[g
],
736 file_version
, ePullOld
, eGeomOld
, dimOld
);
742 static void do_rotgrp(t_fileio
*fio
, t_rotgrp
*rotg
, gmx_bool bRead
)
744 gmx_fio_do_int(fio
, rotg
->eType
);
745 gmx_fio_do_int(fio
, rotg
->bMassW
);
746 gmx_fio_do_int(fio
, rotg
->nat
);
749 snew(rotg
->ind
, rotg
->nat
);
751 gmx_fio_ndo_int(fio
, rotg
->ind
, rotg
->nat
);
754 snew(rotg
->x_ref
, rotg
->nat
);
756 gmx_fio_ndo_rvec(fio
, rotg
->x_ref
, rotg
->nat
);
757 gmx_fio_do_rvec(fio
, rotg
->vec
);
758 gmx_fio_do_rvec(fio
, rotg
->pivot
);
759 gmx_fio_do_real(fio
, rotg
->rate
);
760 gmx_fio_do_real(fio
, rotg
->k
);
761 gmx_fio_do_real(fio
, rotg
->slab_dist
);
762 gmx_fio_do_real(fio
, rotg
->min_gaussian
);
763 gmx_fio_do_real(fio
, rotg
->eps
);
764 gmx_fio_do_int(fio
, rotg
->eFittype
);
765 gmx_fio_do_int(fio
, rotg
->PotAngle_nstep
);
766 gmx_fio_do_real(fio
, rotg
->PotAngle_step
);
769 static void do_rot(t_fileio
*fio
, t_rot
*rot
, gmx_bool bRead
)
773 gmx_fio_do_int(fio
, rot
->ngrp
);
774 gmx_fio_do_int(fio
, rot
->nstrout
);
775 gmx_fio_do_int(fio
, rot
->nstsout
);
778 snew(rot
->grp
, rot
->ngrp
);
780 for (g
= 0; g
< rot
->ngrp
; g
++)
782 do_rotgrp(fio
, &rot
->grp
[g
], bRead
);
787 static void do_swapgroup(t_fileio
*fio
, t_swapGroup
*g
, gmx_bool bRead
)
790 /* Name of the group or molecule */
795 gmx_fio_do_string(fio
, buf
);
796 g
->molname
= gmx_strdup(buf
);
800 gmx_fio_do_string(fio
, g
->molname
);
803 /* Number of atoms in the group */
804 gmx_fio_do_int(fio
, g
->nat
);
806 /* The group's atom indices */
809 snew(g
->ind
, g
->nat
);
811 gmx_fio_ndo_int(fio
, g
->ind
, g
->nat
);
813 /* Requested counts for compartments A and B */
814 gmx_fio_ndo_int(fio
, g
->nmolReq
, eCompNR
);
817 static void do_swapcoords_tpx(t_fileio
*fio
, t_swapcoords
*swap
, gmx_bool bRead
, int file_version
)
819 /* Enums for better readability of the code */
824 eChannel0
= 0, eChannel1
828 if (file_version
>= tpxv_CompElPolyatomicIonsAndMultipleIonTypes
)
830 /* The total number of swap groups is the sum of the fixed groups
831 * (split0, split1, solvent), and the user-defined groups (2+ types of ions)
833 gmx_fio_do_int(fio
, swap
->ngrp
);
836 snew(swap
->grp
, swap
->ngrp
);
838 for (int ig
= 0; ig
< swap
->ngrp
; ig
++)
840 do_swapgroup(fio
, &swap
->grp
[ig
], bRead
);
842 gmx_fio_do_int(fio
, swap
->massw_split
[eChannel0
]);
843 gmx_fio_do_int(fio
, swap
->massw_split
[eChannel1
]);
844 gmx_fio_do_int(fio
, swap
->nstswap
);
845 gmx_fio_do_int(fio
, swap
->nAverage
);
846 gmx_fio_do_real(fio
, swap
->threshold
);
847 gmx_fio_do_real(fio
, swap
->cyl0r
);
848 gmx_fio_do_real(fio
, swap
->cyl0u
);
849 gmx_fio_do_real(fio
, swap
->cyl0l
);
850 gmx_fio_do_real(fio
, swap
->cyl1r
);
851 gmx_fio_do_real(fio
, swap
->cyl1u
);
852 gmx_fio_do_real(fio
, swap
->cyl1l
);
856 /*** Support reading older CompEl .tpr files ***/
858 /* In the original CompEl .tpr files, we always have 5 groups: */
860 snew(swap
->grp
, swap
->ngrp
);
862 swap
->grp
[eGrpSplit0
].molname
= gmx_strdup("split0" ); // group 0: split0
863 swap
->grp
[eGrpSplit1
].molname
= gmx_strdup("split1" ); // group 1: split1
864 swap
->grp
[eGrpSolvent
].molname
= gmx_strdup("solvent"); // group 2: solvent
865 swap
->grp
[3 ].molname
= gmx_strdup("anions" ); // group 3: anions
866 swap
->grp
[4 ].molname
= gmx_strdup("cations"); // group 4: cations
868 gmx_fio_do_int(fio
, swap
->grp
[3].nat
);
869 gmx_fio_do_int(fio
, swap
->grp
[eGrpSolvent
].nat
);
870 gmx_fio_do_int(fio
, swap
->grp
[eGrpSplit0
].nat
);
871 gmx_fio_do_int(fio
, swap
->massw_split
[eChannel0
]);
872 gmx_fio_do_int(fio
, swap
->grp
[eGrpSplit1
].nat
);
873 gmx_fio_do_int(fio
, swap
->massw_split
[eChannel1
]);
874 gmx_fio_do_int(fio
, swap
->nstswap
);
875 gmx_fio_do_int(fio
, swap
->nAverage
);
876 gmx_fio_do_real(fio
, swap
->threshold
);
877 gmx_fio_do_real(fio
, swap
->cyl0r
);
878 gmx_fio_do_real(fio
, swap
->cyl0u
);
879 gmx_fio_do_real(fio
, swap
->cyl0l
);
880 gmx_fio_do_real(fio
, swap
->cyl1r
);
881 gmx_fio_do_real(fio
, swap
->cyl1u
);
882 gmx_fio_do_real(fio
, swap
->cyl1l
);
884 // The order[] array keeps compatibility with older .tpr files
885 // by reading in the groups in the classic order
887 const int order
[4] = {3, eGrpSolvent
, eGrpSplit0
, eGrpSplit1
};
889 for (int ig
= 0; ig
< 4; ig
++)
892 snew(swap
->grp
[g
].ind
, swap
->grp
[g
].nat
);
893 gmx_fio_ndo_int(fio
, swap
->grp
[g
].ind
, swap
->grp
[g
].nat
);
897 for (int j
= eCompA
; j
<= eCompB
; j
++)
899 gmx_fio_do_int(fio
, swap
->grp
[3].nmolReq
[j
]); // group 3 = anions
900 gmx_fio_do_int(fio
, swap
->grp
[4].nmolReq
[j
]); // group 4 = cations
902 } /* End support reading older CompEl .tpr files */
904 if (file_version
>= tpxv_CompElWithSwapLayerOffset
)
906 gmx_fio_do_real(fio
, swap
->bulkOffset
[eCompA
]);
907 gmx_fio_do_real(fio
, swap
->bulkOffset
[eCompB
]);
913 static void do_inputrec(t_fileio
*fio
, t_inputrec
*ir
, gmx_bool bRead
,
914 int file_version
, real
*fudgeQQ
)
916 int i
, j
, k
, *tmp
, idum
= 0;
919 gmx_bool bSimAnn
, bdum
= 0;
920 real zerotemptime
, finish_t
, init_temp
, finish_temp
;
922 if (file_version
!= tpx_version
)
924 /* Give a warning about features that are not accessible */
925 fprintf(stderr
, "Note: file tpx version %d, software tpx version %d\n",
926 file_version
, tpx_version
);
934 if (file_version
== 0)
939 /* Basic inputrec stuff */
940 gmx_fio_do_int(fio
, ir
->eI
);
941 if (file_version
>= 62)
943 gmx_fio_do_int64(fio
, ir
->nsteps
);
947 gmx_fio_do_int(fio
, idum
);
951 if (file_version
> 25)
953 if (file_version
>= 62)
955 gmx_fio_do_int64(fio
, ir
->init_step
);
959 gmx_fio_do_int(fio
, idum
);
960 ir
->init_step
= idum
;
968 if (file_version
>= 58)
970 gmx_fio_do_int(fio
, ir
->simulation_part
);
974 ir
->simulation_part
= 1;
977 if (file_version
>= 67)
979 gmx_fio_do_int(fio
, ir
->nstcalcenergy
);
983 ir
->nstcalcenergy
= 1;
985 if (file_version
< 53)
987 /* The pbc info has been moved out of do_inputrec,
988 * since we always want it, also without reading the inputrec.
990 gmx_fio_do_int(fio
, ir
->ePBC
);
991 if ((file_version
<= 15) && (ir
->ePBC
== 2))
995 if (file_version
>= 45)
997 gmx_fio_do_int(fio
, ir
->bPeriodicMols
);
1004 ir
->bPeriodicMols
= TRUE
;
1008 ir
->bPeriodicMols
= FALSE
;
1012 if (file_version
>= 81)
1014 gmx_fio_do_int(fio
, ir
->cutoff_scheme
);
1015 if (file_version
< 94)
1017 ir
->cutoff_scheme
= 1 - ir
->cutoff_scheme
;
1022 ir
->cutoff_scheme
= ecutsGROUP
;
1024 gmx_fio_do_int(fio
, ir
->ns_type
);
1025 gmx_fio_do_int(fio
, ir
->nstlist
);
1026 gmx_fio_do_int(fio
, idum
); /* used to be ndelta; not used anymore */
1027 if (file_version
< 41)
1029 gmx_fio_do_int(fio
, idum
);
1030 gmx_fio_do_int(fio
, idum
);
1032 if (file_version
>= 45)
1034 gmx_fio_do_real(fio
, ir
->rtpi
);
1040 gmx_fio_do_int(fio
, ir
->nstcomm
);
1041 if (file_version
> 34)
1043 gmx_fio_do_int(fio
, ir
->comm_mode
);
1045 else if (ir
->nstcomm
< 0)
1047 ir
->comm_mode
= ecmANGULAR
;
1051 ir
->comm_mode
= ecmLINEAR
;
1053 ir
->nstcomm
= abs(ir
->nstcomm
);
1055 /* ignore nstcheckpoint */
1056 if (file_version
> 25 && file_version
< tpxv_RemoveObsoleteParameters1
)
1058 gmx_fio_do_int(fio
, idum
);
1061 gmx_fio_do_int(fio
, ir
->nstcgsteep
);
1063 if (file_version
>= 30)
1065 gmx_fio_do_int(fio
, ir
->nbfgscorr
);
1072 gmx_fio_do_int(fio
, ir
->nstlog
);
1073 gmx_fio_do_int(fio
, ir
->nstxout
);
1074 gmx_fio_do_int(fio
, ir
->nstvout
);
1075 gmx_fio_do_int(fio
, ir
->nstfout
);
1076 gmx_fio_do_int(fio
, ir
->nstenergy
);
1077 gmx_fio_do_int(fio
, ir
->nstxout_compressed
);
1078 if (file_version
>= 59)
1080 gmx_fio_do_double(fio
, ir
->init_t
);
1081 gmx_fio_do_double(fio
, ir
->delta_t
);
1085 gmx_fio_do_real(fio
, rdum
);
1087 gmx_fio_do_real(fio
, rdum
);
1090 gmx_fio_do_real(fio
, ir
->x_compression_precision
);
1091 if (file_version
< 19)
1093 gmx_fio_do_int(fio
, idum
);
1094 gmx_fio_do_int(fio
, idum
);
1096 if (file_version
< 18)
1098 gmx_fio_do_int(fio
, idum
);
1100 if (file_version
>= 81)
1102 gmx_fio_do_real(fio
, ir
->verletbuf_tol
);
1106 ir
->verletbuf_tol
= 0;
1108 gmx_fio_do_real(fio
, ir
->rlist
);
1109 if (file_version
>= 67 && file_version
< tpxv_RemoveTwinRange
)
1113 // Reading such a file version could invoke the twin-range
1114 // scheme, about which mdrun should give a fatal error.
1115 real dummy_rlistlong
= -1;
1116 gmx_fio_do_real(fio
, dummy_rlistlong
);
1118 if (ir
->rlist
> 0 && (dummy_rlistlong
== 0 || dummy_rlistlong
> ir
->rlist
))
1120 // Get mdrun to issue an error (regardless of
1121 // ir->cutoff_scheme).
1122 ir
->useTwinRange
= true;
1126 // grompp used to set rlistlong actively. Users were
1127 // probably also confused and set rlistlong == rlist.
1128 // However, in all remaining cases, it is safe to let
1129 // mdrun proceed normally.
1130 ir
->useTwinRange
= false;
1136 // No need to read or write anything
1137 ir
->useTwinRange
= false;
1139 if (file_version
>= 82 && file_version
!= 90)
1141 // Multiple time-stepping is no longer enabled, but the old
1142 // support required the twin-range scheme, for which mdrun
1143 // already emits a fatal error.
1144 int dummy_nstcalclr
= -1;
1145 gmx_fio_do_int(fio
, dummy_nstcalclr
);
1147 gmx_fio_do_int(fio
, ir
->coulombtype
);
1148 if (file_version
< 32 && ir
->coulombtype
== eelRF
)
1150 ir
->coulombtype
= eelRF_NEC_UNSUPPORTED
;
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 if (file_version
>= 37)
1177 gmx_fio_do_real(fio
, ir
->epsilon_rf
);
1181 if (EEL_RF(ir
->coulombtype
))
1183 ir
->epsilon_rf
= ir
->epsilon_r
;
1184 ir
->epsilon_r
= 1.0;
1188 ir
->epsilon_rf
= 1.0;
1191 if (file_version
>= 29)
1193 gmx_fio_do_real(fio
, ir
->tabext
);
1200 if (file_version
> 25)
1202 gmx_fio_do_int(fio
, ir
->gb_algorithm
);
1203 gmx_fio_do_int(fio
, ir
->nstgbradii
);
1204 gmx_fio_do_real(fio
, ir
->rgbradii
);
1205 gmx_fio_do_real(fio
, ir
->gb_saltconc
);
1206 gmx_fio_do_int(fio
, ir
->implicit_solvent
);
1210 ir
->gb_algorithm
= egbSTILL
;
1213 ir
->gb_saltconc
= 0;
1214 ir
->implicit_solvent
= eisNO
;
1216 if (file_version
>= 55)
1218 gmx_fio_do_real(fio
, ir
->gb_epsilon_solvent
);
1219 gmx_fio_do_real(fio
, ir
->gb_obc_alpha
);
1220 gmx_fio_do_real(fio
, ir
->gb_obc_beta
);
1221 gmx_fio_do_real(fio
, ir
->gb_obc_gamma
);
1222 if (file_version
>= 60)
1224 gmx_fio_do_real(fio
, ir
->gb_dielectric_offset
);
1225 gmx_fio_do_int(fio
, ir
->sa_algorithm
);
1229 ir
->gb_dielectric_offset
= 0.009;
1230 ir
->sa_algorithm
= esaAPPROX
;
1232 gmx_fio_do_real(fio
, ir
->sa_surface_tension
);
1234 /* Override sa_surface_tension if it is not changed in the mpd-file */
1235 if (ir
->sa_surface_tension
< 0)
1237 if (ir
->gb_algorithm
== egbSTILL
)
1239 ir
->sa_surface_tension
= 0.0049 * 100 * CAL2JOULE
;
1241 else if (ir
->gb_algorithm
== egbHCT
|| ir
->gb_algorithm
== egbOBC
)
1243 ir
->sa_surface_tension
= 0.0054 * 100 * CAL2JOULE
;
1250 /* Better use sensible values than insane (0.0) ones... */
1251 ir
->gb_epsilon_solvent
= 80;
1252 ir
->gb_obc_alpha
= 1.0;
1253 ir
->gb_obc_beta
= 0.8;
1254 ir
->gb_obc_gamma
= 4.85;
1255 ir
->sa_surface_tension
= 2.092;
1259 if (file_version
>= 81)
1261 gmx_fio_do_real(fio
, ir
->fourier_spacing
);
1265 ir
->fourier_spacing
= 0.0;
1267 gmx_fio_do_int(fio
, ir
->nkx
);
1268 gmx_fio_do_int(fio
, ir
->nky
);
1269 gmx_fio_do_int(fio
, ir
->nkz
);
1270 gmx_fio_do_int(fio
, ir
->pme_order
);
1271 gmx_fio_do_real(fio
, ir
->ewald_rtol
);
1273 if (file_version
>= 93)
1275 gmx_fio_do_real(fio
, ir
->ewald_rtol_lj
);
1279 ir
->ewald_rtol_lj
= ir
->ewald_rtol
;
1282 if (file_version
>= 24)
1284 gmx_fio_do_int(fio
, ir
->ewald_geometry
);
1288 ir
->ewald_geometry
= eewg3D
;
1291 if (file_version
<= 17)
1293 ir
->epsilon_surface
= 0;
1294 if (file_version
== 17)
1296 gmx_fio_do_int(fio
, idum
);
1301 gmx_fio_do_real(fio
, ir
->epsilon_surface
);
1304 /* ignore bOptFFT */
1305 if (file_version
< tpxv_RemoveObsoleteParameters1
)
1307 gmx_fio_do_gmx_bool(fio
, bdum
);
1310 if (file_version
>= 93)
1312 gmx_fio_do_int(fio
, ir
->ljpme_combination_rule
);
1314 gmx_fio_do_gmx_bool(fio
, ir
->bContinuation
);
1315 gmx_fio_do_int(fio
, ir
->etc
);
1316 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1317 * but the values 0 and 1 still mean no and
1318 * berendsen temperature coupling, respectively.
1320 if (file_version
>= 79)
1322 gmx_fio_do_gmx_bool(fio
, ir
->bPrintNHChains
);
1324 if (file_version
>= 71)
1326 gmx_fio_do_int(fio
, ir
->nsttcouple
);
1330 ir
->nsttcouple
= ir
->nstcalcenergy
;
1332 if (file_version
<= 15)
1334 gmx_fio_do_int(fio
, idum
);
1336 if (file_version
<= 17)
1338 gmx_fio_do_int(fio
, ir
->epct
);
1339 if (file_version
<= 15)
1343 ir
->epct
= epctSURFACETENSION
;
1345 gmx_fio_do_int(fio
, idum
);
1348 /* we have removed the NO alternative at the beginning */
1352 ir
->epct
= epctISOTROPIC
;
1356 ir
->epc
= epcBERENDSEN
;
1361 gmx_fio_do_int(fio
, ir
->epc
);
1362 gmx_fio_do_int(fio
, ir
->epct
);
1364 if (file_version
>= 71)
1366 gmx_fio_do_int(fio
, ir
->nstpcouple
);
1370 ir
->nstpcouple
= ir
->nstcalcenergy
;
1372 gmx_fio_do_real(fio
, ir
->tau_p
);
1373 if (file_version
<= 15)
1375 gmx_fio_do_rvec(fio
, vdum
);
1376 clear_mat(ir
->ref_p
);
1377 for (i
= 0; i
< DIM
; i
++)
1379 ir
->ref_p
[i
][i
] = vdum
[i
];
1384 gmx_fio_do_rvec(fio
, ir
->ref_p
[XX
]);
1385 gmx_fio_do_rvec(fio
, ir
->ref_p
[YY
]);
1386 gmx_fio_do_rvec(fio
, ir
->ref_p
[ZZ
]);
1388 if (file_version
<= 15)
1390 gmx_fio_do_rvec(fio
, vdum
);
1391 clear_mat(ir
->compress
);
1392 for (i
= 0; i
< DIM
; i
++)
1394 ir
->compress
[i
][i
] = vdum
[i
];
1399 gmx_fio_do_rvec(fio
, ir
->compress
[XX
]);
1400 gmx_fio_do_rvec(fio
, ir
->compress
[YY
]);
1401 gmx_fio_do_rvec(fio
, ir
->compress
[ZZ
]);
1403 if (file_version
>= 47)
1405 gmx_fio_do_int(fio
, ir
->refcoord_scaling
);
1406 gmx_fio_do_rvec(fio
, ir
->posres_com
);
1407 gmx_fio_do_rvec(fio
, ir
->posres_comB
);
1411 ir
->refcoord_scaling
= erscNO
;
1412 clear_rvec(ir
->posres_com
);
1413 clear_rvec(ir
->posres_comB
);
1415 if ((file_version
> 25) && (file_version
< 79))
1417 gmx_fio_do_int(fio
, ir
->andersen_seed
);
1421 ir
->andersen_seed
= 0;
1423 if (file_version
< 26)
1425 gmx_fio_do_gmx_bool(fio
, bSimAnn
);
1426 gmx_fio_do_real(fio
, zerotemptime
);
1429 if (file_version
< 37)
1431 gmx_fio_do_real(fio
, rdum
);
1434 gmx_fio_do_real(fio
, ir
->shake_tol
);
1435 if (file_version
< 54)
1437 gmx_fio_do_real(fio
, *fudgeQQ
);
1440 gmx_fio_do_int(fio
, ir
->efep
);
1441 if (file_version
<= 14 && ir
->efep
!= efepNO
)
1445 do_fepvals(fio
, ir
->fepvals
, bRead
, file_version
);
1447 if (file_version
>= 79)
1449 gmx_fio_do_gmx_bool(fio
, ir
->bSimTemp
);
1452 ir
->bSimTemp
= TRUE
;
1457 ir
->bSimTemp
= FALSE
;
1461 do_simtempvals(fio
, ir
->simtempvals
, ir
->fepvals
->n_lambda
, bRead
, file_version
);
1464 if (file_version
>= 79)
1466 gmx_fio_do_gmx_bool(fio
, ir
->bExpanded
);
1469 ir
->bExpanded
= TRUE
;
1473 ir
->bExpanded
= FALSE
;
1478 do_expandedvals(fio
, ir
->expandedvals
, ir
->fepvals
, bRead
, file_version
);
1480 if (file_version
>= 57)
1482 gmx_fio_do_int(fio
, ir
->eDisre
);
1484 gmx_fio_do_int(fio
, ir
->eDisreWeighting
);
1485 if (file_version
< 22)
1487 if (ir
->eDisreWeighting
== 0)
1489 ir
->eDisreWeighting
= edrwEqual
;
1493 ir
->eDisreWeighting
= edrwConservative
;
1496 gmx_fio_do_gmx_bool(fio
, ir
->bDisreMixed
);
1497 gmx_fio_do_real(fio
, ir
->dr_fc
);
1498 gmx_fio_do_real(fio
, ir
->dr_tau
);
1499 gmx_fio_do_int(fio
, ir
->nstdisreout
);
1500 if (file_version
>= 22)
1502 gmx_fio_do_real(fio
, ir
->orires_fc
);
1503 gmx_fio_do_real(fio
, ir
->orires_tau
);
1504 gmx_fio_do_int(fio
, ir
->nstorireout
);
1510 ir
->nstorireout
= 0;
1513 /* ignore dihre_fc */
1514 if (file_version
>= 26 && file_version
< 79)
1516 gmx_fio_do_real(fio
, rdum
);
1517 if (file_version
< 56)
1519 gmx_fio_do_real(fio
, rdum
);
1520 gmx_fio_do_int(fio
, idum
);
1524 gmx_fio_do_real(fio
, ir
->em_stepsize
);
1525 gmx_fio_do_real(fio
, ir
->em_tol
);
1526 if (file_version
>= 22)
1528 gmx_fio_do_gmx_bool(fio
, ir
->bShakeSOR
);
1532 ir
->bShakeSOR
= TRUE
;
1534 if (file_version
>= 11)
1536 gmx_fio_do_int(fio
, ir
->niter
);
1541 fprintf(stderr
, "Note: niter not in run input file, setting it to %d\n",
1544 if (file_version
>= 21)
1546 gmx_fio_do_real(fio
, ir
->fc_stepsize
);
1550 ir
->fc_stepsize
= 0;
1552 gmx_fio_do_int(fio
, ir
->eConstrAlg
);
1553 gmx_fio_do_int(fio
, ir
->nProjOrder
);
1554 gmx_fio_do_real(fio
, ir
->LincsWarnAngle
);
1555 if (file_version
<= 14)
1557 gmx_fio_do_int(fio
, idum
);
1559 if (file_version
>= 26)
1561 gmx_fio_do_int(fio
, ir
->nLincsIter
);
1566 fprintf(stderr
, "Note: nLincsIter not in run input file, setting it to %d\n",
1569 if (file_version
< 33)
1571 gmx_fio_do_real(fio
, bd_temp
);
1573 gmx_fio_do_real(fio
, ir
->bd_fric
);
1574 if (file_version
>= tpxv_Use64BitRandomSeed
)
1576 gmx_fio_do_int64(fio
, ir
->ld_seed
);
1580 gmx_fio_do_int(fio
, idum
);
1583 if (file_version
>= 33)
1585 for (i
= 0; i
< DIM
; i
++)
1587 gmx_fio_do_rvec(fio
, ir
->deform
[i
]);
1592 for (i
= 0; i
< DIM
; i
++)
1594 clear_rvec(ir
->deform
[i
]);
1597 if (file_version
>= 14)
1599 gmx_fio_do_real(fio
, ir
->cos_accel
);
1605 gmx_fio_do_int(fio
, ir
->userint1
);
1606 gmx_fio_do_int(fio
, ir
->userint2
);
1607 gmx_fio_do_int(fio
, ir
->userint3
);
1608 gmx_fio_do_int(fio
, ir
->userint4
);
1609 gmx_fio_do_real(fio
, ir
->userreal1
);
1610 gmx_fio_do_real(fio
, ir
->userreal2
);
1611 gmx_fio_do_real(fio
, ir
->userreal3
);
1612 gmx_fio_do_real(fio
, ir
->userreal4
);
1614 /* AdResS is removed, but we need to be able to read old files,
1615 and let mdrun refuse to run them */
1616 if (file_version
>= 77 && file_version
< tpxv_RemoveAdress
)
1618 gmx_fio_do_gmx_bool(fio
, ir
->bAdress
);
1621 int idum
, numThermoForceGroups
, numEnergyGroups
;
1624 gmx_fio_do_int(fio
, idum
);
1625 gmx_fio_do_real(fio
, rdum
);
1626 gmx_fio_do_real(fio
, rdum
);
1627 gmx_fio_do_real(fio
, rdum
);
1628 gmx_fio_do_int(fio
, idum
);
1629 gmx_fio_do_int(fio
, idum
);
1630 gmx_fio_do_rvec(fio
, rvecdum
);
1631 gmx_fio_do_int(fio
, numThermoForceGroups
);
1632 gmx_fio_do_real(fio
, rdum
);
1633 gmx_fio_do_int(fio
, numEnergyGroups
);
1634 gmx_fio_do_int(fio
, idum
);
1636 if (numThermoForceGroups
> 0)
1638 std::vector
<int> idumn(numThermoForceGroups
);
1639 gmx_fio_ndo_int(fio
, idumn
.data(), idumn
.size());
1641 if (numEnergyGroups
> 0)
1643 std::vector
<int> idumn(numEnergyGroups
);
1644 gmx_fio_ndo_int(fio
, idumn
.data(), idumn
.size());
1650 ir
->bAdress
= FALSE
;
1654 if (file_version
>= 48)
1658 if (file_version
>= tpxv_PullCoordTypeGeom
)
1660 gmx_fio_do_gmx_bool(fio
, ir
->bPull
);
1664 gmx_fio_do_int(fio
, ePullOld
);
1665 ir
->bPull
= (ePullOld
> 0);
1666 /* We removed the first ePull=ePullNo for the enum */
1675 do_pull(fio
, ir
->pull
, bRead
, file_version
, ePullOld
);
1683 /* Enforced rotation */
1684 if (file_version
>= 74)
1686 gmx_fio_do_int(fio
, ir
->bRot
);
1687 if (ir
->bRot
== TRUE
)
1693 do_rot(fio
, ir
->rot
, bRead
);
1701 /* Interactive molecular dynamics */
1702 if (file_version
>= tpxv_InteractiveMolecularDynamics
)
1704 gmx_fio_do_int(fio
, ir
->bIMD
);
1705 if (TRUE
== ir
->bIMD
)
1711 do_imd(fio
, ir
->imd
, bRead
);
1716 /* We don't support IMD sessions for old .tpr files */
1721 gmx_fio_do_int(fio
, ir
->opts
.ngtc
);
1722 if (file_version
>= 69)
1724 gmx_fio_do_int(fio
, ir
->opts
.nhchainlength
);
1728 ir
->opts
.nhchainlength
= 1;
1730 gmx_fio_do_int(fio
, ir
->opts
.ngacc
);
1731 gmx_fio_do_int(fio
, ir
->opts
.ngfrz
);
1732 gmx_fio_do_int(fio
, ir
->opts
.ngener
);
1736 snew(ir
->opts
.nrdf
, ir
->opts
.ngtc
);
1737 snew(ir
->opts
.ref_t
, ir
->opts
.ngtc
);
1738 snew(ir
->opts
.annealing
, ir
->opts
.ngtc
);
1739 snew(ir
->opts
.anneal_npoints
, ir
->opts
.ngtc
);
1740 snew(ir
->opts
.anneal_time
, ir
->opts
.ngtc
);
1741 snew(ir
->opts
.anneal_temp
, ir
->opts
.ngtc
);
1742 snew(ir
->opts
.tau_t
, ir
->opts
.ngtc
);
1743 snew(ir
->opts
.nFreeze
, ir
->opts
.ngfrz
);
1744 snew(ir
->opts
.acc
, ir
->opts
.ngacc
);
1745 snew(ir
->opts
.egp_flags
, ir
->opts
.ngener
*ir
->opts
.ngener
);
1747 if (ir
->opts
.ngtc
> 0)
1749 if (bRead
&& file_version
< 13)
1751 snew(tmp
, ir
->opts
.ngtc
);
1752 gmx_fio_ndo_int(fio
, tmp
, ir
->opts
.ngtc
);
1753 for (i
= 0; i
< ir
->opts
.ngtc
; i
++)
1755 ir
->opts
.nrdf
[i
] = tmp
[i
];
1761 gmx_fio_ndo_real(fio
, ir
->opts
.nrdf
, ir
->opts
.ngtc
);
1763 gmx_fio_ndo_real(fio
, ir
->opts
.ref_t
, ir
->opts
.ngtc
);
1764 gmx_fio_ndo_real(fio
, ir
->opts
.tau_t
, ir
->opts
.ngtc
);
1765 if (file_version
< 33 && ir
->eI
== eiBD
)
1767 for (i
= 0; i
< ir
->opts
.ngtc
; i
++)
1769 ir
->opts
.tau_t
[i
] = bd_temp
;
1773 if (ir
->opts
.ngfrz
> 0)
1775 gmx_fio_ndo_ivec(fio
, ir
->opts
.nFreeze
, ir
->opts
.ngfrz
);
1777 if (ir
->opts
.ngacc
> 0)
1779 gmx_fio_ndo_rvec(fio
, ir
->opts
.acc
, ir
->opts
.ngacc
);
1781 if (file_version
>= 12)
1783 gmx_fio_ndo_int(fio
, ir
->opts
.egp_flags
,
1784 ir
->opts
.ngener
*ir
->opts
.ngener
);
1787 if (bRead
&& file_version
< 26)
1789 for (i
= 0; i
< ir
->opts
.ngtc
; i
++)
1793 ir
->opts
.annealing
[i
] = eannSINGLE
;
1794 ir
->opts
.anneal_npoints
[i
] = 2;
1795 snew(ir
->opts
.anneal_time
[i
], 2);
1796 snew(ir
->opts
.anneal_temp
[i
], 2);
1797 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1798 finish_t
= ir
->init_t
+ ir
->nsteps
* ir
->delta_t
;
1799 init_temp
= ir
->opts
.ref_t
[i
]*(1-ir
->init_t
/zerotemptime
);
1800 finish_temp
= ir
->opts
.ref_t
[i
]*(1-finish_t
/zerotemptime
);
1801 ir
->opts
.anneal_time
[i
][0] = ir
->init_t
;
1802 ir
->opts
.anneal_time
[i
][1] = finish_t
;
1803 ir
->opts
.anneal_temp
[i
][0] = init_temp
;
1804 ir
->opts
.anneal_temp
[i
][1] = finish_temp
;
1808 ir
->opts
.annealing
[i
] = eannNO
;
1809 ir
->opts
.anneal_npoints
[i
] = 0;
1815 /* file version 26 or later */
1816 /* First read the lists with annealing and npoints for each group */
1817 gmx_fio_ndo_int(fio
, ir
->opts
.annealing
, ir
->opts
.ngtc
);
1818 gmx_fio_ndo_int(fio
, ir
->opts
.anneal_npoints
, ir
->opts
.ngtc
);
1819 for (j
= 0; j
< (ir
->opts
.ngtc
); j
++)
1821 k
= ir
->opts
.anneal_npoints
[j
];
1824 snew(ir
->opts
.anneal_time
[j
], k
);
1825 snew(ir
->opts
.anneal_temp
[j
], k
);
1827 gmx_fio_ndo_real(fio
, ir
->opts
.anneal_time
[j
], k
);
1828 gmx_fio_ndo_real(fio
, ir
->opts
.anneal_temp
[j
], k
);
1832 if (file_version
>= 45)
1834 gmx_fio_do_int(fio
, ir
->nwall
);
1835 gmx_fio_do_int(fio
, ir
->wall_type
);
1836 if (file_version
>= 50)
1838 gmx_fio_do_real(fio
, ir
->wall_r_linpot
);
1842 ir
->wall_r_linpot
= -1;
1844 gmx_fio_do_int(fio
, ir
->wall_atomtype
[0]);
1845 gmx_fio_do_int(fio
, ir
->wall_atomtype
[1]);
1846 gmx_fio_do_real(fio
, ir
->wall_density
[0]);
1847 gmx_fio_do_real(fio
, ir
->wall_density
[1]);
1848 gmx_fio_do_real(fio
, ir
->wall_ewald_zfac
);
1854 ir
->wall_atomtype
[0] = -1;
1855 ir
->wall_atomtype
[1] = -1;
1856 ir
->wall_density
[0] = 0;
1857 ir
->wall_density
[1] = 0;
1858 ir
->wall_ewald_zfac
= 3;
1860 /* Cosine stuff for electric fields */
1861 for (j
= 0; (j
< DIM
); j
++)
1863 gmx_fio_do_int(fio
, ir
->ex
[j
].n
);
1864 gmx_fio_do_int(fio
, ir
->et
[j
].n
);
1867 snew(ir
->ex
[j
].a
, ir
->ex
[j
].n
);
1868 snew(ir
->ex
[j
].phi
, ir
->ex
[j
].n
);
1869 snew(ir
->et
[j
].a
, ir
->et
[j
].n
);
1870 snew(ir
->et
[j
].phi
, ir
->et
[j
].n
);
1872 gmx_fio_ndo_real(fio
, ir
->ex
[j
].a
, ir
->ex
[j
].n
);
1873 gmx_fio_ndo_real(fio
, ir
->ex
[j
].phi
, ir
->ex
[j
].n
);
1874 gmx_fio_ndo_real(fio
, ir
->et
[j
].a
, ir
->et
[j
].n
);
1875 gmx_fio_ndo_real(fio
, ir
->et
[j
].phi
, ir
->et
[j
].n
);
1879 if (file_version
>= tpxv_ComputationalElectrophysiology
)
1881 gmx_fio_do_int(fio
, ir
->eSwapCoords
);
1882 if (ir
->eSwapCoords
!= eswapNO
)
1888 do_swapcoords_tpx(fio
, ir
->swap
, bRead
, file_version
);
1893 if (file_version
>= 39)
1895 gmx_fio_do_gmx_bool(fio
, ir
->bQMMM
);
1896 gmx_fio_do_int(fio
, ir
->QMMMscheme
);
1897 gmx_fio_do_real(fio
, ir
->scalefactor
);
1898 gmx_fio_do_int(fio
, ir
->opts
.ngQM
);
1901 snew(ir
->opts
.QMmethod
, ir
->opts
.ngQM
);
1902 snew(ir
->opts
.QMbasis
, ir
->opts
.ngQM
);
1903 snew(ir
->opts
.QMcharge
, ir
->opts
.ngQM
);
1904 snew(ir
->opts
.QMmult
, ir
->opts
.ngQM
);
1905 snew(ir
->opts
.bSH
, ir
->opts
.ngQM
);
1906 snew(ir
->opts
.CASorbitals
, ir
->opts
.ngQM
);
1907 snew(ir
->opts
.CASelectrons
, ir
->opts
.ngQM
);
1908 snew(ir
->opts
.SAon
, ir
->opts
.ngQM
);
1909 snew(ir
->opts
.SAoff
, ir
->opts
.ngQM
);
1910 snew(ir
->opts
.SAsteps
, ir
->opts
.ngQM
);
1911 snew(ir
->opts
.bOPT
, ir
->opts
.ngQM
);
1912 snew(ir
->opts
.bTS
, ir
->opts
.ngQM
);
1914 if (ir
->opts
.ngQM
> 0)
1916 gmx_fio_ndo_int(fio
, ir
->opts
.QMmethod
, ir
->opts
.ngQM
);
1917 gmx_fio_ndo_int(fio
, ir
->opts
.QMbasis
, ir
->opts
.ngQM
);
1918 gmx_fio_ndo_int(fio
, ir
->opts
.QMcharge
, ir
->opts
.ngQM
);
1919 gmx_fio_ndo_int(fio
, ir
->opts
.QMmult
, ir
->opts
.ngQM
);
1920 gmx_fio_ndo_gmx_bool(fio
, ir
->opts
.bSH
, ir
->opts
.ngQM
);
1921 gmx_fio_ndo_int(fio
, ir
->opts
.CASorbitals
, ir
->opts
.ngQM
);
1922 gmx_fio_ndo_int(fio
, ir
->opts
.CASelectrons
, ir
->opts
.ngQM
);
1923 gmx_fio_ndo_real(fio
, ir
->opts
.SAon
, ir
->opts
.ngQM
);
1924 gmx_fio_ndo_real(fio
, ir
->opts
.SAoff
, ir
->opts
.ngQM
);
1925 gmx_fio_ndo_int(fio
, ir
->opts
.SAsteps
, ir
->opts
.ngQM
);
1926 gmx_fio_ndo_gmx_bool(fio
, ir
->opts
.bOPT
, ir
->opts
.ngQM
);
1927 gmx_fio_ndo_gmx_bool(fio
, ir
->opts
.bTS
, ir
->opts
.ngQM
);
1929 /* end of QMMM stuff */
1934 static void do_harm(t_fileio
*fio
, t_iparams
*iparams
)
1936 gmx_fio_do_real(fio
, iparams
->harmonic
.rA
);
1937 gmx_fio_do_real(fio
, iparams
->harmonic
.krA
);
1938 gmx_fio_do_real(fio
, iparams
->harmonic
.rB
);
1939 gmx_fio_do_real(fio
, iparams
->harmonic
.krB
);
1942 void do_iparams(t_fileio
*fio
, t_functype ftype
, t_iparams
*iparams
,
1943 gmx_bool bRead
, int file_version
)
1956 do_harm(fio
, iparams
);
1957 if ((ftype
== F_ANGRES
|| ftype
== F_ANGRESZ
) && bRead
)
1959 /* Correct incorrect storage of parameters */
1960 iparams
->pdihs
.phiB
= iparams
->pdihs
.phiA
;
1961 iparams
->pdihs
.cpB
= iparams
->pdihs
.cpA
;
1965 gmx_fio_do_real(fio
, iparams
->harmonic
.rA
);
1966 gmx_fio_do_real(fio
, iparams
->harmonic
.krA
);
1968 case F_LINEAR_ANGLES
:
1969 gmx_fio_do_real(fio
, iparams
->linangle
.klinA
);
1970 gmx_fio_do_real(fio
, iparams
->linangle
.aA
);
1971 gmx_fio_do_real(fio
, iparams
->linangle
.klinB
);
1972 gmx_fio_do_real(fio
, iparams
->linangle
.aB
);
1975 gmx_fio_do_real(fio
, iparams
->fene
.bm
);
1976 gmx_fio_do_real(fio
, iparams
->fene
.kb
);
1980 gmx_fio_do_real(fio
, iparams
->restraint
.lowA
);
1981 gmx_fio_do_real(fio
, iparams
->restraint
.up1A
);
1982 gmx_fio_do_real(fio
, iparams
->restraint
.up2A
);
1983 gmx_fio_do_real(fio
, iparams
->restraint
.kA
);
1984 gmx_fio_do_real(fio
, iparams
->restraint
.lowB
);
1985 gmx_fio_do_real(fio
, iparams
->restraint
.up1B
);
1986 gmx_fio_do_real(fio
, iparams
->restraint
.up2B
);
1987 gmx_fio_do_real(fio
, iparams
->restraint
.kB
);
1993 gmx_fio_do_real(fio
, iparams
->tab
.kA
);
1994 gmx_fio_do_int(fio
, iparams
->tab
.table
);
1995 gmx_fio_do_real(fio
, iparams
->tab
.kB
);
1997 case F_CROSS_BOND_BONDS
:
1998 gmx_fio_do_real(fio
, iparams
->cross_bb
.r1e
);
1999 gmx_fio_do_real(fio
, iparams
->cross_bb
.r2e
);
2000 gmx_fio_do_real(fio
, iparams
->cross_bb
.krr
);
2002 case F_CROSS_BOND_ANGLES
:
2003 gmx_fio_do_real(fio
, iparams
->cross_ba
.r1e
);
2004 gmx_fio_do_real(fio
, iparams
->cross_ba
.r2e
);
2005 gmx_fio_do_real(fio
, iparams
->cross_ba
.r3e
);
2006 gmx_fio_do_real(fio
, iparams
->cross_ba
.krt
);
2008 case F_UREY_BRADLEY
:
2009 gmx_fio_do_real(fio
, iparams
->u_b
.thetaA
);
2010 gmx_fio_do_real(fio
, iparams
->u_b
.kthetaA
);
2011 gmx_fio_do_real(fio
, iparams
->u_b
.r13A
);
2012 gmx_fio_do_real(fio
, iparams
->u_b
.kUBA
);
2013 if (file_version
>= 79)
2015 gmx_fio_do_real(fio
, iparams
->u_b
.thetaB
);
2016 gmx_fio_do_real(fio
, iparams
->u_b
.kthetaB
);
2017 gmx_fio_do_real(fio
, iparams
->u_b
.r13B
);
2018 gmx_fio_do_real(fio
, iparams
->u_b
.kUBB
);
2022 iparams
->u_b
.thetaB
= iparams
->u_b
.thetaA
;
2023 iparams
->u_b
.kthetaB
= iparams
->u_b
.kthetaA
;
2024 iparams
->u_b
.r13B
= iparams
->u_b
.r13A
;
2025 iparams
->u_b
.kUBB
= iparams
->u_b
.kUBA
;
2028 case F_QUARTIC_ANGLES
:
2029 gmx_fio_do_real(fio
, iparams
->qangle
.theta
);
2030 gmx_fio_ndo_real(fio
, iparams
->qangle
.c
, 5);
2033 gmx_fio_do_real(fio
, iparams
->bham
.a
);
2034 gmx_fio_do_real(fio
, iparams
->bham
.b
);
2035 gmx_fio_do_real(fio
, iparams
->bham
.c
);
2038 gmx_fio_do_real(fio
, iparams
->morse
.b0A
);
2039 gmx_fio_do_real(fio
, iparams
->morse
.cbA
);
2040 gmx_fio_do_real(fio
, iparams
->morse
.betaA
);
2041 if (file_version
>= 79)
2043 gmx_fio_do_real(fio
, iparams
->morse
.b0B
);
2044 gmx_fio_do_real(fio
, iparams
->morse
.cbB
);
2045 gmx_fio_do_real(fio
, iparams
->morse
.betaB
);
2049 iparams
->morse
.b0B
= iparams
->morse
.b0A
;
2050 iparams
->morse
.cbB
= iparams
->morse
.cbA
;
2051 iparams
->morse
.betaB
= iparams
->morse
.betaA
;
2055 gmx_fio_do_real(fio
, iparams
->cubic
.b0
);
2056 gmx_fio_do_real(fio
, iparams
->cubic
.kb
);
2057 gmx_fio_do_real(fio
, iparams
->cubic
.kcub
);
2061 case F_POLARIZATION
:
2062 gmx_fio_do_real(fio
, iparams
->polarize
.alpha
);
2065 gmx_fio_do_real(fio
, iparams
->anharm_polarize
.alpha
);
2066 gmx_fio_do_real(fio
, iparams
->anharm_polarize
.drcut
);
2067 gmx_fio_do_real(fio
, iparams
->anharm_polarize
.khyp
);
2070 if (file_version
< 31)
2072 gmx_fatal(FARGS
, "Old tpr files with water_polarization not supported. Make a new.");
2074 gmx_fio_do_real(fio
, iparams
->wpol
.al_x
);
2075 gmx_fio_do_real(fio
, iparams
->wpol
.al_y
);
2076 gmx_fio_do_real(fio
, iparams
->wpol
.al_z
);
2077 gmx_fio_do_real(fio
, iparams
->wpol
.rOH
);
2078 gmx_fio_do_real(fio
, iparams
->wpol
.rHH
);
2079 gmx_fio_do_real(fio
, iparams
->wpol
.rOD
);
2082 gmx_fio_do_real(fio
, iparams
->thole
.a
);
2083 gmx_fio_do_real(fio
, iparams
->thole
.alpha1
);
2084 gmx_fio_do_real(fio
, iparams
->thole
.alpha2
);
2085 gmx_fio_do_real(fio
, iparams
->thole
.rfac
);
2088 gmx_fio_do_real(fio
, iparams
->lj
.c6
);
2089 gmx_fio_do_real(fio
, iparams
->lj
.c12
);
2092 gmx_fio_do_real(fio
, iparams
->lj14
.c6A
);
2093 gmx_fio_do_real(fio
, iparams
->lj14
.c12A
);
2094 gmx_fio_do_real(fio
, iparams
->lj14
.c6B
);
2095 gmx_fio_do_real(fio
, iparams
->lj14
.c12B
);
2098 gmx_fio_do_real(fio
, iparams
->ljc14
.fqq
);
2099 gmx_fio_do_real(fio
, iparams
->ljc14
.qi
);
2100 gmx_fio_do_real(fio
, iparams
->ljc14
.qj
);
2101 gmx_fio_do_real(fio
, iparams
->ljc14
.c6
);
2102 gmx_fio_do_real(fio
, iparams
->ljc14
.c12
);
2104 case F_LJC_PAIRS_NB
:
2105 gmx_fio_do_real(fio
, iparams
->ljcnb
.qi
);
2106 gmx_fio_do_real(fio
, iparams
->ljcnb
.qj
);
2107 gmx_fio_do_real(fio
, iparams
->ljcnb
.c6
);
2108 gmx_fio_do_real(fio
, iparams
->ljcnb
.c12
);
2114 gmx_fio_do_real(fio
, iparams
->pdihs
.phiA
);
2115 gmx_fio_do_real(fio
, iparams
->pdihs
.cpA
);
2116 if ((ftype
== F_ANGRES
|| ftype
== F_ANGRESZ
) && file_version
< 42)
2118 /* Read the incorrectly stored multiplicity */
2119 gmx_fio_do_real(fio
, iparams
->harmonic
.rB
);
2120 gmx_fio_do_real(fio
, iparams
->harmonic
.krB
);
2121 iparams
->pdihs
.phiB
= iparams
->pdihs
.phiA
;
2122 iparams
->pdihs
.cpB
= iparams
->pdihs
.cpA
;
2126 gmx_fio_do_real(fio
, iparams
->pdihs
.phiB
);
2127 gmx_fio_do_real(fio
, iparams
->pdihs
.cpB
);
2128 gmx_fio_do_int(fio
, iparams
->pdihs
.mult
);
2132 gmx_fio_do_real(fio
, iparams
->pdihs
.phiA
);
2133 gmx_fio_do_real(fio
, iparams
->pdihs
.cpA
);
2136 gmx_fio_do_int(fio
, iparams
->disres
.label
);
2137 gmx_fio_do_int(fio
, iparams
->disres
.type
);
2138 gmx_fio_do_real(fio
, iparams
->disres
.low
);
2139 gmx_fio_do_real(fio
, iparams
->disres
.up1
);
2140 gmx_fio_do_real(fio
, iparams
->disres
.up2
);
2141 gmx_fio_do_real(fio
, iparams
->disres
.kfac
);
2144 gmx_fio_do_int(fio
, iparams
->orires
.ex
);
2145 gmx_fio_do_int(fio
, iparams
->orires
.label
);
2146 gmx_fio_do_int(fio
, iparams
->orires
.power
);
2147 gmx_fio_do_real(fio
, iparams
->orires
.c
);
2148 gmx_fio_do_real(fio
, iparams
->orires
.obs
);
2149 gmx_fio_do_real(fio
, iparams
->orires
.kfac
);
2152 if (file_version
< 82)
2154 gmx_fio_do_int(fio
, idum
);
2155 gmx_fio_do_int(fio
, idum
);
2157 gmx_fio_do_real(fio
, iparams
->dihres
.phiA
);
2158 gmx_fio_do_real(fio
, iparams
->dihres
.dphiA
);
2159 gmx_fio_do_real(fio
, iparams
->dihres
.kfacA
);
2160 if (file_version
>= 82)
2162 gmx_fio_do_real(fio
, iparams
->dihres
.phiB
);
2163 gmx_fio_do_real(fio
, iparams
->dihres
.dphiB
);
2164 gmx_fio_do_real(fio
, iparams
->dihres
.kfacB
);
2168 iparams
->dihres
.phiB
= iparams
->dihres
.phiA
;
2169 iparams
->dihres
.dphiB
= iparams
->dihres
.dphiA
;
2170 iparams
->dihres
.kfacB
= iparams
->dihres
.kfacA
;
2174 gmx_fio_do_rvec(fio
, iparams
->posres
.pos0A
);
2175 gmx_fio_do_rvec(fio
, iparams
->posres
.fcA
);
2176 if (bRead
&& file_version
< 27)
2178 copy_rvec(iparams
->posres
.pos0A
, iparams
->posres
.pos0B
);
2179 copy_rvec(iparams
->posres
.fcA
, iparams
->posres
.fcB
);
2183 gmx_fio_do_rvec(fio
, iparams
->posres
.pos0B
);
2184 gmx_fio_do_rvec(fio
, iparams
->posres
.fcB
);
2188 gmx_fio_do_int(fio
, iparams
->fbposres
.geom
);
2189 gmx_fio_do_rvec(fio
, iparams
->fbposres
.pos0
);
2190 gmx_fio_do_real(fio
, iparams
->fbposres
.r
);
2191 gmx_fio_do_real(fio
, iparams
->fbposres
.k
);
2194 gmx_fio_ndo_real(fio
, iparams
->cbtdihs
.cbtcA
, NR_CBTDIHS
);
2197 gmx_fio_ndo_real(fio
, iparams
->rbdihs
.rbcA
, NR_RBDIHS
);
2198 if (file_version
>= 25)
2200 gmx_fio_ndo_real(fio
, iparams
->rbdihs
.rbcB
, NR_RBDIHS
);
2204 /* Fourier dihedrals are internally represented
2205 * as Ryckaert-Bellemans since those are faster to compute.
2207 gmx_fio_ndo_real(fio
, iparams
->rbdihs
.rbcA
, NR_RBDIHS
);
2208 gmx_fio_ndo_real(fio
, iparams
->rbdihs
.rbcB
, NR_RBDIHS
);
2212 gmx_fio_do_real(fio
, iparams
->constr
.dA
);
2213 gmx_fio_do_real(fio
, iparams
->constr
.dB
);
2216 gmx_fio_do_real(fio
, iparams
->settle
.doh
);
2217 gmx_fio_do_real(fio
, iparams
->settle
.dhh
);
2220 gmx_fio_do_real(fio
, iparams
->vsite
.a
);
2225 gmx_fio_do_real(fio
, iparams
->vsite
.a
);
2226 gmx_fio_do_real(fio
, iparams
->vsite
.b
);
2231 gmx_fio_do_real(fio
, iparams
->vsite
.a
);
2232 gmx_fio_do_real(fio
, iparams
->vsite
.b
);
2233 gmx_fio_do_real(fio
, iparams
->vsite
.c
);
2236 gmx_fio_do_int(fio
, iparams
->vsiten
.n
);
2237 gmx_fio_do_real(fio
, iparams
->vsiten
.a
);
2242 /* We got rid of some parameters in version 68 */
2243 if (bRead
&& file_version
< 68)
2245 gmx_fio_do_real(fio
, rdum
);
2246 gmx_fio_do_real(fio
, rdum
);
2247 gmx_fio_do_real(fio
, rdum
);
2248 gmx_fio_do_real(fio
, rdum
);
2250 gmx_fio_do_real(fio
, iparams
->gb
.sar
);
2251 gmx_fio_do_real(fio
, iparams
->gb
.st
);
2252 gmx_fio_do_real(fio
, iparams
->gb
.pi
);
2253 gmx_fio_do_real(fio
, iparams
->gb
.gbr
);
2254 gmx_fio_do_real(fio
, iparams
->gb
.bmlt
);
2257 gmx_fio_do_int(fio
, iparams
->cmap
.cmapA
);
2258 gmx_fio_do_int(fio
, iparams
->cmap
.cmapB
);
2261 gmx_fatal(FARGS
, "unknown function type %d (%s) in %s line %d",
2262 ftype
, interaction_function
[ftype
].name
, __FILE__
, __LINE__
);
2266 static void do_ilist(t_fileio
*fio
, t_ilist
*ilist
, gmx_bool bRead
, int file_version
)
2270 if (file_version
< 44)
2272 for (i
= 0; i
< MAXNODES
; i
++)
2274 gmx_fio_do_int(fio
, idum
);
2277 gmx_fio_do_int(fio
, ilist
->nr
);
2280 snew(ilist
->iatoms
, ilist
->nr
);
2282 gmx_fio_ndo_int(fio
, ilist
->iatoms
, ilist
->nr
);
2285 static void do_ffparams(t_fileio
*fio
, gmx_ffparams_t
*ffparams
,
2286 gmx_bool bRead
, int file_version
)
2291 gmx_fio_do_int(fio
, ffparams
->atnr
);
2292 if (file_version
< 57)
2294 gmx_fio_do_int(fio
, idum
);
2296 gmx_fio_do_int(fio
, ffparams
->ntypes
);
2299 snew(ffparams
->functype
, ffparams
->ntypes
);
2300 snew(ffparams
->iparams
, ffparams
->ntypes
);
2302 /* Read/write all the function types */
2303 gmx_fio_ndo_int(fio
, ffparams
->functype
, ffparams
->ntypes
);
2305 if (file_version
>= 66)
2307 gmx_fio_do_double(fio
, ffparams
->reppow
);
2311 ffparams
->reppow
= 12.0;
2314 if (file_version
>= 57)
2316 gmx_fio_do_real(fio
, ffparams
->fudgeQQ
);
2319 /* Check whether all these function types are supported by the code.
2320 * In practice the code is backwards compatible, which means that the
2321 * numbering may have to be altered from old numbering to new numbering
2323 for (i
= 0; (i
< ffparams
->ntypes
); i
++)
2327 /* Loop over file versions */
2328 for (k
= 0; (k
< NFTUPD
); k
++)
2330 /* Compare the read file_version to the update table */
2331 if ((file_version
< ftupd
[k
].fvnr
) &&
2332 (ffparams
->functype
[i
] >= ftupd
[k
].ftype
))
2334 ffparams
->functype
[i
] += 1;
2339 do_iparams(fio
, ffparams
->functype
[i
], &ffparams
->iparams
[i
], bRead
,
2344 static void add_settle_atoms(t_ilist
*ilist
)
2348 /* Settle used to only store the first atom: add the other two */
2349 srenew(ilist
->iatoms
, 2*ilist
->nr
);
2350 for (i
= ilist
->nr
/2-1; i
>= 0; i
--)
2352 ilist
->iatoms
[4*i
+0] = ilist
->iatoms
[2*i
+0];
2353 ilist
->iatoms
[4*i
+1] = ilist
->iatoms
[2*i
+1];
2354 ilist
->iatoms
[4*i
+2] = ilist
->iatoms
[2*i
+1] + 1;
2355 ilist
->iatoms
[4*i
+3] = ilist
->iatoms
[2*i
+1] + 2;
2357 ilist
->nr
= 2*ilist
->nr
;
2360 static void do_ilists(t_fileio
*fio
, t_ilist
*ilist
, gmx_bool bRead
,
2367 for (j
= 0; (j
< F_NRE
); j
++)
2372 for (k
= 0; k
< NFTUPD
; k
++)
2374 if ((file_version
< ftupd
[k
].fvnr
) && (j
== ftupd
[k
].ftype
))
2383 ilist
[j
].iatoms
= NULL
;
2387 do_ilist(fio
, &ilist
[j
], bRead
, file_version
);
2388 if (file_version
< 78 && j
== F_SETTLE
&& ilist
[j
].nr
> 0)
2390 add_settle_atoms(&ilist
[j
]);
2396 static void do_idef(t_fileio
*fio
, gmx_ffparams_t
*ffparams
, gmx_moltype_t
*molt
,
2397 gmx_bool bRead
, int file_version
)
2399 do_ffparams(fio
, ffparams
, bRead
, file_version
);
2401 if (file_version
>= 54)
2403 gmx_fio_do_real(fio
, ffparams
->fudgeQQ
);
2406 do_ilists(fio
, molt
->ilist
, bRead
, file_version
);
2409 static void do_block(t_fileio
*fio
, t_block
*block
, gmx_bool bRead
, int file_version
)
2411 int i
, idum
, dum_nra
, *dum_a
;
2413 if (file_version
< 44)
2415 for (i
= 0; i
< MAXNODES
; i
++)
2417 gmx_fio_do_int(fio
, idum
);
2420 gmx_fio_do_int(fio
, block
->nr
);
2421 if (file_version
< 51)
2423 gmx_fio_do_int(fio
, dum_nra
);
2427 if ((block
->nalloc_index
> 0) && (NULL
!= block
->index
))
2429 sfree(block
->index
);
2431 block
->nalloc_index
= block
->nr
+1;
2432 snew(block
->index
, block
->nalloc_index
);
2434 gmx_fio_ndo_int(fio
, block
->index
, block
->nr
+1);
2436 if (file_version
< 51 && dum_nra
> 0)
2438 snew(dum_a
, dum_nra
);
2439 gmx_fio_ndo_int(fio
, dum_a
, dum_nra
);
2444 static void do_blocka(t_fileio
*fio
, t_blocka
*block
, gmx_bool bRead
,
2449 if (file_version
< 44)
2451 for (i
= 0; i
< MAXNODES
; i
++)
2453 gmx_fio_do_int(fio
, idum
);
2456 gmx_fio_do_int(fio
, block
->nr
);
2457 gmx_fio_do_int(fio
, block
->nra
);
2460 block
->nalloc_index
= block
->nr
+1;
2461 snew(block
->index
, block
->nalloc_index
);
2462 block
->nalloc_a
= block
->nra
;
2463 snew(block
->a
, block
->nalloc_a
);
2465 gmx_fio_ndo_int(fio
, block
->index
, block
->nr
+1);
2466 gmx_fio_ndo_int(fio
, block
->a
, block
->nra
);
2469 /* This is a primitive routine to make it possible to translate atomic numbers
2470 * to element names when reading TPR files, without making the Gromacs library
2471 * directory a dependency on mdrun (which is the case if we need elements.dat).
2474 atomicnumber_to_element(int atomicnumber
)
2478 /* This does not have to be complete, so we only include elements likely
2479 * to occur in PDB files.
2481 switch (atomicnumber
)
2483 case 1: p
= "H"; break;
2484 case 5: p
= "B"; break;
2485 case 6: p
= "C"; break;
2486 case 7: p
= "N"; break;
2487 case 8: p
= "O"; break;
2488 case 9: p
= "F"; break;
2489 case 11: p
= "Na"; break;
2490 case 12: p
= "Mg"; break;
2491 case 15: p
= "P"; break;
2492 case 16: p
= "S"; break;
2493 case 17: p
= "Cl"; break;
2494 case 18: p
= "Ar"; break;
2495 case 19: p
= "K"; break;
2496 case 20: p
= "Ca"; break;
2497 case 25: p
= "Mn"; break;
2498 case 26: p
= "Fe"; break;
2499 case 28: p
= "Ni"; break;
2500 case 29: p
= "Cu"; break;
2501 case 30: p
= "Zn"; break;
2502 case 35: p
= "Br"; break;
2503 case 47: p
= "Ag"; break;
2504 default: p
= ""; break;
2510 static void do_atom(t_fileio
*fio
, t_atom
*atom
, int ngrp
, gmx_bool bRead
,
2511 int file_version
, gmx_groups_t
*groups
, int atnr
)
2515 gmx_fio_do_real(fio
, atom
->m
);
2516 gmx_fio_do_real(fio
, atom
->q
);
2517 gmx_fio_do_real(fio
, atom
->mB
);
2518 gmx_fio_do_real(fio
, atom
->qB
);
2519 gmx_fio_do_ushort(fio
, atom
->type
);
2520 gmx_fio_do_ushort(fio
, atom
->typeB
);
2521 gmx_fio_do_int(fio
, atom
->ptype
);
2522 gmx_fio_do_int(fio
, atom
->resind
);
2523 if (file_version
>= 52)
2525 gmx_fio_do_int(fio
, atom
->atomnumber
);
2528 /* Set element string from atomic number if present.
2529 * This routine returns an empty string if the name is not found.
2531 std::strncpy(atom
->elem
, atomicnumber_to_element(atom
->atomnumber
), 4);
2532 /* avoid warnings about potentially unterminated string */
2533 atom
->elem
[3] = '\0';
2538 atom
->atomnumber
= 0;
2540 if (file_version
< 23)
2544 else if (file_version
< 39)
2553 if (file_version
< 57)
2555 unsigned char uchar
[egcNR
];
2556 gmx_fio_ndo_uchar(fio
, uchar
, myngrp
);
2557 for (i
= myngrp
; (i
< ngrp
); i
++)
2561 /* Copy the old data format to the groups struct */
2562 for (i
= 0; i
< ngrp
; i
++)
2564 groups
->grpnr
[i
][atnr
] = uchar
[i
];
2569 static void do_grps(t_fileio
*fio
, int ngrp
, t_grps grps
[], gmx_bool bRead
,
2574 if (file_version
< 23)
2578 else if (file_version
< 39)
2587 for (j
= 0; (j
< ngrp
); j
++)
2591 gmx_fio_do_int(fio
, grps
[j
].nr
);
2594 snew(grps
[j
].nm_ind
, grps
[j
].nr
);
2596 gmx_fio_ndo_int(fio
, grps
[j
].nm_ind
, grps
[j
].nr
);
2601 snew(grps
[j
].nm_ind
, grps
[j
].nr
);
2606 static void do_symstr(t_fileio
*fio
, char ***nm
, gmx_bool bRead
, t_symtab
*symtab
)
2612 gmx_fio_do_int(fio
, ls
);
2613 *nm
= get_symtab_handle(symtab
, ls
);
2617 ls
= lookup_symtab(symtab
, *nm
);
2618 gmx_fio_do_int(fio
, ls
);
2622 static void do_strstr(t_fileio
*fio
, int nstr
, char ***nm
, gmx_bool bRead
,
2627 for (j
= 0; (j
< nstr
); j
++)
2629 do_symstr(fio
, &(nm
[j
]), bRead
, symtab
);
2633 static void do_resinfo(t_fileio
*fio
, int n
, t_resinfo
*ri
, gmx_bool bRead
,
2634 t_symtab
*symtab
, int file_version
)
2638 for (j
= 0; (j
< n
); j
++)
2640 do_symstr(fio
, &(ri
[j
].name
), bRead
, symtab
);
2641 if (file_version
>= 63)
2643 gmx_fio_do_int(fio
, ri
[j
].nr
);
2644 gmx_fio_do_uchar(fio
, ri
[j
].ic
);
2654 static void do_atoms(t_fileio
*fio
, t_atoms
*atoms
, gmx_bool bRead
, t_symtab
*symtab
,
2656 gmx_groups_t
*groups
)
2660 gmx_fio_do_int(fio
, atoms
->nr
);
2661 gmx_fio_do_int(fio
, atoms
->nres
);
2662 if (file_version
< 57)
2664 gmx_fio_do_int(fio
, groups
->ngrpname
);
2665 for (i
= 0; i
< egcNR
; i
++)
2667 groups
->ngrpnr
[i
] = atoms
->nr
;
2668 snew(groups
->grpnr
[i
], groups
->ngrpnr
[i
]);
2673 snew(atoms
->atom
, atoms
->nr
);
2674 snew(atoms
->atomname
, atoms
->nr
);
2675 snew(atoms
->atomtype
, atoms
->nr
);
2676 snew(atoms
->atomtypeB
, atoms
->nr
);
2677 snew(atoms
->resinfo
, atoms
->nres
);
2678 if (file_version
< 57)
2680 snew(groups
->grpname
, groups
->ngrpname
);
2682 atoms
->pdbinfo
= NULL
;
2684 for (i
= 0; (i
< atoms
->nr
); i
++)
2686 do_atom(fio
, &atoms
->atom
[i
], egcNR
, bRead
, file_version
, groups
, i
);
2688 do_strstr(fio
, atoms
->nr
, atoms
->atomname
, bRead
, symtab
);
2689 if (bRead
&& (file_version
<= 20))
2691 for (i
= 0; i
< atoms
->nr
; i
++)
2693 atoms
->atomtype
[i
] = put_symtab(symtab
, "?");
2694 atoms
->atomtypeB
[i
] = put_symtab(symtab
, "?");
2699 do_strstr(fio
, atoms
->nr
, atoms
->atomtype
, bRead
, symtab
);
2700 do_strstr(fio
, atoms
->nr
, atoms
->atomtypeB
, bRead
, symtab
);
2702 do_resinfo(fio
, atoms
->nres
, atoms
->resinfo
, bRead
, symtab
, file_version
);
2704 if (file_version
< 57)
2706 do_strstr(fio
, groups
->ngrpname
, groups
->grpname
, bRead
, symtab
);
2708 do_grps(fio
, egcNR
, groups
->grps
, bRead
, file_version
);
2712 static void do_groups(t_fileio
*fio
, gmx_groups_t
*groups
,
2713 gmx_bool bRead
, t_symtab
*symtab
,
2718 do_grps(fio
, egcNR
, groups
->grps
, bRead
, file_version
);
2719 gmx_fio_do_int(fio
, groups
->ngrpname
);
2722 snew(groups
->grpname
, groups
->ngrpname
);
2724 do_strstr(fio
, groups
->ngrpname
, groups
->grpname
, bRead
, symtab
);
2725 for (g
= 0; g
< egcNR
; g
++)
2727 gmx_fio_do_int(fio
, groups
->ngrpnr
[g
]);
2728 if (groups
->ngrpnr
[g
] == 0)
2732 groups
->grpnr
[g
] = NULL
;
2739 snew(groups
->grpnr
[g
], groups
->ngrpnr
[g
]);
2741 gmx_fio_ndo_uchar(fio
, groups
->grpnr
[g
], groups
->ngrpnr
[g
]);
2746 static void do_atomtypes(t_fileio
*fio
, t_atomtypes
*atomtypes
, gmx_bool bRead
,
2751 if (file_version
> 25)
2753 gmx_fio_do_int(fio
, atomtypes
->nr
);
2757 snew(atomtypes
->radius
, j
);
2758 snew(atomtypes
->vol
, j
);
2759 snew(atomtypes
->surftens
, j
);
2760 snew(atomtypes
->atomnumber
, j
);
2761 snew(atomtypes
->gb_radius
, j
);
2762 snew(atomtypes
->S_hct
, j
);
2764 gmx_fio_ndo_real(fio
, atomtypes
->radius
, j
);
2765 gmx_fio_ndo_real(fio
, atomtypes
->vol
, j
);
2766 gmx_fio_ndo_real(fio
, atomtypes
->surftens
, j
);
2767 if (file_version
>= 40)
2769 gmx_fio_ndo_int(fio
, atomtypes
->atomnumber
, j
);
2771 if (file_version
>= 60)
2773 gmx_fio_ndo_real(fio
, atomtypes
->gb_radius
, j
);
2774 gmx_fio_ndo_real(fio
, atomtypes
->S_hct
, j
);
2779 /* File versions prior to 26 cannot do GBSA,
2780 * so they dont use this structure
2783 atomtypes
->radius
= NULL
;
2784 atomtypes
->vol
= NULL
;
2785 atomtypes
->surftens
= NULL
;
2786 atomtypes
->atomnumber
= NULL
;
2787 atomtypes
->gb_radius
= NULL
;
2788 atomtypes
->S_hct
= NULL
;
2792 static void do_symtab(t_fileio
*fio
, t_symtab
*symtab
, gmx_bool bRead
)
2798 gmx_fio_do_int(fio
, symtab
->nr
);
2802 snew(symtab
->symbuf
, 1);
2803 symbuf
= symtab
->symbuf
;
2804 symbuf
->bufsize
= nr
;
2805 snew(symbuf
->buf
, nr
);
2806 for (i
= 0; (i
< nr
); i
++)
2808 gmx_fio_do_string(fio
, buf
);
2809 symbuf
->buf
[i
] = gmx_strdup(buf
);
2814 symbuf
= symtab
->symbuf
;
2815 while (symbuf
!= NULL
)
2817 for (i
= 0; (i
< symbuf
->bufsize
) && (i
< nr
); i
++)
2819 gmx_fio_do_string(fio
, symbuf
->buf
[i
]);
2822 symbuf
= symbuf
->next
;
2826 gmx_fatal(FARGS
, "nr of symtab strings left: %d", nr
);
2831 static void do_cmap(t_fileio
*fio
, gmx_cmap_t
*cmap_grid
, gmx_bool bRead
)
2833 int i
, j
, ngrid
, gs
, nelem
;
2835 gmx_fio_do_int(fio
, cmap_grid
->ngrid
);
2836 gmx_fio_do_int(fio
, cmap_grid
->grid_spacing
);
2838 ngrid
= cmap_grid
->ngrid
;
2839 gs
= cmap_grid
->grid_spacing
;
2844 snew(cmap_grid
->cmapdata
, ngrid
);
2846 for (i
= 0; i
< cmap_grid
->ngrid
; i
++)
2848 snew(cmap_grid
->cmapdata
[i
].cmap
, 4*nelem
);
2852 for (i
= 0; i
< cmap_grid
->ngrid
; i
++)
2854 for (j
= 0; j
< nelem
; j
++)
2856 gmx_fio_do_real(fio
, cmap_grid
->cmapdata
[i
].cmap
[j
*4]);
2857 gmx_fio_do_real(fio
, cmap_grid
->cmapdata
[i
].cmap
[j
*4+1]);
2858 gmx_fio_do_real(fio
, cmap_grid
->cmapdata
[i
].cmap
[j
*4+2]);
2859 gmx_fio_do_real(fio
, cmap_grid
->cmapdata
[i
].cmap
[j
*4+3]);
2865 static void do_moltype(t_fileio
*fio
, gmx_moltype_t
*molt
, gmx_bool bRead
,
2866 t_symtab
*symtab
, int file_version
,
2867 gmx_groups_t
*groups
)
2869 if (file_version
>= 57)
2871 do_symstr(fio
, &(molt
->name
), bRead
, symtab
);
2874 do_atoms(fio
, &molt
->atoms
, bRead
, symtab
, file_version
, groups
);
2876 if (file_version
>= 57)
2878 do_ilists(fio
, molt
->ilist
, bRead
, file_version
);
2880 do_block(fio
, &molt
->cgs
, bRead
, file_version
);
2883 /* This used to be in the atoms struct */
2884 do_blocka(fio
, &molt
->excls
, bRead
, file_version
);
2887 static void do_molblock(t_fileio
*fio
, gmx_molblock_t
*molb
, gmx_bool bRead
)
2889 gmx_fio_do_int(fio
, molb
->type
);
2890 gmx_fio_do_int(fio
, molb
->nmol
);
2891 gmx_fio_do_int(fio
, molb
->natoms_mol
);
2892 /* Position restraint coordinates */
2893 gmx_fio_do_int(fio
, molb
->nposres_xA
);
2894 if (molb
->nposres_xA
> 0)
2898 snew(molb
->posres_xA
, molb
->nposres_xA
);
2900 gmx_fio_ndo_rvec(fio
, molb
->posres_xA
, molb
->nposres_xA
);
2902 gmx_fio_do_int(fio
, molb
->nposres_xB
);
2903 if (molb
->nposres_xB
> 0)
2907 snew(molb
->posres_xB
, molb
->nposres_xB
);
2909 gmx_fio_ndo_rvec(fio
, molb
->posres_xB
, molb
->nposres_xB
);
2914 static t_block
mtop_mols(gmx_mtop_t
*mtop
)
2920 for (mb
= 0; mb
< mtop
->nmolblock
; mb
++)
2922 mols
.nr
+= mtop
->molblock
[mb
].nmol
;
2924 mols
.nalloc_index
= mols
.nr
+ 1;
2925 snew(mols
.index
, mols
.nalloc_index
);
2930 for (mb
= 0; mb
< mtop
->nmolblock
; mb
++)
2932 for (mol
= 0; mol
< mtop
->molblock
[mb
].nmol
; mol
++)
2934 a
+= mtop
->molblock
[mb
].natoms_mol
;
2943 static void add_posres_molblock(gmx_mtop_t
*mtop
)
2948 gmx_molblock_t
*molb
;
2951 /* posres reference positions are stored in ip->posres (if present) and
2952 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2953 posres.pos0A are identical to fbposres.pos0. */
2954 il
= &mtop
->moltype
[0].ilist
[F_POSRES
];
2955 ilfb
= &mtop
->moltype
[0].ilist
[F_FBPOSRES
];
2956 if (il
->nr
== 0 && ilfb
->nr
== 0)
2962 for (i
= 0; i
< il
->nr
; i
+= 2)
2964 ip
= &mtop
->ffparams
.iparams
[il
->iatoms
[i
]];
2965 am
= std::max(am
, il
->iatoms
[i
+1]);
2966 if (ip
->posres
.pos0B
[XX
] != ip
->posres
.pos0A
[XX
] ||
2967 ip
->posres
.pos0B
[YY
] != ip
->posres
.pos0A
[YY
] ||
2968 ip
->posres
.pos0B
[ZZ
] != ip
->posres
.pos0A
[ZZ
])
2973 /* This loop is required if we have only flat-bottomed posres:
2975 - bFE == FALSE (no B-state for flat-bottomed posres) */
2978 for (i
= 0; i
< ilfb
->nr
; i
+= 2)
2980 am
= std::max(am
, ilfb
->iatoms
[i
+1]);
2983 /* Make the posres coordinate block end at a molecule end */
2985 while (am
>= mtop
->mols
.index
[mol
+1])
2989 molb
= &mtop
->molblock
[0];
2990 molb
->nposres_xA
= mtop
->mols
.index
[mol
+1];
2991 snew(molb
->posres_xA
, molb
->nposres_xA
);
2994 molb
->nposres_xB
= molb
->nposres_xA
;
2995 snew(molb
->posres_xB
, molb
->nposres_xB
);
2999 molb
->nposres_xB
= 0;
3001 for (i
= 0; i
< il
->nr
; i
+= 2)
3003 ip
= &mtop
->ffparams
.iparams
[il
->iatoms
[i
]];
3004 a
= il
->iatoms
[i
+1];
3005 molb
->posres_xA
[a
][XX
] = ip
->posres
.pos0A
[XX
];
3006 molb
->posres_xA
[a
][YY
] = ip
->posres
.pos0A
[YY
];
3007 molb
->posres_xA
[a
][ZZ
] = ip
->posres
.pos0A
[ZZ
];
3010 molb
->posres_xB
[a
][XX
] = ip
->posres
.pos0B
[XX
];
3011 molb
->posres_xB
[a
][YY
] = ip
->posres
.pos0B
[YY
];
3012 molb
->posres_xB
[a
][ZZ
] = ip
->posres
.pos0B
[ZZ
];
3017 /* If only flat-bottomed posres are present, take reference pos from them.
3018 Here: bFE == FALSE */
3019 for (i
= 0; i
< ilfb
->nr
; i
+= 2)
3021 ip
= &mtop
->ffparams
.iparams
[ilfb
->iatoms
[i
]];
3022 a
= ilfb
->iatoms
[i
+1];
3023 molb
->posres_xA
[a
][XX
] = ip
->fbposres
.pos0
[XX
];
3024 molb
->posres_xA
[a
][YY
] = ip
->fbposres
.pos0
[YY
];
3025 molb
->posres_xA
[a
][ZZ
] = ip
->fbposres
.pos0
[ZZ
];
3030 static void set_disres_npair(gmx_mtop_t
*mtop
)
3033 gmx_mtop_ilistloop_t iloop
;
3034 t_ilist
*ilist
, *il
;
3038 ip
= mtop
->ffparams
.iparams
;
3040 iloop
= gmx_mtop_ilistloop_init(mtop
);
3041 while (gmx_mtop_ilistloop_next(iloop
, &ilist
, &nmol
))
3043 il
= &ilist
[F_DISRES
];
3049 for (i
= 0; i
< il
->nr
; i
+= 3)
3052 if (i
+3 == il
->nr
|| ip
[a
[i
]].disres
.label
!= ip
[a
[i
+3]].disres
.label
)
3054 ip
[a
[i
]].disres
.npair
= npair
;
3062 static void do_mtop(t_fileio
*fio
, gmx_mtop_t
*mtop
, gmx_bool bRead
,
3072 do_symtab(fio
, &(mtop
->symtab
), bRead
);
3074 do_symstr(fio
, &(mtop
->name
), bRead
, &(mtop
->symtab
));
3076 if (file_version
>= 57)
3078 do_ffparams(fio
, &mtop
->ffparams
, bRead
, file_version
);
3080 gmx_fio_do_int(fio
, mtop
->nmoltype
);
3088 snew(mtop
->moltype
, mtop
->nmoltype
);
3089 if (file_version
< 57)
3091 mtop
->moltype
[0].name
= mtop
->name
;
3094 for (mt
= 0; mt
< mtop
->nmoltype
; mt
++)
3096 do_moltype(fio
, &mtop
->moltype
[mt
], bRead
, &mtop
->symtab
, file_version
,
3100 if (file_version
>= 57)
3102 gmx_fio_do_int(fio
, mtop
->nmolblock
);
3106 mtop
->nmolblock
= 1;
3110 snew(mtop
->molblock
, mtop
->nmolblock
);
3112 if (file_version
>= 57)
3114 for (mb
= 0; mb
< mtop
->nmolblock
; mb
++)
3116 do_molblock(fio
, &mtop
->molblock
[mb
], bRead
);
3118 gmx_fio_do_int(fio
, mtop
->natoms
);
3122 mtop
->molblock
[0].type
= 0;
3123 mtop
->molblock
[0].nmol
= 1;
3124 mtop
->molblock
[0].natoms_mol
= mtop
->moltype
[0].atoms
.nr
;
3125 mtop
->molblock
[0].nposres_xA
= 0;
3126 mtop
->molblock
[0].nposres_xB
= 0;
3129 if (file_version
>= tpxv_IntermolecularBondeds
)
3131 gmx_fio_do_gmx_bool(fio
, mtop
->bIntermolecularInteractions
);
3132 if (mtop
->bIntermolecularInteractions
)
3136 snew(mtop
->intermolecular_ilist
, F_NRE
);
3138 do_ilists(fio
, mtop
->intermolecular_ilist
, bRead
, file_version
);
3143 mtop
->bIntermolecularInteractions
= FALSE
;
3146 do_atomtypes (fio
, &(mtop
->atomtypes
), bRead
, file_version
);
3148 if (file_version
< 57)
3150 do_idef (fio
, &mtop
->ffparams
, &mtop
->moltype
[0], bRead
, file_version
);
3151 mtop
->natoms
= mtop
->moltype
[0].atoms
.nr
;
3154 if (file_version
>= 65)
3156 do_cmap(fio
, &mtop
->ffparams
.cmap_grid
, bRead
);
3160 mtop
->ffparams
.cmap_grid
.ngrid
= 0;
3161 mtop
->ffparams
.cmap_grid
.grid_spacing
= 0;
3162 mtop
->ffparams
.cmap_grid
.cmapdata
= NULL
;
3165 if (file_version
>= 57)
3167 do_groups(fio
, &mtop
->groups
, bRead
, &(mtop
->symtab
), file_version
);
3170 if (file_version
< 57)
3172 do_block(fio
, &mtop
->moltype
[0].cgs
, bRead
, file_version
);
3173 do_block(fio
, &mtop
->mols
, bRead
, file_version
);
3174 /* Add the posres coordinates to the molblock */
3175 add_posres_molblock(mtop
);
3179 if (file_version
>= 57)
3181 done_block(&mtop
->mols
);
3182 mtop
->mols
= mtop_mols(mtop
);
3186 if (file_version
< 51)
3188 /* Here used to be the shake blocks */
3189 do_blocka(fio
, &dumb
, bRead
, file_version
);
3202 close_symtab(&(mtop
->symtab
));
3206 /* If TopOnlyOK is TRUE then we can read even future versions
3207 * of tpx files, provided the file_generation hasn't changed.
3208 * If it is FALSE, we need the inputrecord too, and bail out
3209 * if the file is newer than the program.
3211 * The version and generation if the topology (see top of this file)
3212 * are returned in the two last arguments.
3214 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
3216 static void do_tpxheader(t_fileio
*fio
, gmx_bool bRead
, t_tpxheader
*tpx
,
3217 gmx_bool TopOnlyOK
, int *file_version
,
3218 int *file_generation
)
3221 char file_tag
[STRLEN
];
3228 /* XDR binary topology file */
3229 precision
= sizeof(real
);
3232 gmx_fio_do_string(fio
, buf
);
3233 if (std::strncmp(buf
, "VERSION", 7))
3235 gmx_fatal(FARGS
, "Can not read file %s,\n"
3236 " this file is from a GROMACS version which is older than 2.0\n"
3237 " Make a new one with grompp or use a gro or pdb file, if possible",
3238 gmx_fio_getname(fio
));
3240 gmx_fio_do_int(fio
, precision
);
3241 bDouble
= (precision
== sizeof(double));
3242 if ((precision
!= sizeof(float)) && !bDouble
)
3244 gmx_fatal(FARGS
, "Unknown precision in file %s: real is %d bytes "
3245 "instead of %d or %d",
3246 gmx_fio_getname(fio
), precision
, sizeof(float), sizeof(double));
3248 gmx_fio_setprecision(fio
, bDouble
);
3249 fprintf(stderr
, "Reading file %s, %s (%s precision)\n",
3250 gmx_fio_getname(fio
), buf
, bDouble
? "double" : "single");
3254 gmx_fio_write_string(fio
, gmx_version());
3255 bDouble
= (precision
== sizeof(double));
3256 gmx_fio_setprecision(fio
, bDouble
);
3257 gmx_fio_do_int(fio
, precision
);
3259 sprintf(file_tag
, "%s", tpx_tag
);
3260 fgen
= tpx_generation
;
3263 /* Check versions! */
3264 gmx_fio_do_int(fio
, fver
);
3266 /* This is for backward compatibility with development versions 77-79
3267 * where the tag was, mistakenly, placed before the generation,
3268 * which would cause a segv instead of a proper error message
3269 * when reading the topology only from tpx with <77 code.
3271 if (fver
>= 77 && fver
<= 79)
3273 gmx_fio_do_string(fio
, file_tag
);
3278 gmx_fio_do_int(fio
, fgen
);
3287 gmx_fio_do_string(fio
, file_tag
);
3293 /* Versions before 77 don't have the tag, set it to release */
3294 sprintf(file_tag
, "%s", TPX_TAG_RELEASE
);
3297 if (std::strcmp(file_tag
, tpx_tag
) != 0)
3299 fprintf(stderr
, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3302 /* We only support reading tpx files with the same tag as the code
3303 * or tpx files with the release tag and with lower version number.
3305 if (std::strcmp(file_tag
, TPX_TAG_RELEASE
) != 0 && fver
< tpx_version
)
3307 gmx_fatal(FARGS
, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3308 gmx_fio_getname(fio
), fver
, file_tag
,
3309 tpx_version
, tpx_tag
);
3314 if (file_version
!= NULL
)
3316 *file_version
= fver
;
3318 if (file_generation
!= NULL
)
3320 *file_generation
= fgen
;
3324 if ((fver
<= tpx_incompatible_version
) ||
3325 ((fver
> tpx_version
) && !TopOnlyOK
) ||
3326 (fgen
> tpx_generation
) ||
3327 tpx_version
== 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3329 gmx_fatal(FARGS
, "reading tpx file (%s) version %d with version %d program",
3330 gmx_fio_getname(fio
), fver
, tpx_version
);
3333 gmx_fio_do_int(fio
, tpx
->natoms
);
3336 gmx_fio_do_int(fio
, tpx
->ngtc
);
3344 gmx_fio_do_int(fio
, idum
);
3345 gmx_fio_do_real(fio
, rdum
);
3349 gmx_fio_do_int(fio
, tpx
->fep_state
);
3351 gmx_fio_do_real(fio
, tpx
->lambda
);
3352 gmx_fio_do_int(fio
, tpx
->bIr
);
3353 gmx_fio_do_int(fio
, tpx
->bTop
);
3354 gmx_fio_do_int(fio
, tpx
->bX
);
3355 gmx_fio_do_int(fio
, tpx
->bV
);
3356 gmx_fio_do_int(fio
, tpx
->bF
);
3357 gmx_fio_do_int(fio
, tpx
->bBox
);
3359 if ((fgen
> tpx_generation
))
3361 /* This can only happen if TopOnlyOK=TRUE */
3366 static int do_tpx(t_fileio
*fio
, gmx_bool bRead
,
3367 t_inputrec
*ir
, t_state
*state
, gmx_mtop_t
*mtop
,
3368 gmx_bool bXVallocated
)
3374 int file_version
, file_generation
;
3377 gmx_bool bPeriodicMols
;
3381 tpx
.natoms
= state
->natoms
;
3382 tpx
.ngtc
= state
->ngtc
;
3383 tpx
.fep_state
= state
->fep_state
;
3384 tpx
.lambda
= state
->lambda
[efptFEP
];
3385 tpx
.bIr
= (ir
!= NULL
);
3386 tpx
.bTop
= (mtop
!= NULL
);
3387 tpx
.bX
= (state
->x
!= NULL
);
3388 tpx
.bV
= (state
->v
!= NULL
);
3393 TopOnlyOK
= (ir
== NULL
);
3395 do_tpxheader(fio
, bRead
, &tpx
, TopOnlyOK
, &file_version
, &file_generation
);
3404 init_state(state
, 0, tpx
.ngtc
, 0, 0, 0);
3405 state
->natoms
= tpx
.natoms
;
3406 state
->nalloc
= tpx
.natoms
;
3412 init_state(state
, tpx
.natoms
, tpx
.ngtc
, 0, 0, 0);
3416 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3418 do_test(fio
, tpx
.bBox
, state
->box
);
3421 gmx_fio_ndo_rvec(fio
, state
->box
, DIM
);
3422 if (file_version
>= 51)
3424 gmx_fio_ndo_rvec(fio
, state
->box_rel
, DIM
);
3428 /* We initialize box_rel after reading the inputrec */
3429 clear_mat(state
->box_rel
);
3431 if (file_version
>= 28)
3433 gmx_fio_ndo_rvec(fio
, state
->boxv
, DIM
);
3434 if (file_version
< 56)
3437 gmx_fio_ndo_rvec(fio
, mdum
, DIM
);
3442 if (state
->ngtc
> 0 && file_version
>= 28)
3445 snew(dumv
, state
->ngtc
);
3446 if (file_version
< 69)
3448 gmx_fio_ndo_real(fio
, dumv
, state
->ngtc
);
3450 /* These used to be the Berendsen tcoupl_lambda's */
3451 gmx_fio_ndo_real(fio
, dumv
, state
->ngtc
);
3455 /* Prior to tpx version 26, the inputrec was here.
3456 * I moved it to enable partial forward-compatibility
3457 * for analysis/viewer programs.
3459 if (file_version
< 26)
3461 do_test(fio
, tpx
.bIr
, ir
);
3466 do_inputrec(fio
, ir
, bRead
, file_version
,
3467 mtop
? &mtop
->ffparams
.fudgeQQ
: NULL
);
3471 do_inputrec(fio
, &dum_ir
, bRead
, file_version
,
3472 mtop
? &mtop
->ffparams
.fudgeQQ
: NULL
);
3473 done_inputrec(&dum_ir
);
3479 do_test(fio
, tpx
.bTop
, mtop
);
3484 do_mtop(fio
, mtop
, bRead
, file_version
);
3488 do_mtop(fio
, &dum_top
, bRead
, file_version
);
3489 done_mtop(&dum_top
, TRUE
);
3492 do_test(fio
, tpx
.bX
, state
->x
);
3497 state
->flags
|= (1<<estX
);
3499 gmx_fio_ndo_rvec(fio
, state
->x
, state
->natoms
);
3502 do_test(fio
, tpx
.bV
, state
->v
);
3507 state
->flags
|= (1<<estV
);
3509 gmx_fio_ndo_rvec(fio
, state
->v
, state
->natoms
);
3512 // No need to run do_test when the last argument is NULL
3516 snew(dummyForces
, state
->natoms
);
3517 gmx_fio_ndo_rvec(fio
, dummyForces
, state
->natoms
);
3521 /* Starting with tpx version 26, we have the inputrec
3522 * at the end of the file, so we can ignore it
3523 * if the file is never than the software (but still the
3524 * same generation - see comments at the top of this file.
3529 bPeriodicMols
= FALSE
;
3530 if (file_version
>= 26)
3532 do_test(fio
, tpx
.bIr
, ir
);
3535 if (file_version
>= 53)
3537 /* Removed the pbc info from do_inputrec, since we always want it */
3541 bPeriodicMols
= ir
->bPeriodicMols
;
3543 gmx_fio_do_int(fio
, ePBC
);
3544 gmx_fio_do_gmx_bool(fio
, bPeriodicMols
);
3546 if (file_generation
<= tpx_generation
&& ir
)
3548 do_inputrec(fio
, ir
, bRead
, file_version
, mtop
? &mtop
->ffparams
.fudgeQQ
: NULL
);
3549 if (file_version
< 51)
3551 set_box_rel(ir
, state
);
3553 if (file_version
< 53)
3556 bPeriodicMols
= ir
->bPeriodicMols
;
3559 if (bRead
&& ir
&& file_version
>= 53)
3561 /* We need to do this after do_inputrec, since that initializes ir */
3563 ir
->bPeriodicMols
= bPeriodicMols
;
3572 if (state
->ngtc
== 0)
3574 /* Reading old version without tcoupl state data: set it */
3575 init_gtc_state(state
, ir
->opts
.ngtc
, 0, ir
->opts
.nhchainlength
);
3577 if (tpx
.bTop
&& mtop
)
3579 if (file_version
< 57)
3581 if (mtop
->moltype
[0].ilist
[F_DISRES
].nr
> 0)
3583 ir
->eDisre
= edrSimple
;
3587 ir
->eDisre
= edrNone
;
3590 set_disres_npair(mtop
);
3594 if (tpx
.bTop
&& mtop
)
3596 gmx_mtop_finalize(mtop
);
3603 static t_fileio
*open_tpx(const char *fn
, const char *mode
)
3605 return gmx_fio_open(fn
, mode
);
3608 static void close_tpx(t_fileio
*fio
)
3613 /************************************************************
3615 * The following routines are the exported ones
3617 ************************************************************/
3619 void read_tpxheader(const char *fn
, t_tpxheader
*tpx
, gmx_bool TopOnlyOK
,
3620 int *file_version
, int *file_generation
)
3624 fio
= open_tpx(fn
, "r");
3625 do_tpxheader(fio
, TRUE
, tpx
, TopOnlyOK
, file_version
, file_generation
);
3629 void write_tpx_state(const char *fn
,
3630 t_inputrec
*ir
, t_state
*state
, gmx_mtop_t
*mtop
)
3634 fio
= open_tpx(fn
, "w");
3635 do_tpx(fio
, FALSE
, ir
, state
, mtop
, FALSE
);
3639 void read_tpx_state(const char *fn
,
3640 t_inputrec
*ir
, t_state
*state
, gmx_mtop_t
*mtop
)
3644 fio
= open_tpx(fn
, "r");
3645 do_tpx(fio
, TRUE
, ir
, state
, mtop
, FALSE
);
3649 int read_tpx(const char *fn
,
3650 t_inputrec
*ir
, matrix box
, int *natoms
,
3651 rvec
*x
, rvec
*v
, gmx_mtop_t
*mtop
)
3659 fio
= open_tpx(fn
, "r");
3660 ePBC
= do_tpx(fio
, TRUE
, ir
, &state
, mtop
, TRUE
);
3662 *natoms
= state
.natoms
;
3665 copy_mat(state
.box
, box
);
3674 int read_tpx_top(const char *fn
,
3675 t_inputrec
*ir
, matrix box
, int *natoms
,
3676 rvec
*x
, rvec
*v
, t_topology
*top
)
3681 ePBC
= read_tpx(fn
, ir
, box
, natoms
, x
, v
, &mtop
);
3683 *top
= gmx_mtop_t_to_t_topology(&mtop
);
3688 gmx_bool
fn2bTPX(const char *file
)
3690 return (efTPR
== fn2ftp(file
));
3693 void pr_tpxheader(FILE *fp
, int indent
, const char *title
, const t_tpxheader
*sh
)
3695 if (available(fp
, sh
, indent
, title
))
3697 indent
= pr_title(fp
, indent
, title
);
3698 pr_indent(fp
, indent
);
3699 fprintf(fp
, "bIr = %spresent\n", sh
->bIr
? "" : "not ");
3700 pr_indent(fp
, indent
);
3701 fprintf(fp
, "bBox = %spresent\n", sh
->bBox
? "" : "not ");
3702 pr_indent(fp
, indent
);
3703 fprintf(fp
, "bTop = %spresent\n", sh
->bTop
? "" : "not ");
3704 pr_indent(fp
, indent
);
3705 fprintf(fp
, "bX = %spresent\n", sh
->bX
? "" : "not ");
3706 pr_indent(fp
, indent
);
3707 fprintf(fp
, "bV = %spresent\n", sh
->bV
? "" : "not ");
3708 pr_indent(fp
, indent
);
3709 fprintf(fp
, "bF = %spresent\n", sh
->bF
? "" : "not ");
3711 pr_indent(fp
, indent
);
3712 fprintf(fp
, "natoms = %d\n", sh
->natoms
);
3713 pr_indent(fp
, indent
);
3714 fprintf(fp
, "lambda = %e\n", sh
->lambda
);