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! */
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/fileio/filenm.h"
50 #include "gromacs/fileio/gmxfio.h"
51 #include "gromacs/fileio/gmxfio-xdr.h"
52 #include "gromacs/legacyheaders/copyrite.h"
53 #include "gromacs/legacyheaders/macros.h"
54 #include "gromacs/legacyheaders/names.h"
55 #include "gromacs/legacyheaders/txtdump.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/topology/atomprop.h"
58 #include "gromacs/topology/block.h"
59 #include "gromacs/topology/mtop_util.h"
60 #include "gromacs/topology/symtab.h"
61 #include "gromacs/topology/topology.h"
62 #include "gromacs/utility/cstringutil.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/futil.h"
65 #include "gromacs/utility/smalloc.h"
67 #define TPX_TAG_RELEASE "release"
69 /*! \brief Tag string for the file format written to run input files
70 * written by this version of the code.
72 * Change this if you want to change the run input format in a feature
73 * branch. This ensures that there will not be different run input
74 * formats around which cannot be distinguished, while not causing
75 * problems rebasing the feature branch onto upstream changes. When
76 * merging with mainstream GROMACS, set this tag string back to
77 * TPX_TAG_RELEASE, and instead add an element to tpxv and set
78 * tpx_version to that.
80 static const char *tpx_tag = TPX_TAG_RELEASE;
82 /*! \brief Enum of values that describe the contents of a tpr file
83 * whose format matches a version number
85 * The enum helps the code be more self-documenting and ensure merges
86 * do not silently resolve when two patches make the same bump. When
87 * adding new functionality, add a new element to the end of this
88 * enumeration, change the definition of tpx_version, and write code
89 * below that does the right thing according to the value of
92 tpxv_ComputationalElectrophysiology = 96, /**< support for ion/water position swaps (computational electrophysiology) */
93 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to gmx_int64_t */
94 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
95 tpxv_InteractiveMolecularDynamics, /**< interactive molecular dynamics (IMD) */
96 tpxv_RemoveObsoleteParameters1, /**< remove optimize_fft, dihre_fc, nstcheckpoint */
97 tpxv_PullCoordTypeGeom, /**< add pull type and geometry per group and flat-bottom */
98 tpxv_PullGeomDirRel, /**< add pull geometry direction-relative */
99 tpxv_IntermolecularBondeds /**< permit inter-molecular bonded interactions in the topology */
102 /*! \brief Version number of the file format written to run input
103 * files by this version of the code.
105 * The tpx_version number should be increased whenever the file format
106 * in the main development branch changes, generally to the highest
107 * value present in tpxv. Backward compatibility for reading old run
108 * input files is maintained by checking this version number against
109 * that of the file and then using the correct code path.
111 * When developing a feature branch that needs to change the run input
112 * file format, change tpx_tag instead. */
113 static const int tpx_version = tpxv_IntermolecularBondeds;
116 /* This number should only be increased when you edit the TOPOLOGY section
117 * or the HEADER of the tpx format.
118 * This way we can maintain forward compatibility too for all analysis tools
119 * and/or external programs that only need to know the atom/residue names,
120 * charges, and bond connectivity.
122 * It first appeared in tpx version 26, when I also moved the inputrecord
123 * to the end of the tpx file, so we can just skip it if we only
126 * In particular, it must be increased when adding new elements to
127 * ftupd, so that old code can read new .tpr files.
129 static const int tpx_generation = 26;
131 /* This number should be the most recent backwards incompatible version
132 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
134 static const int tpx_incompatible_version = 9;
138 /* Struct used to maintain tpx compatibility when function types are added */
140 int fvnr; /* file version number in which the function type first appeared */
141 int ftype; /* function type */
145 * The entries should be ordered in:
146 * 1. ascending file version number
147 * 2. ascending function type number
149 /*static const t_ftupd ftupd[] = {
150 { 20, F_CUBICBONDS },
154 { 22, F_DISRESVIOL },
160 { 26, F_DIHRESVIOL },
161 { 30, F_CROSS_BOND_BONDS },
162 { 30, F_CROSS_BOND_ANGLES },
163 { 30, F_UREY_BRADLEY },
164 { 30, F_POLARIZATION },
168 * The entries should be ordered in:
169 * 1. ascending function type number
170 * 2. ascending file version number
172 /* question; what is the purpose of the commented code above? */
173 static const t_ftupd ftupd[] = {
174 { 20, F_CUBICBONDS },
179 { 43, F_TABBONDSNC },
180 { 70, F_RESTRBONDS },
181 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
182 { 76, F_LINEAR_ANGLES },
183 { 30, F_CROSS_BOND_BONDS },
184 { 30, F_CROSS_BOND_ANGLES },
185 { 30, F_UREY_BRADLEY },
186 { 34, F_QUARTIC_ANGLES },
188 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
189 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
198 { 72, F_NPSOLVATION },
200 { 41, F_LJC_PAIRS_NB },
203 { 32, F_COUL_RECIP },
206 { 30, F_POLARIZATION },
209 { 22, F_DISRESVIOL },
213 { 26, F_DIHRESVIOL },
218 { 46, F_ECONSERVED },
219 { 69, F_VTEMP_NOLONGERUSED},
221 { 54, F_DVDL_CONSTR },
222 { 76, F_ANHARM_POL },
225 { 79, F_DVDL_BONDED, },
226 { 79, F_DVDL_RESTRAINT },
227 { 79, F_DVDL_TEMPERATURE },
229 #define NFTUPD asize(ftupd)
231 /* Needed for backward compatibility */
234 /**************************************************************
236 * Now the higer level routines that do io of the structures and arrays
238 **************************************************************/
239 static void do_pullgrp_tpx_pre95(t_fileio *fio,
247 gmx_fio_do_int(fio, pgrp->nat);
250 snew(pgrp->ind, pgrp->nat);
252 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
253 gmx_fio_do_int(fio, pgrp->nweight);
256 snew(pgrp->weight, pgrp->nweight);
258 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
259 gmx_fio_do_int(fio, pgrp->pbcatom);
260 gmx_fio_do_rvec(fio, pcrd->vec);
261 clear_rvec(pcrd->origin);
262 gmx_fio_do_rvec(fio, tmp);
264 gmx_fio_do_real(fio, pcrd->rate);
265 gmx_fio_do_real(fio, pcrd->k);
266 if (file_version >= 56)
268 gmx_fio_do_real(fio, pcrd->kB);
276 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
278 gmx_fio_do_int(fio, pgrp->nat);
281 snew(pgrp->ind, pgrp->nat);
283 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
284 gmx_fio_do_int(fio, pgrp->nweight);
287 snew(pgrp->weight, pgrp->nweight);
289 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
290 gmx_fio_do_int(fio, pgrp->pbcatom);
293 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd, int file_version,
294 int ePullOld, int eGeomOld, ivec dimOld)
296 gmx_fio_do_int(fio, pcrd->group[0]);
297 gmx_fio_do_int(fio, pcrd->group[1]);
298 if (file_version >= tpxv_PullCoordTypeGeom)
300 gmx_fio_do_int(fio, pcrd->eType);
301 gmx_fio_do_int(fio, pcrd->eGeom);
302 if (pcrd->eGeom == epullgDIRRELATIVE)
304 gmx_fio_do_int(fio, pcrd->group[2]);
305 gmx_fio_do_int(fio, pcrd->group[3]);
307 gmx_fio_do_ivec(fio, pcrd->dim);
311 pcrd->eType = ePullOld;
312 pcrd->eGeom = eGeomOld;
313 copy_ivec(dimOld, pcrd->dim);
315 gmx_fio_do_rvec(fio, pcrd->origin);
316 gmx_fio_do_rvec(fio, pcrd->vec);
317 if (file_version >= tpxv_PullCoordTypeGeom)
319 gmx_fio_do_gmx_bool(fio, pcrd->bStart);
323 /* This parameter is only printed, but not actually used by mdrun */
324 pcrd->bStart = FALSE;
326 gmx_fio_do_real(fio, pcrd->init);
327 gmx_fio_do_real(fio, pcrd->rate);
328 gmx_fio_do_real(fio, pcrd->k);
329 gmx_fio_do_real(fio, pcrd->kB);
332 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
334 int n_lambda = fepvals->n_lambda;
336 /* reset the lambda calculation window */
337 fepvals->lambda_start_n = 0;
338 fepvals->lambda_stop_n = n_lambda;
339 if (file_version >= 79)
345 snew(expand->init_lambda_weights, n_lambda);
347 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
348 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
351 gmx_fio_do_int(fio, expand->nstexpanded);
352 gmx_fio_do_int(fio, expand->elmcmove);
353 gmx_fio_do_int(fio, expand->elamstats);
354 gmx_fio_do_int(fio, expand->lmc_repeats);
355 gmx_fio_do_int(fio, expand->gibbsdeltalam);
356 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
357 gmx_fio_do_int(fio, expand->lmc_seed);
358 gmx_fio_do_real(fio, expand->mc_temp);
359 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
360 gmx_fio_do_int(fio, expand->nstTij);
361 gmx_fio_do_int(fio, expand->minvarmin);
362 gmx_fio_do_int(fio, expand->c_range);
363 gmx_fio_do_real(fio, expand->wl_scale);
364 gmx_fio_do_real(fio, expand->wl_ratio);
365 gmx_fio_do_real(fio, expand->init_wl_delta);
366 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
367 gmx_fio_do_int(fio, expand->elmceq);
368 gmx_fio_do_int(fio, expand->equil_steps);
369 gmx_fio_do_int(fio, expand->equil_samples);
370 gmx_fio_do_int(fio, expand->equil_n_at_lam);
371 gmx_fio_do_real(fio, expand->equil_wl_delta);
372 gmx_fio_do_real(fio, expand->equil_ratio);
376 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
379 if (file_version >= 79)
381 gmx_fio_do_int(fio, simtemp->eSimTempScale);
382 gmx_fio_do_real(fio, simtemp->simtemp_high);
383 gmx_fio_do_real(fio, simtemp->simtemp_low);
388 snew(simtemp->temperatures, n_lambda);
390 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
395 static void do_imd(t_fileio *fio, t_IMD *imd, gmx_bool bRead)
397 gmx_fio_do_int(fio, imd->nat);
400 snew(imd->ind, imd->nat);
402 gmx_fio_ndo_int(fio, imd->ind, imd->nat);
405 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
407 /* i is defined in the ndo_double macro; use g to iterate. */
411 /* free energy values */
413 if (file_version >= 79)
415 gmx_fio_do_int(fio, fepvals->init_fep_state);
416 gmx_fio_do_double(fio, fepvals->init_lambda);
417 gmx_fio_do_double(fio, fepvals->delta_lambda);
419 else if (file_version >= 59)
421 gmx_fio_do_double(fio, fepvals->init_lambda);
422 gmx_fio_do_double(fio, fepvals->delta_lambda);
426 gmx_fio_do_real(fio, rdum);
427 fepvals->init_lambda = rdum;
428 gmx_fio_do_real(fio, rdum);
429 fepvals->delta_lambda = rdum;
431 if (file_version >= 79)
433 gmx_fio_do_int(fio, fepvals->n_lambda);
436 snew(fepvals->all_lambda, efptNR);
438 for (g = 0; g < efptNR; g++)
440 if (fepvals->n_lambda > 0)
444 snew(fepvals->all_lambda[g], fepvals->n_lambda);
446 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
447 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
449 else if (fepvals->init_lambda >= 0)
451 fepvals->separate_dvdl[efptFEP] = TRUE;
455 else if (file_version >= 64)
457 gmx_fio_do_int(fio, fepvals->n_lambda);
462 snew(fepvals->all_lambda, efptNR);
463 /* still allocate the all_lambda array's contents. */
464 for (g = 0; g < efptNR; g++)
466 if (fepvals->n_lambda > 0)
468 snew(fepvals->all_lambda[g], fepvals->n_lambda);
472 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
474 if (fepvals->init_lambda >= 0)
478 fepvals->separate_dvdl[efptFEP] = TRUE;
482 /* copy the contents of the efptFEP lambda component to all
483 the other components */
484 for (g = 0; g < efptNR; g++)
486 for (h = 0; h < fepvals->n_lambda; h++)
490 fepvals->all_lambda[g][h] =
491 fepvals->all_lambda[efptFEP][h];
500 fepvals->n_lambda = 0;
501 fepvals->all_lambda = NULL;
502 if (fepvals->init_lambda >= 0)
504 fepvals->separate_dvdl[efptFEP] = TRUE;
507 if (file_version >= 13)
509 gmx_fio_do_real(fio, fepvals->sc_alpha);
513 fepvals->sc_alpha = 0;
515 if (file_version >= 38)
517 gmx_fio_do_int(fio, fepvals->sc_power);
521 fepvals->sc_power = 2;
523 if (file_version >= 79)
525 gmx_fio_do_real(fio, fepvals->sc_r_power);
529 fepvals->sc_r_power = 6.0;
531 if (file_version >= 15)
533 gmx_fio_do_real(fio, fepvals->sc_sigma);
537 fepvals->sc_sigma = 0.3;
541 if (file_version >= 71)
543 fepvals->sc_sigma_min = fepvals->sc_sigma;
547 fepvals->sc_sigma_min = 0;
550 if (file_version >= 79)
552 gmx_fio_do_int(fio, fepvals->bScCoul);
556 fepvals->bScCoul = TRUE;
558 if (file_version >= 64)
560 gmx_fio_do_int(fio, fepvals->nstdhdl);
564 fepvals->nstdhdl = 1;
567 if (file_version >= 73)
569 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
570 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
574 fepvals->separate_dhdl_file = esepdhdlfileYES;
575 fepvals->dhdl_derivatives = edhdlderivativesYES;
577 if (file_version >= 71)
579 gmx_fio_do_int(fio, fepvals->dh_hist_size);
580 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
584 fepvals->dh_hist_size = 0;
585 fepvals->dh_hist_spacing = 0.1;
587 if (file_version >= 79)
589 gmx_fio_do_int(fio, fepvals->edHdLPrintEnergy);
593 fepvals->edHdLPrintEnergy = edHdLPrintEnergyNO;
596 /* handle lambda_neighbors */
597 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
599 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
600 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
601 (fepvals->init_lambda < 0) )
603 fepvals->lambda_start_n = (fepvals->init_fep_state -
604 fepvals->lambda_neighbors);
605 fepvals->lambda_stop_n = (fepvals->init_fep_state +
606 fepvals->lambda_neighbors + 1);
607 if (fepvals->lambda_start_n < 0)
609 fepvals->lambda_start_n = 0;;
611 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
613 fepvals->lambda_stop_n = fepvals->n_lambda;
618 fepvals->lambda_start_n = 0;
619 fepvals->lambda_stop_n = fepvals->n_lambda;
624 fepvals->lambda_start_n = 0;
625 fepvals->lambda_stop_n = fepvals->n_lambda;
629 static void do_pull(t_fileio *fio, pull_params_t *pull, gmx_bool bRead,
630 int file_version, int ePullOld)
636 if (file_version >= 95)
638 gmx_fio_do_int(fio, pull->ngroup);
640 gmx_fio_do_int(fio, pull->ncoord);
641 if (file_version < 95)
643 pull->ngroup = pull->ncoord + 1;
645 if (file_version < tpxv_PullCoordTypeGeom)
649 gmx_fio_do_int(fio, eGeomOld);
650 gmx_fio_do_ivec(fio, dimOld);
651 /* The inner cylinder radius, now removed */
652 gmx_fio_do_real(fio, dum);
654 gmx_fio_do_real(fio, pull->cylinder_r);
655 gmx_fio_do_real(fio, pull->constr_tol);
656 if (file_version >= 95)
658 gmx_fio_do_int(fio, pull->bPrintCOM1);
659 /* With file_version < 95 this value is set below */
661 if (file_version >= tpxv_PullCoordTypeGeom)
663 gmx_fio_do_int(fio, pull->bPrintCOM2);
664 gmx_fio_do_int(fio, pull->bPrintRefValue);
665 gmx_fio_do_int(fio, pull->bPrintComp);
669 pull->bPrintCOM2 = FALSE;
670 pull->bPrintRefValue = FALSE;
671 pull->bPrintComp = TRUE;
673 gmx_fio_do_int(fio, pull->nstxout);
674 gmx_fio_do_int(fio, pull->nstfout);
677 snew(pull->group, pull->ngroup);
678 snew(pull->coord, pull->ncoord);
680 if (file_version < 95)
682 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
683 if (eGeomOld == epullgDIRPBC)
685 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
687 if (eGeomOld > epullgDIRPBC)
692 for (g = 0; g < pull->ngroup; g++)
694 /* We read and ignore a pull coordinate for group 0 */
695 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[std::max(g-1, 0)],
696 bRead, file_version);
699 pull->coord[g-1].group[0] = 0;
700 pull->coord[g-1].group[1] = g;
704 pull->bPrintCOM1 = (pull->group[0].nat > 0);
708 for (g = 0; g < pull->ngroup; g++)
710 do_pull_group(fio, &pull->group[g], bRead);
712 for (g = 0; g < pull->ncoord; g++)
714 do_pull_coord(fio, &pull->coord[g],
715 file_version, ePullOld, eGeomOld, dimOld);
721 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
723 gmx_fio_do_int(fio, rotg->eType);
724 gmx_fio_do_int(fio, rotg->bMassW);
725 gmx_fio_do_int(fio, rotg->nat);
728 snew(rotg->ind, rotg->nat);
730 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
733 snew(rotg->x_ref, rotg->nat);
735 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
736 gmx_fio_do_rvec(fio, rotg->vec);
737 gmx_fio_do_rvec(fio, rotg->pivot);
738 gmx_fio_do_real(fio, rotg->rate);
739 gmx_fio_do_real(fio, rotg->k);
740 gmx_fio_do_real(fio, rotg->slab_dist);
741 gmx_fio_do_real(fio, rotg->min_gaussian);
742 gmx_fio_do_real(fio, rotg->eps);
743 gmx_fio_do_int(fio, rotg->eFittype);
744 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
745 gmx_fio_do_real(fio, rotg->PotAngle_step);
748 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
752 gmx_fio_do_int(fio, rot->ngrp);
753 gmx_fio_do_int(fio, rot->nstrout);
754 gmx_fio_do_int(fio, rot->nstsout);
757 snew(rot->grp, rot->ngrp);
759 for (g = 0; g < rot->ngrp; g++)
761 do_rotgrp(fio, &rot->grp[g], bRead);
766 static void do_swapcoords(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead)
770 gmx_fio_do_int(fio, swap->nat);
771 gmx_fio_do_int(fio, swap->nat_sol);
772 for (j = 0; j < 2; j++)
774 gmx_fio_do_int(fio, swap->nat_split[j]);
775 gmx_fio_do_int(fio, swap->massw_split[j]);
777 gmx_fio_do_int(fio, swap->nstswap);
778 gmx_fio_do_int(fio, swap->nAverage);
779 gmx_fio_do_real(fio, swap->threshold);
780 gmx_fio_do_real(fio, swap->cyl0r);
781 gmx_fio_do_real(fio, swap->cyl0u);
782 gmx_fio_do_real(fio, swap->cyl0l);
783 gmx_fio_do_real(fio, swap->cyl1r);
784 gmx_fio_do_real(fio, swap->cyl1u);
785 gmx_fio_do_real(fio, swap->cyl1l);
789 snew(swap->ind, swap->nat);
790 snew(swap->ind_sol, swap->nat_sol);
791 for (j = 0; j < 2; j++)
793 snew(swap->ind_split[j], swap->nat_split[j]);
797 gmx_fio_ndo_int(fio, swap->ind, swap->nat);
798 gmx_fio_ndo_int(fio, swap->ind_sol, swap->nat_sol);
799 for (j = 0; j < 2; j++)
801 gmx_fio_ndo_int(fio, swap->ind_split[j], swap->nat_split[j]);
804 for (j = 0; j < eCompNR; j++)
806 gmx_fio_do_int(fio, swap->nanions[j]);
807 gmx_fio_do_int(fio, swap->ncations[j]);
813 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
814 int file_version, real *fudgeQQ)
816 int i, j, k, *tmp, idum = 0;
819 gmx_bool bSimAnn, bdum = 0;
820 real zerotemptime, finish_t, init_temp, finish_temp;
822 if (file_version != tpx_version)
824 /* Give a warning about features that are not accessible */
825 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
826 file_version, tpx_version);
834 if (file_version == 0)
839 /* Basic inputrec stuff */
840 gmx_fio_do_int(fio, ir->eI);
841 if (file_version >= 62)
843 gmx_fio_do_int64(fio, ir->nsteps);
847 gmx_fio_do_int(fio, idum);
851 if (file_version > 25)
853 if (file_version >= 62)
855 gmx_fio_do_int64(fio, ir->init_step);
859 gmx_fio_do_int(fio, idum);
860 ir->init_step = idum;
868 if (file_version >= 58)
870 gmx_fio_do_int(fio, ir->simulation_part);
874 ir->simulation_part = 1;
877 if (file_version >= 67)
879 gmx_fio_do_int(fio, ir->nstcalcenergy);
883 ir->nstcalcenergy = 1;
885 if (file_version < 53)
887 /* The pbc info has been moved out of do_inputrec,
888 * since we always want it, also without reading the inputrec.
890 gmx_fio_do_int(fio, ir->ePBC);
891 if ((file_version <= 15) && (ir->ePBC == 2))
895 if (file_version >= 45)
897 gmx_fio_do_int(fio, ir->bPeriodicMols);
904 ir->bPeriodicMols = TRUE;
908 ir->bPeriodicMols = FALSE;
912 if (file_version >= 81)
914 gmx_fio_do_int(fio, ir->cutoff_scheme);
915 if (file_version < 94)
917 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
922 ir->cutoff_scheme = ecutsGROUP;
924 gmx_fio_do_int(fio, ir->ns_type);
925 gmx_fio_do_int(fio, ir->nstlist);
926 gmx_fio_do_int(fio, idum); /* used to be ndelta; not used anymore */
927 if (file_version < 41)
929 gmx_fio_do_int(fio, idum);
930 gmx_fio_do_int(fio, idum);
932 if (file_version >= 45)
934 gmx_fio_do_real(fio, ir->rtpi);
940 gmx_fio_do_int(fio, ir->nstcomm);
941 if (file_version > 34)
943 gmx_fio_do_int(fio, ir->comm_mode);
945 else if (ir->nstcomm < 0)
947 ir->comm_mode = ecmANGULAR;
951 ir->comm_mode = ecmLINEAR;
953 ir->nstcomm = abs(ir->nstcomm);
955 /* ignore nstcheckpoint */
956 if (file_version > 25 && file_version < tpxv_RemoveObsoleteParameters1)
958 gmx_fio_do_int(fio, idum);
961 gmx_fio_do_int(fio, ir->nstcgsteep);
963 if (file_version >= 30)
965 gmx_fio_do_int(fio, ir->nbfgscorr);
972 gmx_fio_do_int(fio, ir->nstlog);
973 gmx_fio_do_int(fio, ir->nstxout);
974 gmx_fio_do_int(fio, ir->nstvout);
975 gmx_fio_do_int(fio, ir->nstfout);
976 gmx_fio_do_int(fio, ir->nstenergy);
977 gmx_fio_do_int(fio, ir->nstxout_compressed);
978 if (file_version >= 59)
980 gmx_fio_do_double(fio, ir->init_t);
981 gmx_fio_do_double(fio, ir->delta_t);
985 gmx_fio_do_real(fio, rdum);
987 gmx_fio_do_real(fio, rdum);
990 gmx_fio_do_real(fio, ir->x_compression_precision);
991 if (file_version < 19)
993 gmx_fio_do_int(fio, idum);
994 gmx_fio_do_int(fio, idum);
996 if (file_version < 18)
998 gmx_fio_do_int(fio, idum);
1000 if (file_version >= 81)
1002 gmx_fio_do_real(fio, ir->verletbuf_tol);
1006 ir->verletbuf_tol = 0;
1008 gmx_fio_do_real(fio, ir->rlist);
1009 if (file_version >= 67)
1011 gmx_fio_do_real(fio, ir->rlistlong);
1013 if (file_version >= 82 && file_version != 90)
1015 gmx_fio_do_int(fio, ir->nstcalclr);
1019 /* Calculate at NS steps */
1020 ir->nstcalclr = ir->nstlist;
1022 gmx_fio_do_int(fio, ir->coulombtype);
1023 if (file_version < 32 && ir->coulombtype == eelRF)
1025 ir->coulombtype = eelRF_NEC;
1027 if (file_version >= 81)
1029 gmx_fio_do_int(fio, ir->coulomb_modifier);
1033 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1035 gmx_fio_do_real(fio, ir->rcoulomb_switch);
1036 gmx_fio_do_real(fio, ir->rcoulomb);
1037 gmx_fio_do_int(fio, ir->vdwtype);
1038 if (file_version >= 81)
1040 gmx_fio_do_int(fio, ir->vdw_modifier);
1044 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1046 gmx_fio_do_real(fio, ir->rvdw_switch);
1047 gmx_fio_do_real(fio, ir->rvdw);
1048 if (file_version < 67)
1050 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1052 gmx_fio_do_int(fio, ir->eDispCorr);
1053 gmx_fio_do_real(fio, ir->epsilon_r);
1054 if (file_version >= 37)
1056 gmx_fio_do_real(fio, ir->epsilon_rf);
1060 if (EEL_RF(ir->coulombtype))
1062 ir->epsilon_rf = ir->epsilon_r;
1063 ir->epsilon_r = 1.0;
1067 ir->epsilon_rf = 1.0;
1070 if (file_version >= 29)
1072 gmx_fio_do_real(fio, ir->tabext);
1079 if (file_version > 25)
1081 gmx_fio_do_int(fio, ir->gb_algorithm);
1082 gmx_fio_do_int(fio, ir->nstgbradii);
1083 gmx_fio_do_real(fio, ir->rgbradii);
1084 gmx_fio_do_real(fio, ir->gb_saltconc);
1085 gmx_fio_do_int(fio, ir->implicit_solvent);
1089 ir->gb_algorithm = egbSTILL;
1092 ir->gb_saltconc = 0;
1093 ir->implicit_solvent = eisNO;
1095 if (file_version >= 55)
1097 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1098 gmx_fio_do_real(fio, ir->gb_obc_alpha);
1099 gmx_fio_do_real(fio, ir->gb_obc_beta);
1100 gmx_fio_do_real(fio, ir->gb_obc_gamma);
1101 if (file_version >= 60)
1103 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1104 gmx_fio_do_int(fio, ir->sa_algorithm);
1108 ir->gb_dielectric_offset = 0.009;
1109 ir->sa_algorithm = esaAPPROX;
1111 gmx_fio_do_real(fio, ir->sa_surface_tension);
1113 /* Override sa_surface_tension if it is not changed in the mpd-file */
1114 if (ir->sa_surface_tension < 0)
1116 if (ir->gb_algorithm == egbSTILL)
1118 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1120 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1122 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1129 /* Better use sensible values than insane (0.0) ones... */
1130 ir->gb_epsilon_solvent = 80;
1131 ir->gb_obc_alpha = 1.0;
1132 ir->gb_obc_beta = 0.8;
1133 ir->gb_obc_gamma = 4.85;
1134 ir->sa_surface_tension = 2.092;
1138 if (file_version >= 81)
1140 gmx_fio_do_real(fio, ir->fourier_spacing);
1144 ir->fourier_spacing = 0.0;
1146 gmx_fio_do_int(fio, ir->nkx);
1147 gmx_fio_do_int(fio, ir->nky);
1148 gmx_fio_do_int(fio, ir->nkz);
1149 gmx_fio_do_int(fio, ir->pme_order);
1150 gmx_fio_do_real(fio, ir->ewald_rtol);
1152 if (file_version >= 93)
1154 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1158 ir->ewald_rtol_lj = ir->ewald_rtol;
1161 if (file_version >= 24)
1163 gmx_fio_do_int(fio, ir->ewald_geometry);
1167 ir->ewald_geometry = eewg3D;
1170 if (file_version <= 17)
1172 ir->epsilon_surface = 0;
1173 if (file_version == 17)
1175 gmx_fio_do_int(fio, idum);
1180 gmx_fio_do_real(fio, ir->epsilon_surface);
1183 /* ignore bOptFFT */
1184 if (file_version < tpxv_RemoveObsoleteParameters1)
1186 gmx_fio_do_gmx_bool(fio, bdum);
1189 if (file_version >= 93)
1191 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1193 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1194 gmx_fio_do_int(fio, ir->etc);
1195 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1196 * but the values 0 and 1 still mean no and
1197 * berendsen temperature coupling, respectively.
1199 if (file_version >= 79)
1201 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1203 if (file_version >= 71)
1205 gmx_fio_do_int(fio, ir->nsttcouple);
1209 ir->nsttcouple = ir->nstcalcenergy;
1211 if (file_version <= 15)
1213 gmx_fio_do_int(fio, idum);
1215 if (file_version <= 17)
1217 gmx_fio_do_int(fio, ir->epct);
1218 if (file_version <= 15)
1222 ir->epct = epctSURFACETENSION;
1224 gmx_fio_do_int(fio, idum);
1227 /* we have removed the NO alternative at the beginning */
1231 ir->epct = epctISOTROPIC;
1235 ir->epc = epcBERENDSEN;
1240 gmx_fio_do_int(fio, ir->epc);
1241 gmx_fio_do_int(fio, ir->epct);
1243 if (file_version >= 71)
1245 gmx_fio_do_int(fio, ir->nstpcouple);
1249 ir->nstpcouple = ir->nstcalcenergy;
1251 gmx_fio_do_real(fio, ir->tau_p);
1252 if (file_version <= 15)
1254 gmx_fio_do_rvec(fio, vdum);
1255 clear_mat(ir->ref_p);
1256 for (i = 0; i < DIM; i++)
1258 ir->ref_p[i][i] = vdum[i];
1263 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1264 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1265 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1267 if (file_version <= 15)
1269 gmx_fio_do_rvec(fio, vdum);
1270 clear_mat(ir->compress);
1271 for (i = 0; i < DIM; i++)
1273 ir->compress[i][i] = vdum[i];
1278 gmx_fio_do_rvec(fio, ir->compress[XX]);
1279 gmx_fio_do_rvec(fio, ir->compress[YY]);
1280 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1282 if (file_version >= 47)
1284 gmx_fio_do_int(fio, ir->refcoord_scaling);
1285 gmx_fio_do_rvec(fio, ir->posres_com);
1286 gmx_fio_do_rvec(fio, ir->posres_comB);
1290 ir->refcoord_scaling = erscNO;
1291 clear_rvec(ir->posres_com);
1292 clear_rvec(ir->posres_comB);
1294 if ((file_version > 25) && (file_version < 79))
1296 gmx_fio_do_int(fio, ir->andersen_seed);
1300 ir->andersen_seed = 0;
1302 if (file_version < 26)
1304 gmx_fio_do_gmx_bool(fio, bSimAnn);
1305 gmx_fio_do_real(fio, zerotemptime);
1308 if (file_version < 37)
1310 gmx_fio_do_real(fio, rdum);
1313 gmx_fio_do_real(fio, ir->shake_tol);
1314 if (file_version < 54)
1316 gmx_fio_do_real(fio, *fudgeQQ);
1319 gmx_fio_do_int(fio, ir->efep);
1320 if (file_version <= 14 && ir->efep != efepNO)
1324 do_fepvals(fio, ir->fepvals, bRead, file_version);
1326 if (file_version >= 79)
1328 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1331 ir->bSimTemp = TRUE;
1336 ir->bSimTemp = FALSE;
1340 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1343 if (file_version >= 79)
1345 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1348 ir->bExpanded = TRUE;
1352 ir->bExpanded = FALSE;
1357 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1359 if (file_version >= 57)
1361 gmx_fio_do_int(fio, ir->eDisre);
1363 gmx_fio_do_int(fio, ir->eDisreWeighting);
1364 if (file_version < 22)
1366 if (ir->eDisreWeighting == 0)
1368 ir->eDisreWeighting = edrwEqual;
1372 ir->eDisreWeighting = edrwConservative;
1375 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1376 gmx_fio_do_real(fio, ir->dr_fc);
1377 gmx_fio_do_real(fio, ir->dr_tau);
1378 gmx_fio_do_int(fio, ir->nstdisreout);
1379 if (file_version >= 22)
1381 gmx_fio_do_real(fio, ir->orires_fc);
1382 gmx_fio_do_real(fio, ir->orires_tau);
1383 gmx_fio_do_int(fio, ir->nstorireout);
1389 ir->nstorireout = 0;
1392 /* ignore dihre_fc */
1393 if (file_version >= 26 && file_version < 79)
1395 gmx_fio_do_real(fio, rdum);
1396 if (file_version < 56)
1398 gmx_fio_do_real(fio, rdum);
1399 gmx_fio_do_int(fio, idum);
1403 gmx_fio_do_real(fio, ir->em_stepsize);
1404 gmx_fio_do_real(fio, ir->em_tol);
1405 if (file_version >= 22)
1407 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1411 ir->bShakeSOR = TRUE;
1413 if (file_version >= 11)
1415 gmx_fio_do_int(fio, ir->niter);
1420 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1423 if (file_version >= 21)
1425 gmx_fio_do_real(fio, ir->fc_stepsize);
1429 ir->fc_stepsize = 0;
1431 gmx_fio_do_int(fio, ir->eConstrAlg);
1432 gmx_fio_do_int(fio, ir->nProjOrder);
1433 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1434 if (file_version <= 14)
1436 gmx_fio_do_int(fio, idum);
1438 if (file_version >= 26)
1440 gmx_fio_do_int(fio, ir->nLincsIter);
1445 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1448 if (file_version < 33)
1450 gmx_fio_do_real(fio, bd_temp);
1452 gmx_fio_do_real(fio, ir->bd_fric);
1453 if (file_version >= tpxv_Use64BitRandomSeed)
1455 gmx_fio_do_int64(fio, ir->ld_seed);
1459 gmx_fio_do_int(fio, idum);
1462 if (file_version >= 33)
1464 for (i = 0; i < DIM; i++)
1466 gmx_fio_do_rvec(fio, ir->deform[i]);
1471 for (i = 0; i < DIM; i++)
1473 clear_rvec(ir->deform[i]);
1476 if (file_version >= 14)
1478 gmx_fio_do_real(fio, ir->cos_accel);
1484 gmx_fio_do_int(fio, ir->userint1);
1485 gmx_fio_do_int(fio, ir->userint2);
1486 gmx_fio_do_int(fio, ir->userint3);
1487 gmx_fio_do_int(fio, ir->userint4);
1488 gmx_fio_do_real(fio, ir->userreal1);
1489 gmx_fio_do_real(fio, ir->userreal2);
1490 gmx_fio_do_real(fio, ir->userreal3);
1491 gmx_fio_do_real(fio, ir->userreal4);
1494 if (file_version >= 77)
1496 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1501 snew(ir->adress, 1);
1503 gmx_fio_do_int(fio, ir->adress->type);
1504 gmx_fio_do_real(fio, ir->adress->const_wf);
1505 gmx_fio_do_real(fio, ir->adress->ex_width);
1506 gmx_fio_do_real(fio, ir->adress->hy_width);
1507 gmx_fio_do_int(fio, ir->adress->icor);
1508 gmx_fio_do_int(fio, ir->adress->site);
1509 gmx_fio_do_rvec(fio, ir->adress->refs);
1510 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1511 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1512 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1513 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1517 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1519 if (ir->adress->n_tf_grps > 0)
1521 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1525 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1527 if (ir->adress->n_energy_grps > 0)
1529 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1535 ir->bAdress = FALSE;
1539 if (file_version >= 48)
1543 if (file_version >= tpxv_PullCoordTypeGeom)
1545 gmx_fio_do_gmx_bool(fio, ir->bPull);
1549 gmx_fio_do_int(fio, ePullOld);
1550 ir->bPull = (ePullOld > 0);
1551 /* We removed the first ePull=ePullNo for the enum */
1560 do_pull(fio, ir->pull, bRead, file_version, ePullOld);
1568 /* Enforced rotation */
1569 if (file_version >= 74)
1571 gmx_fio_do_int(fio, ir->bRot);
1572 if (ir->bRot == TRUE)
1578 do_rot(fio, ir->rot, bRead);
1586 /* Interactive molecular dynamics */
1587 if (file_version >= tpxv_InteractiveMolecularDynamics)
1589 gmx_fio_do_int(fio, ir->bIMD);
1590 if (TRUE == ir->bIMD)
1596 do_imd(fio, ir->imd, bRead);
1601 /* We don't support IMD sessions for old .tpr files */
1606 gmx_fio_do_int(fio, ir->opts.ngtc);
1607 if (file_version >= 69)
1609 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1613 ir->opts.nhchainlength = 1;
1615 gmx_fio_do_int(fio, ir->opts.ngacc);
1616 gmx_fio_do_int(fio, ir->opts.ngfrz);
1617 gmx_fio_do_int(fio, ir->opts.ngener);
1621 snew(ir->opts.nrdf, ir->opts.ngtc);
1622 snew(ir->opts.ref_t, ir->opts.ngtc);
1623 snew(ir->opts.annealing, ir->opts.ngtc);
1624 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1625 snew(ir->opts.anneal_time, ir->opts.ngtc);
1626 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1627 snew(ir->opts.tau_t, ir->opts.ngtc);
1628 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1629 snew(ir->opts.acc, ir->opts.ngacc);
1630 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1632 if (ir->opts.ngtc > 0)
1634 if (bRead && file_version < 13)
1636 snew(tmp, ir->opts.ngtc);
1637 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1638 for (i = 0; i < ir->opts.ngtc; i++)
1640 ir->opts.nrdf[i] = tmp[i];
1646 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1648 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1649 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1650 if (file_version < 33 && ir->eI == eiBD)
1652 for (i = 0; i < ir->opts.ngtc; i++)
1654 ir->opts.tau_t[i] = bd_temp;
1658 if (ir->opts.ngfrz > 0)
1660 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1662 if (ir->opts.ngacc > 0)
1664 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1666 if (file_version >= 12)
1668 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1669 ir->opts.ngener*ir->opts.ngener);
1672 if (bRead && file_version < 26)
1674 for (i = 0; i < ir->opts.ngtc; i++)
1678 ir->opts.annealing[i] = eannSINGLE;
1679 ir->opts.anneal_npoints[i] = 2;
1680 snew(ir->opts.anneal_time[i], 2);
1681 snew(ir->opts.anneal_temp[i], 2);
1682 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1683 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1684 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1685 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1686 ir->opts.anneal_time[i][0] = ir->init_t;
1687 ir->opts.anneal_time[i][1] = finish_t;
1688 ir->opts.anneal_temp[i][0] = init_temp;
1689 ir->opts.anneal_temp[i][1] = finish_temp;
1693 ir->opts.annealing[i] = eannNO;
1694 ir->opts.anneal_npoints[i] = 0;
1700 /* file version 26 or later */
1701 /* First read the lists with annealing and npoints for each group */
1702 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1703 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1704 for (j = 0; j < (ir->opts.ngtc); j++)
1706 k = ir->opts.anneal_npoints[j];
1709 snew(ir->opts.anneal_time[j], k);
1710 snew(ir->opts.anneal_temp[j], k);
1712 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1713 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1717 if (file_version >= 45)
1719 gmx_fio_do_int(fio, ir->nwall);
1720 gmx_fio_do_int(fio, ir->wall_type);
1721 if (file_version >= 50)
1723 gmx_fio_do_real(fio, ir->wall_r_linpot);
1727 ir->wall_r_linpot = -1;
1729 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1730 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1731 gmx_fio_do_real(fio, ir->wall_density[0]);
1732 gmx_fio_do_real(fio, ir->wall_density[1]);
1733 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1739 ir->wall_atomtype[0] = -1;
1740 ir->wall_atomtype[1] = -1;
1741 ir->wall_density[0] = 0;
1742 ir->wall_density[1] = 0;
1743 ir->wall_ewald_zfac = 3;
1745 /* Cosine stuff for electric fields */
1746 for (j = 0; (j < DIM); j++)
1748 gmx_fio_do_int(fio, ir->ex[j].n);
1749 gmx_fio_do_int(fio, ir->et[j].n);
1752 snew(ir->ex[j].a, ir->ex[j].n);
1753 snew(ir->ex[j].phi, ir->ex[j].n);
1754 snew(ir->et[j].a, ir->et[j].n);
1755 snew(ir->et[j].phi, ir->et[j].n);
1757 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1758 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1759 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1760 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1764 if (file_version >= tpxv_ComputationalElectrophysiology)
1766 gmx_fio_do_int(fio, ir->eSwapCoords);
1767 if (ir->eSwapCoords != eswapNO)
1773 do_swapcoords(fio, ir->swap, bRead);
1778 if (file_version >= 39)
1780 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1781 gmx_fio_do_int(fio, ir->QMMMscheme);
1782 gmx_fio_do_real(fio, ir->scalefactor);
1783 gmx_fio_do_int(fio, ir->opts.ngQM);
1786 snew(ir->opts.QMmethod, ir->opts.ngQM);
1787 snew(ir->opts.QMbasis, ir->opts.ngQM);
1788 snew(ir->opts.QMcharge, ir->opts.ngQM);
1789 snew(ir->opts.QMmult, ir->opts.ngQM);
1790 snew(ir->opts.bSH, ir->opts.ngQM);
1791 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1792 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1793 snew(ir->opts.SAon, ir->opts.ngQM);
1794 snew(ir->opts.SAoff, ir->opts.ngQM);
1795 snew(ir->opts.SAsteps, ir->opts.ngQM);
1796 snew(ir->opts.bOPT, ir->opts.ngQM);
1797 snew(ir->opts.bTS, ir->opts.ngQM);
1799 if (ir->opts.ngQM > 0)
1801 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1802 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1803 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1804 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1805 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1806 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1807 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1808 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1809 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1810 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1811 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1812 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1814 /* end of QMMM stuff */
1819 static void do_harm(t_fileio *fio, t_iparams *iparams)
1821 gmx_fio_do_real(fio, iparams->harmonic.rA);
1822 gmx_fio_do_real(fio, iparams->harmonic.krA);
1823 gmx_fio_do_real(fio, iparams->harmonic.rB);
1824 gmx_fio_do_real(fio, iparams->harmonic.krB);
1827 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1828 gmx_bool bRead, int file_version)
1841 do_harm(fio, iparams);
1842 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1844 /* Correct incorrect storage of parameters */
1845 iparams->pdihs.phiB = iparams->pdihs.phiA;
1846 iparams->pdihs.cpB = iparams->pdihs.cpA;
1850 gmx_fio_do_real(fio, iparams->harmonic.rA);
1851 gmx_fio_do_real(fio, iparams->harmonic.krA);
1853 case F_LINEAR_ANGLES:
1854 gmx_fio_do_real(fio, iparams->linangle.klinA);
1855 gmx_fio_do_real(fio, iparams->linangle.aA);
1856 gmx_fio_do_real(fio, iparams->linangle.klinB);
1857 gmx_fio_do_real(fio, iparams->linangle.aB);
1860 gmx_fio_do_real(fio, iparams->fene.bm);
1861 gmx_fio_do_real(fio, iparams->fene.kb);
1865 gmx_fio_do_real(fio, iparams->restraint.lowA);
1866 gmx_fio_do_real(fio, iparams->restraint.up1A);
1867 gmx_fio_do_real(fio, iparams->restraint.up2A);
1868 gmx_fio_do_real(fio, iparams->restraint.kA);
1869 gmx_fio_do_real(fio, iparams->restraint.lowB);
1870 gmx_fio_do_real(fio, iparams->restraint.up1B);
1871 gmx_fio_do_real(fio, iparams->restraint.up2B);
1872 gmx_fio_do_real(fio, iparams->restraint.kB);
1878 gmx_fio_do_real(fio, iparams->tab.kA);
1879 gmx_fio_do_int(fio, iparams->tab.table);
1880 gmx_fio_do_real(fio, iparams->tab.kB);
1882 case F_CROSS_BOND_BONDS:
1883 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1884 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1885 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1887 case F_CROSS_BOND_ANGLES:
1888 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1889 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1890 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1891 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1893 case F_UREY_BRADLEY:
1894 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1895 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1896 gmx_fio_do_real(fio, iparams->u_b.r13A);
1897 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1898 if (file_version >= 79)
1900 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1901 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1902 gmx_fio_do_real(fio, iparams->u_b.r13B);
1903 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1907 iparams->u_b.thetaB = iparams->u_b.thetaA;
1908 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1909 iparams->u_b.r13B = iparams->u_b.r13A;
1910 iparams->u_b.kUBB = iparams->u_b.kUBA;
1913 case F_QUARTIC_ANGLES:
1914 gmx_fio_do_real(fio, iparams->qangle.theta);
1915 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1918 gmx_fio_do_real(fio, iparams->bham.a);
1919 gmx_fio_do_real(fio, iparams->bham.b);
1920 gmx_fio_do_real(fio, iparams->bham.c);
1923 gmx_fio_do_real(fio, iparams->morse.b0A);
1924 gmx_fio_do_real(fio, iparams->morse.cbA);
1925 gmx_fio_do_real(fio, iparams->morse.betaA);
1926 if (file_version >= 79)
1928 gmx_fio_do_real(fio, iparams->morse.b0B);
1929 gmx_fio_do_real(fio, iparams->morse.cbB);
1930 gmx_fio_do_real(fio, iparams->morse.betaB);
1934 iparams->morse.b0B = iparams->morse.b0A;
1935 iparams->morse.cbB = iparams->morse.cbA;
1936 iparams->morse.betaB = iparams->morse.betaA;
1940 gmx_fio_do_real(fio, iparams->cubic.b0);
1941 gmx_fio_do_real(fio, iparams->cubic.kb);
1942 gmx_fio_do_real(fio, iparams->cubic.kcub);
1946 case F_POLARIZATION:
1947 gmx_fio_do_real(fio, iparams->polarize.alpha);
1950 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1951 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1952 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1955 if (file_version < 31)
1957 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1959 gmx_fio_do_real(fio, iparams->wpol.al_x);
1960 gmx_fio_do_real(fio, iparams->wpol.al_y);
1961 gmx_fio_do_real(fio, iparams->wpol.al_z);
1962 gmx_fio_do_real(fio, iparams->wpol.rOH);
1963 gmx_fio_do_real(fio, iparams->wpol.rHH);
1964 gmx_fio_do_real(fio, iparams->wpol.rOD);
1967 gmx_fio_do_real(fio, iparams->thole.a);
1968 gmx_fio_do_real(fio, iparams->thole.alpha1);
1969 gmx_fio_do_real(fio, iparams->thole.alpha2);
1970 gmx_fio_do_real(fio, iparams->thole.rfac);
1973 gmx_fio_do_real(fio, iparams->lj.c6);
1974 gmx_fio_do_real(fio, iparams->lj.c12);
1977 gmx_fio_do_real(fio, iparams->lj14.c6A);
1978 gmx_fio_do_real(fio, iparams->lj14.c12A);
1979 gmx_fio_do_real(fio, iparams->lj14.c6B);
1980 gmx_fio_do_real(fio, iparams->lj14.c12B);
1983 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1984 gmx_fio_do_real(fio, iparams->ljc14.qi);
1985 gmx_fio_do_real(fio, iparams->ljc14.qj);
1986 gmx_fio_do_real(fio, iparams->ljc14.c6);
1987 gmx_fio_do_real(fio, iparams->ljc14.c12);
1989 case F_LJC_PAIRS_NB:
1990 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1991 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1992 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1993 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1999 gmx_fio_do_real(fio, iparams->pdihs.phiA);
2000 gmx_fio_do_real(fio, iparams->pdihs.cpA);
2001 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
2003 /* Read the incorrectly stored multiplicity */
2004 gmx_fio_do_real(fio, iparams->harmonic.rB);
2005 gmx_fio_do_real(fio, iparams->harmonic.krB);
2006 iparams->pdihs.phiB = iparams->pdihs.phiA;
2007 iparams->pdihs.cpB = iparams->pdihs.cpA;
2011 gmx_fio_do_real(fio, iparams->pdihs.phiB);
2012 gmx_fio_do_real(fio, iparams->pdihs.cpB);
2013 gmx_fio_do_int(fio, iparams->pdihs.mult);
2017 gmx_fio_do_real(fio, iparams->pdihs.phiA);
2018 gmx_fio_do_real(fio, iparams->pdihs.cpA);
2021 gmx_fio_do_int(fio, iparams->disres.label);
2022 gmx_fio_do_int(fio, iparams->disres.type);
2023 gmx_fio_do_real(fio, iparams->disres.low);
2024 gmx_fio_do_real(fio, iparams->disres.up1);
2025 gmx_fio_do_real(fio, iparams->disres.up2);
2026 gmx_fio_do_real(fio, iparams->disres.kfac);
2029 gmx_fio_do_int(fio, iparams->orires.ex);
2030 gmx_fio_do_int(fio, iparams->orires.label);
2031 gmx_fio_do_int(fio, iparams->orires.power);
2032 gmx_fio_do_real(fio, iparams->orires.c);
2033 gmx_fio_do_real(fio, iparams->orires.obs);
2034 gmx_fio_do_real(fio, iparams->orires.kfac);
2037 if (file_version < 82)
2039 gmx_fio_do_int(fio, idum);
2040 gmx_fio_do_int(fio, idum);
2042 gmx_fio_do_real(fio, iparams->dihres.phiA);
2043 gmx_fio_do_real(fio, iparams->dihres.dphiA);
2044 gmx_fio_do_real(fio, iparams->dihres.kfacA);
2045 if (file_version >= 82)
2047 gmx_fio_do_real(fio, iparams->dihres.phiB);
2048 gmx_fio_do_real(fio, iparams->dihres.dphiB);
2049 gmx_fio_do_real(fio, iparams->dihres.kfacB);
2053 iparams->dihres.phiB = iparams->dihres.phiA;
2054 iparams->dihres.dphiB = iparams->dihres.dphiA;
2055 iparams->dihres.kfacB = iparams->dihres.kfacA;
2059 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
2060 gmx_fio_do_rvec(fio, iparams->posres.fcA);
2061 if (bRead && file_version < 27)
2063 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
2064 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
2068 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
2069 gmx_fio_do_rvec(fio, iparams->posres.fcB);
2073 gmx_fio_do_int(fio, iparams->fbposres.geom);
2074 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
2075 gmx_fio_do_real(fio, iparams->fbposres.r);
2076 gmx_fio_do_real(fio, iparams->fbposres.k);
2079 gmx_fio_ndo_real(fio, iparams->cbtdihs.cbtcA, NR_CBTDIHS);
2082 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2083 if (file_version >= 25)
2085 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2089 /* Fourier dihedrals are internally represented
2090 * as Ryckaert-Bellemans since those are faster to compute.
2092 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2093 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2097 gmx_fio_do_real(fio, iparams->constr.dA);
2098 gmx_fio_do_real(fio, iparams->constr.dB);
2101 gmx_fio_do_real(fio, iparams->settle.doh);
2102 gmx_fio_do_real(fio, iparams->settle.dhh);
2105 gmx_fio_do_real(fio, iparams->vsite.a);
2110 gmx_fio_do_real(fio, iparams->vsite.a);
2111 gmx_fio_do_real(fio, iparams->vsite.b);
2116 gmx_fio_do_real(fio, iparams->vsite.a);
2117 gmx_fio_do_real(fio, iparams->vsite.b);
2118 gmx_fio_do_real(fio, iparams->vsite.c);
2121 gmx_fio_do_int(fio, iparams->vsiten.n);
2122 gmx_fio_do_real(fio, iparams->vsiten.a);
2127 /* We got rid of some parameters in version 68 */
2128 if (bRead && file_version < 68)
2130 gmx_fio_do_real(fio, rdum);
2131 gmx_fio_do_real(fio, rdum);
2132 gmx_fio_do_real(fio, rdum);
2133 gmx_fio_do_real(fio, rdum);
2135 gmx_fio_do_real(fio, iparams->gb.sar);
2136 gmx_fio_do_real(fio, iparams->gb.st);
2137 gmx_fio_do_real(fio, iparams->gb.pi);
2138 gmx_fio_do_real(fio, iparams->gb.gbr);
2139 gmx_fio_do_real(fio, iparams->gb.bmlt);
2142 gmx_fio_do_int(fio, iparams->cmap.cmapA);
2143 gmx_fio_do_int(fio, iparams->cmap.cmapB);
2146 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
2147 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
2151 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version)
2155 if (file_version < 44)
2157 for (i = 0; i < MAXNODES; i++)
2159 gmx_fio_do_int(fio, idum);
2162 gmx_fio_do_int(fio, ilist->nr);
2165 snew(ilist->iatoms, ilist->nr);
2167 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2170 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2171 gmx_bool bRead, int file_version)
2176 gmx_fio_do_int(fio, ffparams->atnr);
2177 if (file_version < 57)
2179 gmx_fio_do_int(fio, idum);
2181 gmx_fio_do_int(fio, ffparams->ntypes);
2184 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2185 ffparams->atnr, ffparams->ntypes);
2189 snew(ffparams->functype, ffparams->ntypes);
2190 snew(ffparams->iparams, ffparams->ntypes);
2192 /* Read/write all the function types */
2193 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2196 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2199 if (file_version >= 66)
2201 gmx_fio_do_double(fio, ffparams->reppow);
2205 ffparams->reppow = 12.0;
2208 if (file_version >= 57)
2210 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2213 /* Check whether all these function types are supported by the code.
2214 * In practice the code is backwards compatible, which means that the
2215 * numbering may have to be altered from old numbering to new numbering
2217 for (i = 0; (i < ffparams->ntypes); i++)
2221 /* Loop over file versions */
2222 for (k = 0; (k < NFTUPD); k++)
2224 /* Compare the read file_version to the update table */
2225 if ((file_version < ftupd[k].fvnr) &&
2226 (ffparams->functype[i] >= ftupd[k].ftype))
2228 ffparams->functype[i] += 1;
2231 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2232 i, ffparams->functype[i],
2233 interaction_function[ftupd[k].ftype].longname);
2240 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2244 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2249 static void add_settle_atoms(t_ilist *ilist)
2253 /* Settle used to only store the first atom: add the other two */
2254 srenew(ilist->iatoms, 2*ilist->nr);
2255 for (i = ilist->nr/2-1; i >= 0; i--)
2257 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2258 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2259 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2260 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2262 ilist->nr = 2*ilist->nr;
2265 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2272 for (j = 0; (j < F_NRE); j++)
2277 for (k = 0; k < NFTUPD; k++)
2279 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2288 ilist[j].iatoms = NULL;
2292 do_ilist(fio, &ilist[j], bRead, file_version);
2293 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2295 add_settle_atoms(&ilist[j]);
2299 if (bRead && gmx_debug_at)
2300 pr_ilist(debug,0,interaction_function[j].longname,
2301 functype,&ilist[j],TRUE);
2306 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2307 gmx_bool bRead, int file_version)
2309 do_ffparams(fio, ffparams, bRead, file_version);
2311 if (file_version >= 54)
2313 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2316 do_ilists(fio, molt->ilist, bRead, file_version);
2319 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2321 int i, idum, dum_nra, *dum_a;
2323 if (file_version < 44)
2325 for (i = 0; i < MAXNODES; i++)
2327 gmx_fio_do_int(fio, idum);
2330 gmx_fio_do_int(fio, block->nr);
2331 if (file_version < 51)
2333 gmx_fio_do_int(fio, dum_nra);
2337 if ((block->nalloc_index > 0) && (NULL != block->index))
2339 sfree(block->index);
2341 block->nalloc_index = block->nr+1;
2342 snew(block->index, block->nalloc_index);
2344 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2346 if (file_version < 51 && dum_nra > 0)
2348 snew(dum_a, dum_nra);
2349 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2354 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2359 if (file_version < 44)
2361 for (i = 0; i < MAXNODES; i++)
2363 gmx_fio_do_int(fio, idum);
2366 gmx_fio_do_int(fio, block->nr);
2367 gmx_fio_do_int(fio, block->nra);
2370 block->nalloc_index = block->nr+1;
2371 snew(block->index, block->nalloc_index);
2372 block->nalloc_a = block->nra;
2373 snew(block->a, block->nalloc_a);
2375 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2376 gmx_fio_ndo_int(fio, block->a, block->nra);
2379 /* This is a primitive routine to make it possible to translate atomic numbers
2380 * to element names when reading TPR files, without making the Gromacs library
2381 * directory a dependency on mdrun (which is the case if we need elements.dat).
2384 atomicnumber_to_element(int atomicnumber)
2388 /* This does not have to be complete, so we only include elements likely
2389 * to occur in PDB files.
2391 switch (atomicnumber)
2393 case 1: p = "H"; break;
2394 case 5: p = "B"; break;
2395 case 6: p = "C"; break;
2396 case 7: p = "N"; break;
2397 case 8: p = "O"; break;
2398 case 9: p = "F"; break;
2399 case 11: p = "Na"; break;
2400 case 12: p = "Mg"; break;
2401 case 15: p = "P"; break;
2402 case 16: p = "S"; break;
2403 case 17: p = "Cl"; break;
2404 case 18: p = "Ar"; break;
2405 case 19: p = "K"; break;
2406 case 20: p = "Ca"; break;
2407 case 25: p = "Mn"; break;
2408 case 26: p = "Fe"; break;
2409 case 28: p = "Ni"; break;
2410 case 29: p = "Cu"; break;
2411 case 30: p = "Zn"; break;
2412 case 35: p = "Br"; break;
2413 case 47: p = "Ag"; break;
2414 default: p = ""; break;
2420 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2421 int file_version, gmx_groups_t *groups, int atnr)
2425 gmx_fio_do_real(fio, atom->m);
2426 gmx_fio_do_real(fio, atom->q);
2427 gmx_fio_do_real(fio, atom->mB);
2428 gmx_fio_do_real(fio, atom->qB);
2429 gmx_fio_do_ushort(fio, atom->type);
2430 gmx_fio_do_ushort(fio, atom->typeB);
2431 gmx_fio_do_int(fio, atom->ptype);
2432 gmx_fio_do_int(fio, atom->resind);
2433 if (file_version >= 52)
2435 gmx_fio_do_int(fio, atom->atomnumber);
2438 /* Set element string from atomic number if present.
2439 * This routine returns an empty string if the name is not found.
2441 std::strncpy(atom->elem, atomicnumber_to_element(atom->atomnumber), 4);
2442 /* avoid warnings about potentially unterminated string */
2443 atom->elem[3] = '\0';
2448 atom->atomnumber = NOTSET;
2450 if (file_version < 23)
2454 else if (file_version < 39)
2463 if (file_version < 57)
2465 unsigned char uchar[egcNR];
2466 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2467 for (i = myngrp; (i < ngrp); i++)
2471 /* Copy the old data format to the groups struct */
2472 for (i = 0; i < ngrp; i++)
2474 groups->grpnr[i][atnr] = uchar[i];
2479 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2484 if (file_version < 23)
2488 else if (file_version < 39)
2497 for (j = 0; (j < ngrp); j++)
2501 gmx_fio_do_int(fio, grps[j].nr);
2504 snew(grps[j].nm_ind, grps[j].nr);
2506 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2511 snew(grps[j].nm_ind, grps[j].nr);
2516 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2522 gmx_fio_do_int(fio, ls);
2523 *nm = get_symtab_handle(symtab, ls);
2527 ls = lookup_symtab(symtab, *nm);
2528 gmx_fio_do_int(fio, ls);
2532 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2537 for (j = 0; (j < nstr); j++)
2539 do_symstr(fio, &(nm[j]), bRead, symtab);
2543 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2544 t_symtab *symtab, int file_version)
2548 for (j = 0; (j < n); j++)
2550 do_symstr(fio, &(ri[j].name), bRead, symtab);
2551 if (file_version >= 63)
2553 gmx_fio_do_int(fio, ri[j].nr);
2554 gmx_fio_do_uchar(fio, ri[j].ic);
2564 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2566 gmx_groups_t *groups)
2570 gmx_fio_do_int(fio, atoms->nr);
2571 gmx_fio_do_int(fio, atoms->nres);
2572 if (file_version < 57)
2574 gmx_fio_do_int(fio, groups->ngrpname);
2575 for (i = 0; i < egcNR; i++)
2577 groups->ngrpnr[i] = atoms->nr;
2578 snew(groups->grpnr[i], groups->ngrpnr[i]);
2583 snew(atoms->atom, atoms->nr);
2584 snew(atoms->atomname, atoms->nr);
2585 snew(atoms->atomtype, atoms->nr);
2586 snew(atoms->atomtypeB, atoms->nr);
2587 snew(atoms->resinfo, atoms->nres);
2588 if (file_version < 57)
2590 snew(groups->grpname, groups->ngrpname);
2592 atoms->pdbinfo = NULL;
2594 for (i = 0; (i < atoms->nr); i++)
2596 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2598 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2599 if (bRead && (file_version <= 20))
2601 for (i = 0; i < atoms->nr; i++)
2603 atoms->atomtype[i] = put_symtab(symtab, "?");
2604 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2609 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2610 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2612 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2614 if (file_version < 57)
2616 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2618 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2622 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2623 gmx_bool bRead, t_symtab *symtab,
2628 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2629 gmx_fio_do_int(fio, groups->ngrpname);
2632 snew(groups->grpname, groups->ngrpname);
2634 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2635 for (g = 0; g < egcNR; g++)
2637 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2638 if (groups->ngrpnr[g] == 0)
2642 groups->grpnr[g] = NULL;
2649 snew(groups->grpnr[g], groups->ngrpnr[g]);
2651 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2656 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2661 if (file_version > 25)
2663 gmx_fio_do_int(fio, atomtypes->nr);
2667 snew(atomtypes->radius, j);
2668 snew(atomtypes->vol, j);
2669 snew(atomtypes->surftens, j);
2670 snew(atomtypes->atomnumber, j);
2671 snew(atomtypes->gb_radius, j);
2672 snew(atomtypes->S_hct, j);
2674 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2675 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2676 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2677 if (file_version >= 40)
2679 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2681 if (file_version >= 60)
2683 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2684 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2689 /* File versions prior to 26 cannot do GBSA,
2690 * so they dont use this structure
2693 atomtypes->radius = NULL;
2694 atomtypes->vol = NULL;
2695 atomtypes->surftens = NULL;
2696 atomtypes->atomnumber = NULL;
2697 atomtypes->gb_radius = NULL;
2698 atomtypes->S_hct = NULL;
2702 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2708 gmx_fio_do_int(fio, symtab->nr);
2712 snew(symtab->symbuf, 1);
2713 symbuf = symtab->symbuf;
2714 symbuf->bufsize = nr;
2715 snew(symbuf->buf, nr);
2716 for (i = 0; (i < nr); i++)
2718 gmx_fio_do_string(fio, buf);
2719 symbuf->buf[i] = gmx_strdup(buf);
2724 symbuf = symtab->symbuf;
2725 while (symbuf != NULL)
2727 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2729 gmx_fio_do_string(fio, symbuf->buf[i]);
2732 symbuf = symbuf->next;
2736 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2741 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2743 int i, j, ngrid, gs, nelem;
2745 gmx_fio_do_int(fio, cmap_grid->ngrid);
2746 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2748 ngrid = cmap_grid->ngrid;
2749 gs = cmap_grid->grid_spacing;
2754 snew(cmap_grid->cmapdata, ngrid);
2756 for (i = 0; i < cmap_grid->ngrid; i++)
2758 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2762 for (i = 0; i < cmap_grid->ngrid; i++)
2764 for (j = 0; j < nelem; j++)
2766 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2767 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2768 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2769 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2775 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2777 int m, a, a0, a1, r;
2781 /* We always assign a new chain number, but save the chain id characters
2782 * for larger molecules.
2784 #define CHAIN_MIN_ATOMS 15
2788 for (m = 0; m < mols->nr; m++)
2790 a0 = mols->index[m];
2791 a1 = mols->index[m+1];
2792 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2801 for (a = a0; a < a1; a++)
2803 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2804 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2809 /* Blank out the chain id if there was only one chain */
2812 for (r = 0; r < atoms->nres; r++)
2814 atoms->resinfo[r].chainid = ' ';
2819 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2820 t_symtab *symtab, int file_version,
2821 gmx_groups_t *groups)
2823 if (file_version >= 57)
2825 do_symstr(fio, &(molt->name), bRead, symtab);
2828 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2830 if (bRead && gmx_debug_at)
2832 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2835 if (file_version >= 57)
2837 do_ilists(fio, molt->ilist, bRead, file_version);
2839 do_block(fio, &molt->cgs, bRead, file_version);
2840 if (bRead && gmx_debug_at)
2842 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2846 /* This used to be in the atoms struct */
2847 do_blocka(fio, &molt->excls, bRead, file_version);
2850 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2852 gmx_fio_do_int(fio, molb->type);
2853 gmx_fio_do_int(fio, molb->nmol);
2854 gmx_fio_do_int(fio, molb->natoms_mol);
2855 /* Position restraint coordinates */
2856 gmx_fio_do_int(fio, molb->nposres_xA);
2857 if (molb->nposres_xA > 0)
2861 snew(molb->posres_xA, molb->nposres_xA);
2863 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2865 gmx_fio_do_int(fio, molb->nposres_xB);
2866 if (molb->nposres_xB > 0)
2870 snew(molb->posres_xB, molb->nposres_xB);
2872 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2877 static t_block mtop_mols(gmx_mtop_t *mtop)
2883 for (mb = 0; mb < mtop->nmolblock; mb++)
2885 mols.nr += mtop->molblock[mb].nmol;
2887 mols.nalloc_index = mols.nr + 1;
2888 snew(mols.index, mols.nalloc_index);
2893 for (mb = 0; mb < mtop->nmolblock; mb++)
2895 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2897 a += mtop->molblock[mb].natoms_mol;
2906 static void add_posres_molblock(gmx_mtop_t *mtop)
2911 gmx_molblock_t *molb;
2914 /* posres reference positions are stored in ip->posres (if present) and
2915 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2916 posres.pos0A are identical to fbposres.pos0. */
2917 il = &mtop->moltype[0].ilist[F_POSRES];
2918 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2919 if (il->nr == 0 && ilfb->nr == 0)
2925 for (i = 0; i < il->nr; i += 2)
2927 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2928 am = std::max(am, il->iatoms[i+1]);
2929 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2930 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2931 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2936 /* This loop is required if we have only flat-bottomed posres:
2938 - bFE == FALSE (no B-state for flat-bottomed posres) */
2941 for (i = 0; i < ilfb->nr; i += 2)
2943 am = std::max(am, ilfb->iatoms[i+1]);
2946 /* Make the posres coordinate block end at a molecule end */
2948 while (am >= mtop->mols.index[mol+1])
2952 molb = &mtop->molblock[0];
2953 molb->nposres_xA = mtop->mols.index[mol+1];
2954 snew(molb->posres_xA, molb->nposres_xA);
2957 molb->nposres_xB = molb->nposres_xA;
2958 snew(molb->posres_xB, molb->nposres_xB);
2962 molb->nposres_xB = 0;
2964 for (i = 0; i < il->nr; i += 2)
2966 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2967 a = il->iatoms[i+1];
2968 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2969 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2970 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2973 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2974 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2975 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2980 /* If only flat-bottomed posres are present, take reference pos from them.
2981 Here: bFE == FALSE */
2982 for (i = 0; i < ilfb->nr; i += 2)
2984 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2985 a = ilfb->iatoms[i+1];
2986 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2987 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2988 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2993 static void set_disres_npair(gmx_mtop_t *mtop)
2996 gmx_mtop_ilistloop_t iloop;
2997 t_ilist *ilist, *il;
3001 ip = mtop->ffparams.iparams;
3003 iloop = gmx_mtop_ilistloop_init(mtop);
3004 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
3006 il = &ilist[F_DISRES];
3012 for (i = 0; i < il->nr; i += 3)
3015 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
3017 ip[a[i]].disres.npair = npair;
3025 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
3035 do_symtab(fio, &(mtop->symtab), bRead);
3038 pr_symtab(debug, 0, "symtab", &mtop->symtab);
3041 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
3043 if (file_version >= 57)
3045 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
3047 gmx_fio_do_int(fio, mtop->nmoltype);
3055 snew(mtop->moltype, mtop->nmoltype);
3056 if (file_version < 57)
3058 mtop->moltype[0].name = mtop->name;
3061 for (mt = 0; mt < mtop->nmoltype; mt++)
3063 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
3067 if (file_version >= 57)
3069 gmx_fio_do_int(fio, mtop->nmolblock);
3073 mtop->nmolblock = 1;
3077 snew(mtop->molblock, mtop->nmolblock);
3079 if (file_version >= 57)
3081 for (mb = 0; mb < mtop->nmolblock; mb++)
3083 do_molblock(fio, &mtop->molblock[mb], bRead);
3085 gmx_fio_do_int(fio, mtop->natoms);
3089 mtop->molblock[0].type = 0;
3090 mtop->molblock[0].nmol = 1;
3091 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
3092 mtop->molblock[0].nposres_xA = 0;
3093 mtop->molblock[0].nposres_xB = 0;
3096 if (file_version >= tpxv_IntermolecularBondeds)
3098 gmx_fio_do_gmx_bool(fio, mtop->bIntermolecularInteractions);
3099 if (mtop->bIntermolecularInteractions)
3103 snew(mtop->intermolecular_ilist, F_NRE);
3105 do_ilists(fio, mtop->intermolecular_ilist, bRead, file_version);
3110 mtop->bIntermolecularInteractions = FALSE;
3113 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
3116 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
3119 if (file_version < 57)
3121 /* Debug statements are inside do_idef */
3122 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
3123 mtop->natoms = mtop->moltype[0].atoms.nr;
3126 if (file_version >= 65)
3128 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
3132 mtop->ffparams.cmap_grid.ngrid = 0;
3133 mtop->ffparams.cmap_grid.grid_spacing = 0;
3134 mtop->ffparams.cmap_grid.cmapdata = NULL;
3137 if (file_version >= 57)
3139 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
3142 if (file_version < 57)
3144 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
3145 if (bRead && gmx_debug_at)
3147 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
3149 do_block(fio, &mtop->mols, bRead, file_version);
3150 /* Add the posres coordinates to the molblock */
3151 add_posres_molblock(mtop);
3155 if (file_version >= 57)
3157 done_block(&mtop->mols);
3158 mtop->mols = mtop_mols(mtop);
3162 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
3166 if (file_version < 51)
3168 /* Here used to be the shake blocks */
3169 do_blocka(fio, &dumb, bRead, file_version);
3182 close_symtab(&(mtop->symtab));
3186 /* If TopOnlyOK is TRUE then we can read even future versions
3187 * of tpx files, provided the file_generation hasn't changed.
3188 * If it is FALSE, we need the inputrecord too, and bail out
3189 * if the file is newer than the program.
3191 * The version and generation if the topology (see top of this file)
3192 * are returned in the two last arguments.
3194 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
3196 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
3197 gmx_bool TopOnlyOK, int *file_version,
3198 int *file_generation)
3201 char file_tag[STRLEN];
3208 /* XDR binary topology file */
3209 precision = sizeof(real);
3212 gmx_fio_do_string(fio, buf);
3213 if (std::strncmp(buf, "VERSION", 7))
3215 gmx_fatal(FARGS, "Can not read file %s,\n"
3216 " this file is from a GROMACS version which is older than 2.0\n"
3217 " Make a new one with grompp or use a gro or pdb file, if possible",
3218 gmx_fio_getname(fio));
3220 gmx_fio_do_int(fio, precision);
3221 bDouble = (precision == sizeof(double));
3222 if ((precision != sizeof(float)) && !bDouble)
3224 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3225 "instead of %d or %d",
3226 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3228 gmx_fio_setprecision(fio, bDouble);
3229 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3230 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3234 gmx_fio_write_string(fio, GromacsVersion());
3235 bDouble = (precision == sizeof(double));
3236 gmx_fio_setprecision(fio, bDouble);
3237 gmx_fio_do_int(fio, precision);
3239 sprintf(file_tag, "%s", tpx_tag);
3240 fgen = tpx_generation;
3243 /* Check versions! */
3244 gmx_fio_do_int(fio, fver);
3246 /* This is for backward compatibility with development versions 77-79
3247 * where the tag was, mistakenly, placed before the generation,
3248 * which would cause a segv instead of a proper error message
3249 * when reading the topology only from tpx with <77 code.
3251 if (fver >= 77 && fver <= 79)
3253 gmx_fio_do_string(fio, file_tag);
3258 gmx_fio_do_int(fio, fgen);
3267 gmx_fio_do_string(fio, file_tag);
3273 /* Versions before 77 don't have the tag, set it to release */
3274 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3277 if (std::strcmp(file_tag, tpx_tag) != 0)
3279 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3282 /* We only support reading tpx files with the same tag as the code
3283 * or tpx files with the release tag and with lower version number.
3285 if (std::strcmp(file_tag, TPX_TAG_RELEASE) != 0 && fver < tpx_version)
3287 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3288 gmx_fio_getname(fio), fver, file_tag,
3289 tpx_version, tpx_tag);
3294 if (file_version != NULL)
3296 *file_version = fver;
3298 if (file_generation != NULL)
3300 *file_generation = fgen;
3304 if ((fver <= tpx_incompatible_version) ||
3305 ((fver > tpx_version) && !TopOnlyOK) ||
3306 (fgen > tpx_generation) ||
3307 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3309 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3310 gmx_fio_getname(fio), fver, tpx_version);
3313 gmx_fio_do_int(fio, tpx->natoms);
3316 gmx_fio_do_int(fio, tpx->ngtc);
3324 gmx_fio_do_int(fio, idum);
3325 gmx_fio_do_real(fio, rdum);
3327 /*a better decision will eventually (5.0 or later) need to be made
3328 on how to treat the alchemical state of the system, which can now
3329 vary through a simulation, and cannot be completely described
3330 though a single lambda variable, or even a single state
3331 index. Eventually, should probably be a vector. MRS*/
3334 gmx_fio_do_int(fio, tpx->fep_state);
3336 gmx_fio_do_real(fio, tpx->lambda);
3337 gmx_fio_do_int(fio, tpx->bIr);
3338 gmx_fio_do_int(fio, tpx->bTop);
3339 gmx_fio_do_int(fio, tpx->bX);
3340 gmx_fio_do_int(fio, tpx->bV);
3341 gmx_fio_do_int(fio, tpx->bF);
3342 gmx_fio_do_int(fio, tpx->bBox);
3344 if ((fgen > tpx_generation))
3346 /* This can only happen if TopOnlyOK=TRUE */
3351 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3352 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3353 gmx_bool bXVallocated)
3359 int file_version, file_generation;
3362 gmx_bool bPeriodicMols;
3366 tpx.natoms = state->natoms;
3367 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3368 tpx.fep_state = state->fep_state;
3369 tpx.lambda = state->lambda[efptFEP];
3370 tpx.bIr = (ir != NULL);
3371 tpx.bTop = (mtop != NULL);
3372 tpx.bX = (state->x != NULL);
3373 tpx.bV = (state->v != NULL);
3374 tpx.bF = (f != NULL);
3378 TopOnlyOK = (ir == NULL);
3380 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3385 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3386 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3391 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3392 state->natoms = tpx.natoms;
3393 state->nalloc = tpx.natoms;
3399 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3403 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3405 do_test(fio, tpx.bBox, state->box);
3408 gmx_fio_ndo_rvec(fio, state->box, DIM);
3409 if (file_version >= 51)
3411 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3415 /* We initialize box_rel after reading the inputrec */
3416 clear_mat(state->box_rel);
3418 if (file_version >= 28)
3420 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3421 if (file_version < 56)
3424 gmx_fio_ndo_rvec(fio, mdum, DIM);
3429 if (state->ngtc > 0 && file_version >= 28)
3432 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3433 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3434 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3435 snew(dumv, state->ngtc);
3436 if (file_version < 69)
3438 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3440 /* These used to be the Berendsen tcoupl_lambda's */
3441 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3445 /* Prior to tpx version 26, the inputrec was here.
3446 * I moved it to enable partial forward-compatibility
3447 * for analysis/viewer programs.
3449 if (file_version < 26)
3451 do_test(fio, tpx.bIr, ir);
3456 do_inputrec(fio, ir, bRead, file_version,
3457 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3460 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3465 do_inputrec(fio, &dum_ir, bRead, file_version,
3466 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3469 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3471 done_inputrec(&dum_ir);
3477 do_test(fio, tpx.bTop, mtop);
3482 do_mtop(fio, mtop, bRead, file_version);
3486 do_mtop(fio, &dum_top, bRead, file_version);
3487 done_mtop(&dum_top, TRUE);
3490 do_test(fio, tpx.bX, state->x);
3495 state->flags |= (1<<estX);
3497 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3500 do_test(fio, tpx.bV, state->v);
3505 state->flags |= (1<<estV);
3507 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3510 do_test(fio, tpx.bF, f);
3513 gmx_fio_ndo_rvec(fio, f, state->natoms);
3516 /* Starting with tpx version 26, we have the inputrec
3517 * at the end of the file, so we can ignore it
3518 * if the file is never than the software (but still the
3519 * same generation - see comments at the top of this file.
3524 bPeriodicMols = FALSE;
3525 if (file_version >= 26)
3527 do_test(fio, tpx.bIr, ir);
3530 if (file_version >= 53)
3532 /* Removed the pbc info from do_inputrec, since we always want it */
3536 bPeriodicMols = ir->bPeriodicMols;
3538 gmx_fio_do_int(fio, ePBC);
3539 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3541 if (file_generation <= tpx_generation && ir)
3543 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3546 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3548 if (file_version < 51)
3550 set_box_rel(ir, state);
3552 if (file_version < 53)
3555 bPeriodicMols = ir->bPeriodicMols;
3558 if (bRead && ir && file_version >= 53)
3560 /* We need to do this after do_inputrec, since that initializes ir */
3562 ir->bPeriodicMols = bPeriodicMols;
3571 if (state->ngtc == 0)
3573 /* Reading old version without tcoupl state data: set it */
3574 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3576 if (tpx.bTop && mtop)
3578 if (file_version < 57)
3580 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3582 ir->eDisre = edrSimple;
3586 ir->eDisre = edrNone;
3589 set_disres_npair(mtop);
3593 if (tpx.bTop && mtop)
3595 gmx_mtop_finalize(mtop);
3598 if (file_version >= 57)
3602 env = getenv("GMX_NOCHARGEGROUPS");
3605 sscanf(env, "%d", &ienv);
3606 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3611 "Will make single atomic charge groups in non-solvent%s\n",
3612 ienv > 1 ? " and solvent" : "");
3613 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3615 fprintf(stderr, "\n");
3623 static t_fileio *open_tpx(const char *fn, const char *mode)
3625 return gmx_fio_open(fn, mode);
3628 static void close_tpx(t_fileio *fio)
3633 /************************************************************
3635 * The following routines are the exported ones
3637 ************************************************************/
3639 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3640 int *file_version, int *file_generation)
3644 fio = open_tpx(fn, "r");
3645 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3649 void write_tpx_state(const char *fn,
3650 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3654 fio = open_tpx(fn, "w");
3655 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3659 void read_tpx_state(const char *fn,
3660 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3664 fio = open_tpx(fn, "r");
3665 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3669 int read_tpx(const char *fn,
3670 t_inputrec *ir, matrix box, int *natoms,
3671 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3679 fio = open_tpx(fn, "r");
3680 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3682 *natoms = state.natoms;
3685 copy_mat(state.box, box);
3694 int read_tpx_top(const char *fn,
3695 t_inputrec *ir, matrix box, int *natoms,
3696 rvec *x, rvec *v, rvec *f, t_topology *top)
3701 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3703 *top = gmx_mtop_t_to_t_topology(&mtop);
3708 gmx_bool fn2bTPX(const char *file)
3710 return (efTPR == fn2ftp(file));
3713 static void done_gmx_groups_t(gmx_groups_t *g)
3717 for (i = 0; (i < egcNR); i++)
3719 if (NULL != g->grps[i].nm_ind)
3721 sfree(g->grps[i].nm_ind);
3722 g->grps[i].nm_ind = NULL;
3724 if (NULL != g->grpnr[i])
3730 /* The contents of this array is in symtab, don't free it here */
3734 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3735 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3738 int natoms, i, version, generation;
3739 gmx_bool bTop, bXNULL = FALSE;
3743 bTop = fn2bTPX(infile);
3747 read_tpxheader(infile, &header, TRUE, &version, &generation);
3750 snew(*x, header.natoms);
3754 snew(*v, header.natoms);
3757 *ePBC = read_tpx(infile, NULL, box, &natoms,
3758 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3759 *top = gmx_mtop_t_to_t_topology(mtop);
3760 /* In this case we need to throw away the group data too */
3761 done_gmx_groups_t(&mtop->groups);
3763 std::strcpy(title, *top->name);
3764 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3768 get_stx_coordnum(infile, &natoms);
3769 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3780 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3788 aps = gmx_atomprop_init();
3789 for (i = 0; (i < natoms); i++)
3791 if (!gmx_atomprop_query(aps, epropMass,
3792 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3793 *top->atoms.atomname[i],
3794 &(top->atoms.atom[i].m)))
3798 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3799 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3800 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3801 *top->atoms.atomname[i]);
3805 gmx_atomprop_destroy(aps);
3807 top->idef.ntypes = -1;