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! */
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/fileio/filenm.h"
48 #include "gromacs/fileio/gmxfio.h"
49 #include "gromacs/legacyheaders/copyrite.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gromacs/legacyheaders/txtdump.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/topology/atomprop.h"
55 #include "gromacs/topology/block.h"
56 #include "gromacs/topology/mtop_util.h"
57 #include "gromacs/topology/symtab.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/utility/cstringutil.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
64 #define TPX_TAG_RELEASE "release"
66 /*! \brief Tag string for the file format written to run input files
67 * written by this version of the code.
69 * Change this if you want to change the run input format in a feature
70 * branch. This ensures that there will not be different run input
71 * formats around which cannot be distinguished, while not causing
72 * problems rebasing the feature branch onto upstream changes. When
73 * merging with mainstream GROMACS, set this tag string back to
74 * TPX_TAG_RELEASE, and instead add an element to tpxv and set
75 * tpx_version to that.
77 static const char *tpx_tag = TPX_TAG_RELEASE;
79 /*! \brief Enum of values that describe the contents of a tpr file
80 * whose format matches a version number
82 * The enum helps the code be more self-documenting and ensure merges
83 * do not silently resolve when two patches make the same bump. When
84 * adding new functionality, add a new element to the end of this
85 * enumeration, change the definition of tpx_version, and write code
86 * below that does the right thing according to the value of
89 tpxv_ComputationalElectrophysiology = 96, /**< support for ion/water position swaps (computational electrophysiology) */
90 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to gmx_int64_t */
91 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
92 tpxv_InteractiveMolecularDynamics, /**< interactive molecular dynamics (IMD) */
93 tpxv_RemoveObsoleteParameters1, /**< remove optimize_fft, dihre_fc, nstcheckpoint */
94 tpxv_PullCoordTypeGeom, /**< add pull type and geometry per group and flat-bottom */
95 tpxv_PullGeomDirRel, /**< add pull geometry direction-relative */
96 tpxv_IntermolecularBondeds /**< permit inter-molecular bonded interactions in the topology */
99 /*! \brief Version number of the file format written to run input
100 * files by this version of the code.
102 * The tpx_version number should be increased whenever the file format
103 * in the main development branch changes, generally to the highest
104 * value present in tpxv. Backward compatibility for reading old run
105 * input files is maintained by checking this version number against
106 * that of the file and then using the correct code path.
108 * When developing a feature branch that needs to change the run input
109 * file format, change tpx_tag instead. */
110 static const int tpx_version = tpxv_IntermolecularBondeds;
113 /* This number should only be increased when you edit the TOPOLOGY section
114 * or the HEADER of the tpx format.
115 * This way we can maintain forward compatibility too for all analysis tools
116 * and/or external programs that only need to know the atom/residue names,
117 * charges, and bond connectivity.
119 * It first appeared in tpx version 26, when I also moved the inputrecord
120 * to the end of the tpx file, so we can just skip it if we only
123 * In particular, it must be increased when adding new elements to
124 * ftupd, so that old code can read new .tpr files.
126 static const int tpx_generation = 26;
128 /* This number should be the most recent backwards incompatible version
129 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
131 static const int tpx_incompatible_version = 9;
135 /* Struct used to maintain tpx compatibility when function types are added */
137 int fvnr; /* file version number in which the function type first appeared */
138 int ftype; /* function type */
142 * The entries should be ordered in:
143 * 1. ascending file version number
144 * 2. ascending function type number
146 /*static const t_ftupd ftupd[] = {
147 { 20, F_CUBICBONDS },
151 { 22, F_DISRESVIOL },
157 { 26, F_DIHRESVIOL },
158 { 30, F_CROSS_BOND_BONDS },
159 { 30, F_CROSS_BOND_ANGLES },
160 { 30, F_UREY_BRADLEY },
161 { 30, F_POLARIZATION },
165 * The entries should be ordered in:
166 * 1. ascending function type number
167 * 2. ascending file version number
169 /* question; what is the purpose of the commented code above? */
170 static const t_ftupd ftupd[] = {
171 { 20, F_CUBICBONDS },
176 { 43, F_TABBONDSNC },
177 { 70, F_RESTRBONDS },
178 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
179 { 76, F_LINEAR_ANGLES },
180 { 30, F_CROSS_BOND_BONDS },
181 { 30, F_CROSS_BOND_ANGLES },
182 { 30, F_UREY_BRADLEY },
183 { 34, F_QUARTIC_ANGLES },
185 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
186 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
195 { 72, F_NPSOLVATION },
197 { 41, F_LJC_PAIRS_NB },
200 { 32, F_COUL_RECIP },
203 { 30, F_POLARIZATION },
206 { 22, F_DISRESVIOL },
210 { 26, F_DIHRESVIOL },
215 { 46, F_ECONSERVED },
216 { 69, F_VTEMP_NOLONGERUSED},
218 { 54, F_DVDL_CONSTR },
219 { 76, F_ANHARM_POL },
222 { 79, F_DVDL_BONDED, },
223 { 79, F_DVDL_RESTRAINT },
224 { 79, F_DVDL_TEMPERATURE },
226 #define NFTUPD asize(ftupd)
228 /* Needed for backward compatibility */
231 /**************************************************************
233 * Now the higer level routines that do io of the structures and arrays
235 **************************************************************/
236 static void do_pullgrp_tpx_pre95(t_fileio *fio,
245 gmx_fio_do_int(fio, pgrp->nat);
248 snew(pgrp->ind, pgrp->nat);
250 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
251 gmx_fio_do_int(fio, pgrp->nweight);
254 snew(pgrp->weight, pgrp->nweight);
256 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
257 gmx_fio_do_int(fio, pgrp->pbcatom);
258 gmx_fio_do_rvec(fio, pcrd->vec);
259 clear_rvec(pcrd->origin);
260 gmx_fio_do_rvec(fio, tmp);
262 gmx_fio_do_real(fio, pcrd->rate);
263 gmx_fio_do_real(fio, pcrd->k);
264 if (file_version >= 56)
266 gmx_fio_do_real(fio, pcrd->kB);
274 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)
298 gmx_fio_do_int(fio, pcrd->group[0]);
299 gmx_fio_do_int(fio, pcrd->group[1]);
300 if (file_version >= tpxv_PullCoordTypeGeom)
302 gmx_fio_do_int(fio, pcrd->eType);
303 gmx_fio_do_int(fio, pcrd->eGeom);
304 if (pcrd->eGeom == epullgDIRRELATIVE)
306 gmx_fio_do_int(fio, pcrd->group[2]);
307 gmx_fio_do_int(fio, pcrd->group[3]);
309 gmx_fio_do_ivec(fio, pcrd->dim);
313 pcrd->eType = ePullOld;
314 pcrd->eGeom = eGeomOld;
315 copy_ivec(dimOld, pcrd->dim);
317 gmx_fio_do_rvec(fio, pcrd->origin);
318 gmx_fio_do_rvec(fio, pcrd->vec);
319 if (file_version >= tpxv_PullCoordTypeGeom)
321 gmx_fio_do_gmx_bool(fio, pcrd->bStart);
325 /* This parameter is only printed, but not actually used by mdrun */
326 pcrd->bStart = FALSE;
328 gmx_fio_do_real(fio, pcrd->init);
329 gmx_fio_do_real(fio, pcrd->rate);
330 gmx_fio_do_real(fio, pcrd->k);
331 gmx_fio_do_real(fio, pcrd->kB);
334 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
336 /* i is used in the ndo_double macro*/
340 int n_lambda = fepvals->n_lambda;
342 /* reset the lambda calculation window */
343 fepvals->lambda_start_n = 0;
344 fepvals->lambda_stop_n = n_lambda;
345 if (file_version >= 79)
351 snew(expand->init_lambda_weights, n_lambda);
353 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
354 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
357 gmx_fio_do_int(fio, expand->nstexpanded);
358 gmx_fio_do_int(fio, expand->elmcmove);
359 gmx_fio_do_int(fio, expand->elamstats);
360 gmx_fio_do_int(fio, expand->lmc_repeats);
361 gmx_fio_do_int(fio, expand->gibbsdeltalam);
362 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
363 gmx_fio_do_int(fio, expand->lmc_seed);
364 gmx_fio_do_real(fio, expand->mc_temp);
365 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
366 gmx_fio_do_int(fio, expand->nstTij);
367 gmx_fio_do_int(fio, expand->minvarmin);
368 gmx_fio_do_int(fio, expand->c_range);
369 gmx_fio_do_real(fio, expand->wl_scale);
370 gmx_fio_do_real(fio, expand->wl_ratio);
371 gmx_fio_do_real(fio, expand->init_wl_delta);
372 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
373 gmx_fio_do_int(fio, expand->elmceq);
374 gmx_fio_do_int(fio, expand->equil_steps);
375 gmx_fio_do_int(fio, expand->equil_samples);
376 gmx_fio_do_int(fio, expand->equil_n_at_lam);
377 gmx_fio_do_real(fio, expand->equil_wl_delta);
378 gmx_fio_do_real(fio, expand->equil_ratio);
382 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
385 if (file_version >= 79)
387 gmx_fio_do_int(fio, simtemp->eSimTempScale);
388 gmx_fio_do_real(fio, simtemp->simtemp_high);
389 gmx_fio_do_real(fio, simtemp->simtemp_low);
394 snew(simtemp->temperatures, n_lambda);
396 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
401 static void do_imd(t_fileio *fio, t_IMD *imd, gmx_bool bRead)
403 gmx_fio_do_int(fio, imd->nat);
406 snew(imd->ind, imd->nat);
408 gmx_fio_ndo_int(fio, imd->ind, imd->nat);
411 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
413 /* i is defined in the ndo_double macro; use g to iterate. */
418 /* free energy values */
420 if (file_version >= 79)
422 gmx_fio_do_int(fio, fepvals->init_fep_state);
423 gmx_fio_do_double(fio, fepvals->init_lambda);
424 gmx_fio_do_double(fio, fepvals->delta_lambda);
426 else if (file_version >= 59)
428 gmx_fio_do_double(fio, fepvals->init_lambda);
429 gmx_fio_do_double(fio, fepvals->delta_lambda);
433 gmx_fio_do_real(fio, rdum);
434 fepvals->init_lambda = rdum;
435 gmx_fio_do_real(fio, rdum);
436 fepvals->delta_lambda = rdum;
438 if (file_version >= 79)
440 gmx_fio_do_int(fio, fepvals->n_lambda);
443 snew(fepvals->all_lambda, efptNR);
445 for (g = 0; g < efptNR; g++)
447 if (fepvals->n_lambda > 0)
451 snew(fepvals->all_lambda[g], fepvals->n_lambda);
453 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
454 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
456 else if (fepvals->init_lambda >= 0)
458 fepvals->separate_dvdl[efptFEP] = TRUE;
462 else if (file_version >= 64)
464 gmx_fio_do_int(fio, fepvals->n_lambda);
469 snew(fepvals->all_lambda, efptNR);
470 /* still allocate the all_lambda array's contents. */
471 for (g = 0; g < efptNR; g++)
473 if (fepvals->n_lambda > 0)
475 snew(fepvals->all_lambda[g], fepvals->n_lambda);
479 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
481 if (fepvals->init_lambda >= 0)
485 fepvals->separate_dvdl[efptFEP] = TRUE;
489 /* copy the contents of the efptFEP lambda component to all
490 the other components */
491 for (g = 0; g < efptNR; g++)
493 for (h = 0; h < fepvals->n_lambda; h++)
497 fepvals->all_lambda[g][h] =
498 fepvals->all_lambda[efptFEP][h];
507 fepvals->n_lambda = 0;
508 fepvals->all_lambda = NULL;
509 if (fepvals->init_lambda >= 0)
511 fepvals->separate_dvdl[efptFEP] = TRUE;
514 if (file_version >= 13)
516 gmx_fio_do_real(fio, fepvals->sc_alpha);
520 fepvals->sc_alpha = 0;
522 if (file_version >= 38)
524 gmx_fio_do_int(fio, fepvals->sc_power);
528 fepvals->sc_power = 2;
530 if (file_version >= 79)
532 gmx_fio_do_real(fio, fepvals->sc_r_power);
536 fepvals->sc_r_power = 6.0;
538 if (file_version >= 15)
540 gmx_fio_do_real(fio, fepvals->sc_sigma);
544 fepvals->sc_sigma = 0.3;
548 if (file_version >= 71)
550 fepvals->sc_sigma_min = fepvals->sc_sigma;
554 fepvals->sc_sigma_min = 0;
557 if (file_version >= 79)
559 gmx_fio_do_int(fio, fepvals->bScCoul);
563 fepvals->bScCoul = TRUE;
565 if (file_version >= 64)
567 gmx_fio_do_int(fio, fepvals->nstdhdl);
571 fepvals->nstdhdl = 1;
574 if (file_version >= 73)
576 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
577 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
581 fepvals->separate_dhdl_file = esepdhdlfileYES;
582 fepvals->dhdl_derivatives = edhdlderivativesYES;
584 if (file_version >= 71)
586 gmx_fio_do_int(fio, fepvals->dh_hist_size);
587 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
591 fepvals->dh_hist_size = 0;
592 fepvals->dh_hist_spacing = 0.1;
594 if (file_version >= 79)
596 gmx_fio_do_int(fio, fepvals->edHdLPrintEnergy);
600 fepvals->edHdLPrintEnergy = edHdLPrintEnergyNO;
603 /* handle lambda_neighbors */
604 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
606 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
607 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
608 (fepvals->init_lambda < 0) )
610 fepvals->lambda_start_n = (fepvals->init_fep_state -
611 fepvals->lambda_neighbors);
612 fepvals->lambda_stop_n = (fepvals->init_fep_state +
613 fepvals->lambda_neighbors + 1);
614 if (fepvals->lambda_start_n < 0)
616 fepvals->lambda_start_n = 0;;
618 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
620 fepvals->lambda_stop_n = fepvals->n_lambda;
625 fepvals->lambda_start_n = 0;
626 fepvals->lambda_stop_n = fepvals->n_lambda;
631 fepvals->lambda_start_n = 0;
632 fepvals->lambda_stop_n = fepvals->n_lambda;
636 static void do_pull(t_fileio *fio, pull_params_t *pull, gmx_bool bRead,
637 int file_version, int ePullOld)
643 if (file_version >= 95)
645 gmx_fio_do_int(fio, pull->ngroup);
647 gmx_fio_do_int(fio, pull->ncoord);
648 if (file_version < 95)
650 pull->ngroup = pull->ncoord + 1;
652 if (file_version < tpxv_PullCoordTypeGeom)
656 gmx_fio_do_int(fio, eGeomOld);
657 gmx_fio_do_ivec(fio, dimOld);
658 /* The inner cylinder radius, now removed */
659 gmx_fio_do_real(fio, dum);
661 gmx_fio_do_real(fio, pull->cylinder_r);
662 gmx_fio_do_real(fio, pull->constr_tol);
663 if (file_version >= 95)
665 gmx_fio_do_int(fio, pull->bPrintCOM1);
666 /* With file_version < 95 this value is set below */
668 if (file_version >= tpxv_PullCoordTypeGeom)
670 gmx_fio_do_int(fio, pull->bPrintCOM2);
671 gmx_fio_do_int(fio, pull->bPrintRefValue);
672 gmx_fio_do_int(fio, pull->bPrintComp);
676 pull->bPrintCOM2 = FALSE;
677 pull->bPrintRefValue = FALSE;
678 pull->bPrintComp = TRUE;
680 gmx_fio_do_int(fio, pull->nstxout);
681 gmx_fio_do_int(fio, pull->nstfout);
684 snew(pull->group, pull->ngroup);
685 snew(pull->coord, pull->ncoord);
687 if (file_version < 95)
689 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
690 if (eGeomOld == epullgDIRPBC)
692 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
694 if (eGeomOld > epullgDIRPBC)
699 for (g = 0; g < pull->ngroup; g++)
701 /* We read and ignore a pull coordinate for group 0 */
702 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[max(g-1, 0)],
703 bRead, file_version);
706 pull->coord[g-1].group[0] = 0;
707 pull->coord[g-1].group[1] = g;
711 pull->bPrintCOM1 = (pull->group[0].nat > 0);
715 for (g = 0; g < pull->ngroup; g++)
717 do_pull_group(fio, &pull->group[g], bRead);
719 for (g = 0; g < pull->ncoord; g++)
721 do_pull_coord(fio, &pull->coord[g],
722 file_version, ePullOld, eGeomOld, dimOld);
728 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
732 gmx_fio_do_int(fio, rotg->eType);
733 gmx_fio_do_int(fio, rotg->bMassW);
734 gmx_fio_do_int(fio, rotg->nat);
737 snew(rotg->ind, rotg->nat);
739 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
742 snew(rotg->x_ref, rotg->nat);
744 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
745 gmx_fio_do_rvec(fio, rotg->vec);
746 gmx_fio_do_rvec(fio, rotg->pivot);
747 gmx_fio_do_real(fio, rotg->rate);
748 gmx_fio_do_real(fio, rotg->k);
749 gmx_fio_do_real(fio, rotg->slab_dist);
750 gmx_fio_do_real(fio, rotg->min_gaussian);
751 gmx_fio_do_real(fio, rotg->eps);
752 gmx_fio_do_int(fio, rotg->eFittype);
753 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
754 gmx_fio_do_real(fio, rotg->PotAngle_step);
757 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
761 gmx_fio_do_int(fio, rot->ngrp);
762 gmx_fio_do_int(fio, rot->nstrout);
763 gmx_fio_do_int(fio, rot->nstsout);
766 snew(rot->grp, rot->ngrp);
768 for (g = 0; g < rot->ngrp; g++)
770 do_rotgrp(fio, &rot->grp[g], bRead);
775 static void do_swapcoords(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead)
780 gmx_fio_do_int(fio, swap->nat);
781 gmx_fio_do_int(fio, swap->nat_sol);
782 for (j = 0; j < 2; j++)
784 gmx_fio_do_int(fio, swap->nat_split[j]);
785 gmx_fio_do_int(fio, swap->massw_split[j]);
787 gmx_fio_do_int(fio, swap->nstswap);
788 gmx_fio_do_int(fio, swap->nAverage);
789 gmx_fio_do_real(fio, swap->threshold);
790 gmx_fio_do_real(fio, swap->cyl0r);
791 gmx_fio_do_real(fio, swap->cyl0u);
792 gmx_fio_do_real(fio, swap->cyl0l);
793 gmx_fio_do_real(fio, swap->cyl1r);
794 gmx_fio_do_real(fio, swap->cyl1u);
795 gmx_fio_do_real(fio, swap->cyl1l);
799 snew(swap->ind, swap->nat);
800 snew(swap->ind_sol, swap->nat_sol);
801 for (j = 0; j < 2; j++)
803 snew(swap->ind_split[j], swap->nat_split[j]);
807 gmx_fio_ndo_int(fio, swap->ind, swap->nat);
808 gmx_fio_ndo_int(fio, swap->ind_sol, swap->nat_sol);
809 for (j = 0; j < 2; j++)
811 gmx_fio_ndo_int(fio, swap->ind_split[j], swap->nat_split[j]);
814 for (j = 0; j < eCompNR; j++)
816 gmx_fio_do_int(fio, swap->nanions[j]);
817 gmx_fio_do_int(fio, swap->ncations[j]);
823 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
824 int file_version, real *fudgeQQ)
826 int i, j, k, *tmp, idum = 0;
829 gmx_bool bSimAnn, bdum = 0;
830 real zerotemptime, finish_t, init_temp, finish_temp;
832 if (file_version != tpx_version)
834 /* Give a warning about features that are not accessible */
835 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
836 file_version, tpx_version);
844 if (file_version == 0)
849 /* Basic inputrec stuff */
850 gmx_fio_do_int(fio, ir->eI);
851 if (file_version >= 62)
853 gmx_fio_do_int64(fio, ir->nsteps);
857 gmx_fio_do_int(fio, idum);
861 if (file_version > 25)
863 if (file_version >= 62)
865 gmx_fio_do_int64(fio, ir->init_step);
869 gmx_fio_do_int(fio, idum);
870 ir->init_step = idum;
878 if (file_version >= 58)
880 gmx_fio_do_int(fio, ir->simulation_part);
884 ir->simulation_part = 1;
887 if (file_version >= 67)
889 gmx_fio_do_int(fio, ir->nstcalcenergy);
893 ir->nstcalcenergy = 1;
895 if (file_version < 53)
897 /* The pbc info has been moved out of do_inputrec,
898 * since we always want it, also without reading the inputrec.
900 gmx_fio_do_int(fio, ir->ePBC);
901 if ((file_version <= 15) && (ir->ePBC == 2))
905 if (file_version >= 45)
907 gmx_fio_do_int(fio, ir->bPeriodicMols);
914 ir->bPeriodicMols = TRUE;
918 ir->bPeriodicMols = FALSE;
922 if (file_version >= 81)
924 gmx_fio_do_int(fio, ir->cutoff_scheme);
925 if (file_version < 94)
927 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
932 ir->cutoff_scheme = ecutsGROUP;
934 gmx_fio_do_int(fio, ir->ns_type);
935 gmx_fio_do_int(fio, ir->nstlist);
936 gmx_fio_do_int(fio, idum); /* used to be ndelta; not used anymore */
937 if (file_version < 41)
939 gmx_fio_do_int(fio, idum);
940 gmx_fio_do_int(fio, idum);
942 if (file_version >= 45)
944 gmx_fio_do_real(fio, ir->rtpi);
950 gmx_fio_do_int(fio, ir->nstcomm);
951 if (file_version > 34)
953 gmx_fio_do_int(fio, ir->comm_mode);
955 else if (ir->nstcomm < 0)
957 ir->comm_mode = ecmANGULAR;
961 ir->comm_mode = ecmLINEAR;
963 ir->nstcomm = abs(ir->nstcomm);
965 /* ignore nstcheckpoint */
966 if (file_version > 25 && file_version < tpxv_RemoveObsoleteParameters1)
968 gmx_fio_do_int(fio, idum);
971 gmx_fio_do_int(fio, ir->nstcgsteep);
973 if (file_version >= 30)
975 gmx_fio_do_int(fio, ir->nbfgscorr);
982 gmx_fio_do_int(fio, ir->nstlog);
983 gmx_fio_do_int(fio, ir->nstxout);
984 gmx_fio_do_int(fio, ir->nstvout);
985 gmx_fio_do_int(fio, ir->nstfout);
986 gmx_fio_do_int(fio, ir->nstenergy);
987 gmx_fio_do_int(fio, ir->nstxout_compressed);
988 if (file_version >= 59)
990 gmx_fio_do_double(fio, ir->init_t);
991 gmx_fio_do_double(fio, ir->delta_t);
995 gmx_fio_do_real(fio, rdum);
997 gmx_fio_do_real(fio, rdum);
1000 gmx_fio_do_real(fio, ir->x_compression_precision);
1001 if (file_version < 19)
1003 gmx_fio_do_int(fio, idum);
1004 gmx_fio_do_int(fio, idum);
1006 if (file_version < 18)
1008 gmx_fio_do_int(fio, idum);
1010 if (file_version >= 81)
1012 gmx_fio_do_real(fio, ir->verletbuf_tol);
1016 ir->verletbuf_tol = 0;
1018 gmx_fio_do_real(fio, ir->rlist);
1019 if (file_version >= 67)
1021 gmx_fio_do_real(fio, ir->rlistlong);
1023 if (file_version >= 82 && file_version != 90)
1025 gmx_fio_do_int(fio, ir->nstcalclr);
1029 /* Calculate at NS steps */
1030 ir->nstcalclr = ir->nstlist;
1032 gmx_fio_do_int(fio, ir->coulombtype);
1033 if (file_version < 32 && ir->coulombtype == eelRF)
1035 ir->coulombtype = eelRF_NEC;
1037 if (file_version >= 81)
1039 gmx_fio_do_int(fio, ir->coulomb_modifier);
1043 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1045 gmx_fio_do_real(fio, ir->rcoulomb_switch);
1046 gmx_fio_do_real(fio, ir->rcoulomb);
1047 gmx_fio_do_int(fio, ir->vdwtype);
1048 if (file_version >= 81)
1050 gmx_fio_do_int(fio, ir->vdw_modifier);
1054 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1056 gmx_fio_do_real(fio, ir->rvdw_switch);
1057 gmx_fio_do_real(fio, ir->rvdw);
1058 if (file_version < 67)
1060 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1062 gmx_fio_do_int(fio, ir->eDispCorr);
1063 gmx_fio_do_real(fio, ir->epsilon_r);
1064 if (file_version >= 37)
1066 gmx_fio_do_real(fio, ir->epsilon_rf);
1070 if (EEL_RF(ir->coulombtype))
1072 ir->epsilon_rf = ir->epsilon_r;
1073 ir->epsilon_r = 1.0;
1077 ir->epsilon_rf = 1.0;
1080 if (file_version >= 29)
1082 gmx_fio_do_real(fio, ir->tabext);
1089 if (file_version > 25)
1091 gmx_fio_do_int(fio, ir->gb_algorithm);
1092 gmx_fio_do_int(fio, ir->nstgbradii);
1093 gmx_fio_do_real(fio, ir->rgbradii);
1094 gmx_fio_do_real(fio, ir->gb_saltconc);
1095 gmx_fio_do_int(fio, ir->implicit_solvent);
1099 ir->gb_algorithm = egbSTILL;
1102 ir->gb_saltconc = 0;
1103 ir->implicit_solvent = eisNO;
1105 if (file_version >= 55)
1107 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1108 gmx_fio_do_real(fio, ir->gb_obc_alpha);
1109 gmx_fio_do_real(fio, ir->gb_obc_beta);
1110 gmx_fio_do_real(fio, ir->gb_obc_gamma);
1111 if (file_version >= 60)
1113 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1114 gmx_fio_do_int(fio, ir->sa_algorithm);
1118 ir->gb_dielectric_offset = 0.009;
1119 ir->sa_algorithm = esaAPPROX;
1121 gmx_fio_do_real(fio, ir->sa_surface_tension);
1123 /* Override sa_surface_tension if it is not changed in the mpd-file */
1124 if (ir->sa_surface_tension < 0)
1126 if (ir->gb_algorithm == egbSTILL)
1128 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1130 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1132 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1139 /* Better use sensible values than insane (0.0) ones... */
1140 ir->gb_epsilon_solvent = 80;
1141 ir->gb_obc_alpha = 1.0;
1142 ir->gb_obc_beta = 0.8;
1143 ir->gb_obc_gamma = 4.85;
1144 ir->sa_surface_tension = 2.092;
1148 if (file_version >= 81)
1150 gmx_fio_do_real(fio, ir->fourier_spacing);
1154 ir->fourier_spacing = 0.0;
1156 gmx_fio_do_int(fio, ir->nkx);
1157 gmx_fio_do_int(fio, ir->nky);
1158 gmx_fio_do_int(fio, ir->nkz);
1159 gmx_fio_do_int(fio, ir->pme_order);
1160 gmx_fio_do_real(fio, ir->ewald_rtol);
1162 if (file_version >= 93)
1164 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1168 ir->ewald_rtol_lj = ir->ewald_rtol;
1171 if (file_version >= 24)
1173 gmx_fio_do_int(fio, ir->ewald_geometry);
1177 ir->ewald_geometry = eewg3D;
1180 if (file_version <= 17)
1182 ir->epsilon_surface = 0;
1183 if (file_version == 17)
1185 gmx_fio_do_int(fio, idum);
1190 gmx_fio_do_real(fio, ir->epsilon_surface);
1193 /* ignore bOptFFT */
1194 if (file_version < tpxv_RemoveObsoleteParameters1)
1196 gmx_fio_do_gmx_bool(fio, bdum);
1199 if (file_version >= 93)
1201 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1203 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1204 gmx_fio_do_int(fio, ir->etc);
1205 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1206 * but the values 0 and 1 still mean no and
1207 * berendsen temperature coupling, respectively.
1209 if (file_version >= 79)
1211 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1213 if (file_version >= 71)
1215 gmx_fio_do_int(fio, ir->nsttcouple);
1219 ir->nsttcouple = ir->nstcalcenergy;
1221 if (file_version <= 15)
1223 gmx_fio_do_int(fio, idum);
1225 if (file_version <= 17)
1227 gmx_fio_do_int(fio, ir->epct);
1228 if (file_version <= 15)
1232 ir->epct = epctSURFACETENSION;
1234 gmx_fio_do_int(fio, idum);
1237 /* we have removed the NO alternative at the beginning */
1241 ir->epct = epctISOTROPIC;
1245 ir->epc = epcBERENDSEN;
1250 gmx_fio_do_int(fio, ir->epc);
1251 gmx_fio_do_int(fio, ir->epct);
1253 if (file_version >= 71)
1255 gmx_fio_do_int(fio, ir->nstpcouple);
1259 ir->nstpcouple = ir->nstcalcenergy;
1261 gmx_fio_do_real(fio, ir->tau_p);
1262 if (file_version <= 15)
1264 gmx_fio_do_rvec(fio, vdum);
1265 clear_mat(ir->ref_p);
1266 for (i = 0; i < DIM; i++)
1268 ir->ref_p[i][i] = vdum[i];
1273 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1274 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1275 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1277 if (file_version <= 15)
1279 gmx_fio_do_rvec(fio, vdum);
1280 clear_mat(ir->compress);
1281 for (i = 0; i < DIM; i++)
1283 ir->compress[i][i] = vdum[i];
1288 gmx_fio_do_rvec(fio, ir->compress[XX]);
1289 gmx_fio_do_rvec(fio, ir->compress[YY]);
1290 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1292 if (file_version >= 47)
1294 gmx_fio_do_int(fio, ir->refcoord_scaling);
1295 gmx_fio_do_rvec(fio, ir->posres_com);
1296 gmx_fio_do_rvec(fio, ir->posres_comB);
1300 ir->refcoord_scaling = erscNO;
1301 clear_rvec(ir->posres_com);
1302 clear_rvec(ir->posres_comB);
1304 if ((file_version > 25) && (file_version < 79))
1306 gmx_fio_do_int(fio, ir->andersen_seed);
1310 ir->andersen_seed = 0;
1312 if (file_version < 26)
1314 gmx_fio_do_gmx_bool(fio, bSimAnn);
1315 gmx_fio_do_real(fio, zerotemptime);
1318 if (file_version < 37)
1320 gmx_fio_do_real(fio, rdum);
1323 gmx_fio_do_real(fio, ir->shake_tol);
1324 if (file_version < 54)
1326 gmx_fio_do_real(fio, *fudgeQQ);
1329 gmx_fio_do_int(fio, ir->efep);
1330 if (file_version <= 14 && ir->efep != efepNO)
1334 do_fepvals(fio, ir->fepvals, bRead, file_version);
1336 if (file_version >= 79)
1338 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1341 ir->bSimTemp = TRUE;
1346 ir->bSimTemp = FALSE;
1350 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1353 if (file_version >= 79)
1355 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1358 ir->bExpanded = TRUE;
1362 ir->bExpanded = FALSE;
1367 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1369 if (file_version >= 57)
1371 gmx_fio_do_int(fio, ir->eDisre);
1373 gmx_fio_do_int(fio, ir->eDisreWeighting);
1374 if (file_version < 22)
1376 if (ir->eDisreWeighting == 0)
1378 ir->eDisreWeighting = edrwEqual;
1382 ir->eDisreWeighting = edrwConservative;
1385 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1386 gmx_fio_do_real(fio, ir->dr_fc);
1387 gmx_fio_do_real(fio, ir->dr_tau);
1388 gmx_fio_do_int(fio, ir->nstdisreout);
1389 if (file_version >= 22)
1391 gmx_fio_do_real(fio, ir->orires_fc);
1392 gmx_fio_do_real(fio, ir->orires_tau);
1393 gmx_fio_do_int(fio, ir->nstorireout);
1399 ir->nstorireout = 0;
1402 /* ignore dihre_fc */
1403 if (file_version >= 26 && file_version < 79)
1405 gmx_fio_do_real(fio, rdum);
1406 if (file_version < 56)
1408 gmx_fio_do_real(fio, rdum);
1409 gmx_fio_do_int(fio, idum);
1413 gmx_fio_do_real(fio, ir->em_stepsize);
1414 gmx_fio_do_real(fio, ir->em_tol);
1415 if (file_version >= 22)
1417 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1421 ir->bShakeSOR = TRUE;
1423 if (file_version >= 11)
1425 gmx_fio_do_int(fio, ir->niter);
1430 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1433 if (file_version >= 21)
1435 gmx_fio_do_real(fio, ir->fc_stepsize);
1439 ir->fc_stepsize = 0;
1441 gmx_fio_do_int(fio, ir->eConstrAlg);
1442 gmx_fio_do_int(fio, ir->nProjOrder);
1443 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1444 if (file_version <= 14)
1446 gmx_fio_do_int(fio, idum);
1448 if (file_version >= 26)
1450 gmx_fio_do_int(fio, ir->nLincsIter);
1455 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1458 if (file_version < 33)
1460 gmx_fio_do_real(fio, bd_temp);
1462 gmx_fio_do_real(fio, ir->bd_fric);
1463 if (file_version >= tpxv_Use64BitRandomSeed)
1465 gmx_fio_do_int64(fio, ir->ld_seed);
1469 gmx_fio_do_int(fio, idum);
1472 if (file_version >= 33)
1474 for (i = 0; i < DIM; i++)
1476 gmx_fio_do_rvec(fio, ir->deform[i]);
1481 for (i = 0; i < DIM; i++)
1483 clear_rvec(ir->deform[i]);
1486 if (file_version >= 14)
1488 gmx_fio_do_real(fio, ir->cos_accel);
1494 gmx_fio_do_int(fio, ir->userint1);
1495 gmx_fio_do_int(fio, ir->userint2);
1496 gmx_fio_do_int(fio, ir->userint3);
1497 gmx_fio_do_int(fio, ir->userint4);
1498 gmx_fio_do_real(fio, ir->userreal1);
1499 gmx_fio_do_real(fio, ir->userreal2);
1500 gmx_fio_do_real(fio, ir->userreal3);
1501 gmx_fio_do_real(fio, ir->userreal4);
1504 if (file_version >= 77)
1506 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1511 snew(ir->adress, 1);
1513 gmx_fio_do_int(fio, ir->adress->type);
1514 gmx_fio_do_real(fio, ir->adress->const_wf);
1515 gmx_fio_do_real(fio, ir->adress->ex_width);
1516 gmx_fio_do_real(fio, ir->adress->hy_width);
1517 gmx_fio_do_int(fio, ir->adress->icor);
1518 gmx_fio_do_int(fio, ir->adress->site);
1519 gmx_fio_do_rvec(fio, ir->adress->refs);
1520 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1521 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1522 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1523 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1527 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1529 if (ir->adress->n_tf_grps > 0)
1531 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1535 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1537 if (ir->adress->n_energy_grps > 0)
1539 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1545 ir->bAdress = FALSE;
1549 if (file_version >= 48)
1553 if (file_version >= tpxv_PullCoordTypeGeom)
1555 gmx_fio_do_gmx_bool(fio, ir->bPull);
1559 gmx_fio_do_int(fio, ePullOld);
1560 ir->bPull = (ePullOld > 0);
1561 /* We removed the first ePull=ePullNo for the enum */
1570 do_pull(fio, ir->pull, bRead, file_version, ePullOld);
1578 /* Enforced rotation */
1579 if (file_version >= 74)
1581 gmx_fio_do_int(fio, ir->bRot);
1582 if (ir->bRot == TRUE)
1588 do_rot(fio, ir->rot, bRead);
1596 /* Interactive molecular dynamics */
1597 if (file_version >= tpxv_InteractiveMolecularDynamics)
1599 gmx_fio_do_int(fio, ir->bIMD);
1600 if (TRUE == ir->bIMD)
1606 do_imd(fio, ir->imd, bRead);
1611 /* We don't support IMD sessions for old .tpr files */
1616 gmx_fio_do_int(fio, ir->opts.ngtc);
1617 if (file_version >= 69)
1619 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1623 ir->opts.nhchainlength = 1;
1625 gmx_fio_do_int(fio, ir->opts.ngacc);
1626 gmx_fio_do_int(fio, ir->opts.ngfrz);
1627 gmx_fio_do_int(fio, ir->opts.ngener);
1631 snew(ir->opts.nrdf, ir->opts.ngtc);
1632 snew(ir->opts.ref_t, ir->opts.ngtc);
1633 snew(ir->opts.annealing, ir->opts.ngtc);
1634 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1635 snew(ir->opts.anneal_time, ir->opts.ngtc);
1636 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1637 snew(ir->opts.tau_t, ir->opts.ngtc);
1638 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1639 snew(ir->opts.acc, ir->opts.ngacc);
1640 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1642 if (ir->opts.ngtc > 0)
1644 if (bRead && file_version < 13)
1646 snew(tmp, ir->opts.ngtc);
1647 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1648 for (i = 0; i < ir->opts.ngtc; i++)
1650 ir->opts.nrdf[i] = tmp[i];
1656 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1658 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1659 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1660 if (file_version < 33 && ir->eI == eiBD)
1662 for (i = 0; i < ir->opts.ngtc; i++)
1664 ir->opts.tau_t[i] = bd_temp;
1668 if (ir->opts.ngfrz > 0)
1670 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1672 if (ir->opts.ngacc > 0)
1674 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1676 if (file_version >= 12)
1678 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1679 ir->opts.ngener*ir->opts.ngener);
1682 if (bRead && file_version < 26)
1684 for (i = 0; i < ir->opts.ngtc; i++)
1688 ir->opts.annealing[i] = eannSINGLE;
1689 ir->opts.anneal_npoints[i] = 2;
1690 snew(ir->opts.anneal_time[i], 2);
1691 snew(ir->opts.anneal_temp[i], 2);
1692 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1693 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1694 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1695 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1696 ir->opts.anneal_time[i][0] = ir->init_t;
1697 ir->opts.anneal_time[i][1] = finish_t;
1698 ir->opts.anneal_temp[i][0] = init_temp;
1699 ir->opts.anneal_temp[i][1] = finish_temp;
1703 ir->opts.annealing[i] = eannNO;
1704 ir->opts.anneal_npoints[i] = 0;
1710 /* file version 26 or later */
1711 /* First read the lists with annealing and npoints for each group */
1712 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1713 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1714 for (j = 0; j < (ir->opts.ngtc); j++)
1716 k = ir->opts.anneal_npoints[j];
1719 snew(ir->opts.anneal_time[j], k);
1720 snew(ir->opts.anneal_temp[j], k);
1722 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1723 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1727 if (file_version >= 45)
1729 gmx_fio_do_int(fio, ir->nwall);
1730 gmx_fio_do_int(fio, ir->wall_type);
1731 if (file_version >= 50)
1733 gmx_fio_do_real(fio, ir->wall_r_linpot);
1737 ir->wall_r_linpot = -1;
1739 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1740 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1741 gmx_fio_do_real(fio, ir->wall_density[0]);
1742 gmx_fio_do_real(fio, ir->wall_density[1]);
1743 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1749 ir->wall_atomtype[0] = -1;
1750 ir->wall_atomtype[1] = -1;
1751 ir->wall_density[0] = 0;
1752 ir->wall_density[1] = 0;
1753 ir->wall_ewald_zfac = 3;
1755 /* Cosine stuff for electric fields */
1756 for (j = 0; (j < DIM); j++)
1758 gmx_fio_do_int(fio, ir->ex[j].n);
1759 gmx_fio_do_int(fio, ir->et[j].n);
1762 snew(ir->ex[j].a, ir->ex[j].n);
1763 snew(ir->ex[j].phi, ir->ex[j].n);
1764 snew(ir->et[j].a, ir->et[j].n);
1765 snew(ir->et[j].phi, ir->et[j].n);
1767 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1768 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1769 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1770 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1774 if (file_version >= tpxv_ComputationalElectrophysiology)
1776 gmx_fio_do_int(fio, ir->eSwapCoords);
1777 if (ir->eSwapCoords != eswapNO)
1783 do_swapcoords(fio, ir->swap, bRead);
1788 if (file_version >= 39)
1790 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1791 gmx_fio_do_int(fio, ir->QMMMscheme);
1792 gmx_fio_do_real(fio, ir->scalefactor);
1793 gmx_fio_do_int(fio, ir->opts.ngQM);
1796 snew(ir->opts.QMmethod, ir->opts.ngQM);
1797 snew(ir->opts.QMbasis, ir->opts.ngQM);
1798 snew(ir->opts.QMcharge, ir->opts.ngQM);
1799 snew(ir->opts.QMmult, ir->opts.ngQM);
1800 snew(ir->opts.bSH, ir->opts.ngQM);
1801 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1802 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1803 snew(ir->opts.SAon, ir->opts.ngQM);
1804 snew(ir->opts.SAoff, ir->opts.ngQM);
1805 snew(ir->opts.SAsteps, ir->opts.ngQM);
1806 snew(ir->opts.bOPT, ir->opts.ngQM);
1807 snew(ir->opts.bTS, ir->opts.ngQM);
1809 if (ir->opts.ngQM > 0)
1811 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1812 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1813 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1814 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1815 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1816 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1817 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1818 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1819 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1820 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1821 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1822 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1824 /* end of QMMM stuff */
1829 static void do_harm(t_fileio *fio, t_iparams *iparams)
1831 gmx_fio_do_real(fio, iparams->harmonic.rA);
1832 gmx_fio_do_real(fio, iparams->harmonic.krA);
1833 gmx_fio_do_real(fio, iparams->harmonic.rB);
1834 gmx_fio_do_real(fio, iparams->harmonic.krB);
1837 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1838 gmx_bool bRead, int file_version)
1851 do_harm(fio, iparams);
1852 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1854 /* Correct incorrect storage of parameters */
1855 iparams->pdihs.phiB = iparams->pdihs.phiA;
1856 iparams->pdihs.cpB = iparams->pdihs.cpA;
1860 gmx_fio_do_real(fio, iparams->harmonic.rA);
1861 gmx_fio_do_real(fio, iparams->harmonic.krA);
1863 case F_LINEAR_ANGLES:
1864 gmx_fio_do_real(fio, iparams->linangle.klinA);
1865 gmx_fio_do_real(fio, iparams->linangle.aA);
1866 gmx_fio_do_real(fio, iparams->linangle.klinB);
1867 gmx_fio_do_real(fio, iparams->linangle.aB);
1870 gmx_fio_do_real(fio, iparams->fene.bm);
1871 gmx_fio_do_real(fio, iparams->fene.kb);
1875 gmx_fio_do_real(fio, iparams->restraint.lowA);
1876 gmx_fio_do_real(fio, iparams->restraint.up1A);
1877 gmx_fio_do_real(fio, iparams->restraint.up2A);
1878 gmx_fio_do_real(fio, iparams->restraint.kA);
1879 gmx_fio_do_real(fio, iparams->restraint.lowB);
1880 gmx_fio_do_real(fio, iparams->restraint.up1B);
1881 gmx_fio_do_real(fio, iparams->restraint.up2B);
1882 gmx_fio_do_real(fio, iparams->restraint.kB);
1888 gmx_fio_do_real(fio, iparams->tab.kA);
1889 gmx_fio_do_int(fio, iparams->tab.table);
1890 gmx_fio_do_real(fio, iparams->tab.kB);
1892 case F_CROSS_BOND_BONDS:
1893 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1894 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1895 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1897 case F_CROSS_BOND_ANGLES:
1898 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1899 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1900 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1901 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1903 case F_UREY_BRADLEY:
1904 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1905 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1906 gmx_fio_do_real(fio, iparams->u_b.r13A);
1907 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1908 if (file_version >= 79)
1910 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1911 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1912 gmx_fio_do_real(fio, iparams->u_b.r13B);
1913 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1917 iparams->u_b.thetaB = iparams->u_b.thetaA;
1918 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1919 iparams->u_b.r13B = iparams->u_b.r13A;
1920 iparams->u_b.kUBB = iparams->u_b.kUBA;
1923 case F_QUARTIC_ANGLES:
1924 gmx_fio_do_real(fio, iparams->qangle.theta);
1925 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1928 gmx_fio_do_real(fio, iparams->bham.a);
1929 gmx_fio_do_real(fio, iparams->bham.b);
1930 gmx_fio_do_real(fio, iparams->bham.c);
1933 gmx_fio_do_real(fio, iparams->morse.b0A);
1934 gmx_fio_do_real(fio, iparams->morse.cbA);
1935 gmx_fio_do_real(fio, iparams->morse.betaA);
1936 if (file_version >= 79)
1938 gmx_fio_do_real(fio, iparams->morse.b0B);
1939 gmx_fio_do_real(fio, iparams->morse.cbB);
1940 gmx_fio_do_real(fio, iparams->morse.betaB);
1944 iparams->morse.b0B = iparams->morse.b0A;
1945 iparams->morse.cbB = iparams->morse.cbA;
1946 iparams->morse.betaB = iparams->morse.betaA;
1950 gmx_fio_do_real(fio, iparams->cubic.b0);
1951 gmx_fio_do_real(fio, iparams->cubic.kb);
1952 gmx_fio_do_real(fio, iparams->cubic.kcub);
1956 case F_POLARIZATION:
1957 gmx_fio_do_real(fio, iparams->polarize.alpha);
1960 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1961 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1962 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1965 if (file_version < 31)
1967 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1969 gmx_fio_do_real(fio, iparams->wpol.al_x);
1970 gmx_fio_do_real(fio, iparams->wpol.al_y);
1971 gmx_fio_do_real(fio, iparams->wpol.al_z);
1972 gmx_fio_do_real(fio, iparams->wpol.rOH);
1973 gmx_fio_do_real(fio, iparams->wpol.rHH);
1974 gmx_fio_do_real(fio, iparams->wpol.rOD);
1977 gmx_fio_do_real(fio, iparams->thole.a);
1978 gmx_fio_do_real(fio, iparams->thole.alpha1);
1979 gmx_fio_do_real(fio, iparams->thole.alpha2);
1980 gmx_fio_do_real(fio, iparams->thole.rfac);
1983 gmx_fio_do_real(fio, iparams->lj.c6);
1984 gmx_fio_do_real(fio, iparams->lj.c12);
1987 gmx_fio_do_real(fio, iparams->lj14.c6A);
1988 gmx_fio_do_real(fio, iparams->lj14.c12A);
1989 gmx_fio_do_real(fio, iparams->lj14.c6B);
1990 gmx_fio_do_real(fio, iparams->lj14.c12B);
1993 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1994 gmx_fio_do_real(fio, iparams->ljc14.qi);
1995 gmx_fio_do_real(fio, iparams->ljc14.qj);
1996 gmx_fio_do_real(fio, iparams->ljc14.c6);
1997 gmx_fio_do_real(fio, iparams->ljc14.c12);
1999 case F_LJC_PAIRS_NB:
2000 gmx_fio_do_real(fio, iparams->ljcnb.qi);
2001 gmx_fio_do_real(fio, iparams->ljcnb.qj);
2002 gmx_fio_do_real(fio, iparams->ljcnb.c6);
2003 gmx_fio_do_real(fio, iparams->ljcnb.c12);
2009 gmx_fio_do_real(fio, iparams->pdihs.phiA);
2010 gmx_fio_do_real(fio, iparams->pdihs.cpA);
2011 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
2013 /* Read the incorrectly stored multiplicity */
2014 gmx_fio_do_real(fio, iparams->harmonic.rB);
2015 gmx_fio_do_real(fio, iparams->harmonic.krB);
2016 iparams->pdihs.phiB = iparams->pdihs.phiA;
2017 iparams->pdihs.cpB = iparams->pdihs.cpA;
2021 gmx_fio_do_real(fio, iparams->pdihs.phiB);
2022 gmx_fio_do_real(fio, iparams->pdihs.cpB);
2023 gmx_fio_do_int(fio, iparams->pdihs.mult);
2027 gmx_fio_do_real(fio, iparams->pdihs.phiA);
2028 gmx_fio_do_real(fio, iparams->pdihs.cpA);
2031 gmx_fio_do_int(fio, iparams->disres.label);
2032 gmx_fio_do_int(fio, iparams->disres.type);
2033 gmx_fio_do_real(fio, iparams->disres.low);
2034 gmx_fio_do_real(fio, iparams->disres.up1);
2035 gmx_fio_do_real(fio, iparams->disres.up2);
2036 gmx_fio_do_real(fio, iparams->disres.kfac);
2039 gmx_fio_do_int(fio, iparams->orires.ex);
2040 gmx_fio_do_int(fio, iparams->orires.label);
2041 gmx_fio_do_int(fio, iparams->orires.power);
2042 gmx_fio_do_real(fio, iparams->orires.c);
2043 gmx_fio_do_real(fio, iparams->orires.obs);
2044 gmx_fio_do_real(fio, iparams->orires.kfac);
2047 if (file_version < 82)
2049 gmx_fio_do_int(fio, idum);
2050 gmx_fio_do_int(fio, idum);
2052 gmx_fio_do_real(fio, iparams->dihres.phiA);
2053 gmx_fio_do_real(fio, iparams->dihres.dphiA);
2054 gmx_fio_do_real(fio, iparams->dihres.kfacA);
2055 if (file_version >= 82)
2057 gmx_fio_do_real(fio, iparams->dihres.phiB);
2058 gmx_fio_do_real(fio, iparams->dihres.dphiB);
2059 gmx_fio_do_real(fio, iparams->dihres.kfacB);
2063 iparams->dihres.phiB = iparams->dihres.phiA;
2064 iparams->dihres.dphiB = iparams->dihres.dphiA;
2065 iparams->dihres.kfacB = iparams->dihres.kfacA;
2069 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
2070 gmx_fio_do_rvec(fio, iparams->posres.fcA);
2071 if (bRead && file_version < 27)
2073 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
2074 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
2078 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
2079 gmx_fio_do_rvec(fio, iparams->posres.fcB);
2083 gmx_fio_do_int(fio, iparams->fbposres.geom);
2084 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
2085 gmx_fio_do_real(fio, iparams->fbposres.r);
2086 gmx_fio_do_real(fio, iparams->fbposres.k);
2089 gmx_fio_ndo_real(fio, iparams->cbtdihs.cbtcA, NR_CBTDIHS);
2092 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2093 if (file_version >= 25)
2095 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2099 /* Fourier dihedrals are internally represented
2100 * as Ryckaert-Bellemans since those are faster to compute.
2102 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2103 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2107 gmx_fio_do_real(fio, iparams->constr.dA);
2108 gmx_fio_do_real(fio, iparams->constr.dB);
2111 gmx_fio_do_real(fio, iparams->settle.doh);
2112 gmx_fio_do_real(fio, iparams->settle.dhh);
2115 gmx_fio_do_real(fio, iparams->vsite.a);
2120 gmx_fio_do_real(fio, iparams->vsite.a);
2121 gmx_fio_do_real(fio, iparams->vsite.b);
2126 gmx_fio_do_real(fio, iparams->vsite.a);
2127 gmx_fio_do_real(fio, iparams->vsite.b);
2128 gmx_fio_do_real(fio, iparams->vsite.c);
2131 gmx_fio_do_int(fio, iparams->vsiten.n);
2132 gmx_fio_do_real(fio, iparams->vsiten.a);
2137 /* We got rid of some parameters in version 68 */
2138 if (bRead && file_version < 68)
2140 gmx_fio_do_real(fio, rdum);
2141 gmx_fio_do_real(fio, rdum);
2142 gmx_fio_do_real(fio, rdum);
2143 gmx_fio_do_real(fio, rdum);
2145 gmx_fio_do_real(fio, iparams->gb.sar);
2146 gmx_fio_do_real(fio, iparams->gb.st);
2147 gmx_fio_do_real(fio, iparams->gb.pi);
2148 gmx_fio_do_real(fio, iparams->gb.gbr);
2149 gmx_fio_do_real(fio, iparams->gb.bmlt);
2152 gmx_fio_do_int(fio, iparams->cmap.cmapA);
2153 gmx_fio_do_int(fio, iparams->cmap.cmapB);
2156 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
2157 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
2161 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version)
2165 if (file_version < 44)
2167 for (i = 0; i < MAXNODES; i++)
2169 gmx_fio_do_int(fio, idum);
2172 gmx_fio_do_int(fio, ilist->nr);
2175 snew(ilist->iatoms, ilist->nr);
2177 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2180 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2181 gmx_bool bRead, int file_version)
2186 gmx_fio_do_int(fio, ffparams->atnr);
2187 if (file_version < 57)
2189 gmx_fio_do_int(fio, idum);
2191 gmx_fio_do_int(fio, ffparams->ntypes);
2194 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2195 ffparams->atnr, ffparams->ntypes);
2199 snew(ffparams->functype, ffparams->ntypes);
2200 snew(ffparams->iparams, ffparams->ntypes);
2202 /* Read/write all the function types */
2203 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2206 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2209 if (file_version >= 66)
2211 gmx_fio_do_double(fio, ffparams->reppow);
2215 ffparams->reppow = 12.0;
2218 if (file_version >= 57)
2220 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2223 /* Check whether all these function types are supported by the code.
2224 * In practice the code is backwards compatible, which means that the
2225 * numbering may have to be altered from old numbering to new numbering
2227 for (i = 0; (i < ffparams->ntypes); i++)
2231 /* Loop over file versions */
2232 for (k = 0; (k < NFTUPD); k++)
2234 /* Compare the read file_version to the update table */
2235 if ((file_version < ftupd[k].fvnr) &&
2236 (ffparams->functype[i] >= ftupd[k].ftype))
2238 ffparams->functype[i] += 1;
2241 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2242 i, ffparams->functype[i],
2243 interaction_function[ftupd[k].ftype].longname);
2250 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2254 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2259 static void add_settle_atoms(t_ilist *ilist)
2263 /* Settle used to only store the first atom: add the other two */
2264 srenew(ilist->iatoms, 2*ilist->nr);
2265 for (i = ilist->nr/2-1; i >= 0; i--)
2267 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2268 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2269 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2270 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2272 ilist->nr = 2*ilist->nr;
2275 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2278 int i, j, renum[F_NRE];
2282 for (j = 0; (j < F_NRE); j++)
2287 for (k = 0; k < NFTUPD; k++)
2289 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2298 ilist[j].iatoms = NULL;
2302 do_ilist(fio, &ilist[j], bRead, file_version);
2303 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2305 add_settle_atoms(&ilist[j]);
2309 if (bRead && gmx_debug_at)
2310 pr_ilist(debug,0,interaction_function[j].longname,
2311 functype,&ilist[j],TRUE);
2316 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2317 gmx_bool bRead, int file_version)
2319 do_ffparams(fio, ffparams, bRead, file_version);
2321 if (file_version >= 54)
2323 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2326 do_ilists(fio, molt->ilist, bRead, file_version);
2329 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2331 int i, idum, dum_nra, *dum_a;
2333 if (file_version < 44)
2335 for (i = 0; i < MAXNODES; i++)
2337 gmx_fio_do_int(fio, idum);
2340 gmx_fio_do_int(fio, block->nr);
2341 if (file_version < 51)
2343 gmx_fio_do_int(fio, dum_nra);
2347 if ((block->nalloc_index > 0) && (NULL != block->index))
2349 sfree(block->index);
2351 block->nalloc_index = block->nr+1;
2352 snew(block->index, block->nalloc_index);
2354 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2356 if (file_version < 51 && dum_nra > 0)
2358 snew(dum_a, dum_nra);
2359 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2364 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2369 if (file_version < 44)
2371 for (i = 0; i < MAXNODES; i++)
2373 gmx_fio_do_int(fio, idum);
2376 gmx_fio_do_int(fio, block->nr);
2377 gmx_fio_do_int(fio, block->nra);
2380 block->nalloc_index = block->nr+1;
2381 snew(block->index, block->nalloc_index);
2382 block->nalloc_a = block->nra;
2383 snew(block->a, block->nalloc_a);
2385 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2386 gmx_fio_ndo_int(fio, block->a, block->nra);
2389 /* This is a primitive routine to make it possible to translate atomic numbers
2390 * to element names when reading TPR files, without making the Gromacs library
2391 * directory a dependency on mdrun (which is the case if we need elements.dat).
2394 atomicnumber_to_element(int atomicnumber)
2398 /* This does not have to be complete, so we only include elements likely
2399 * to occur in PDB files.
2401 switch (atomicnumber)
2403 case 1: p = "H"; break;
2404 case 5: p = "B"; break;
2405 case 6: p = "C"; break;
2406 case 7: p = "N"; break;
2407 case 8: p = "O"; break;
2408 case 9: p = "F"; break;
2409 case 11: p = "Na"; break;
2410 case 12: p = "Mg"; break;
2411 case 15: p = "P"; break;
2412 case 16: p = "S"; break;
2413 case 17: p = "Cl"; break;
2414 case 18: p = "Ar"; break;
2415 case 19: p = "K"; break;
2416 case 20: p = "Ca"; break;
2417 case 25: p = "Mn"; break;
2418 case 26: p = "Fe"; break;
2419 case 28: p = "Ni"; break;
2420 case 29: p = "Cu"; break;
2421 case 30: p = "Zn"; break;
2422 case 35: p = "Br"; break;
2423 case 47: p = "Ag"; break;
2424 default: p = ""; break;
2430 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2431 int file_version, gmx_groups_t *groups, int atnr)
2436 gmx_fio_do_real(fio, atom->m);
2437 gmx_fio_do_real(fio, atom->q);
2438 gmx_fio_do_real(fio, atom->mB);
2439 gmx_fio_do_real(fio, atom->qB);
2440 gmx_fio_do_ushort(fio, atom->type);
2441 gmx_fio_do_ushort(fio, atom->typeB);
2442 gmx_fio_do_int(fio, atom->ptype);
2443 gmx_fio_do_int(fio, atom->resind);
2444 if (file_version >= 52)
2446 gmx_fio_do_int(fio, atom->atomnumber);
2449 /* Set element string from atomic number if present.
2450 * This routine returns an empty string if the name is not found.
2452 strncpy(atom->elem, atomicnumber_to_element(atom->atomnumber), 4);
2453 /* avoid warnings about potentially unterminated string */
2454 atom->elem[3] = '\0';
2459 atom->atomnumber = NOTSET;
2461 if (file_version < 23)
2465 else if (file_version < 39)
2474 if (file_version < 57)
2476 unsigned char uchar[egcNR];
2477 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2478 for (i = myngrp; (i < ngrp); i++)
2482 /* Copy the old data format to the groups struct */
2483 for (i = 0; i < ngrp; i++)
2485 groups->grpnr[i][atnr] = uchar[i];
2490 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2495 if (file_version < 23)
2499 else if (file_version < 39)
2508 for (j = 0; (j < ngrp); j++)
2512 gmx_fio_do_int(fio, grps[j].nr);
2515 snew(grps[j].nm_ind, grps[j].nr);
2517 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2522 snew(grps[j].nm_ind, grps[j].nr);
2527 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2533 gmx_fio_do_int(fio, ls);
2534 *nm = get_symtab_handle(symtab, ls);
2538 ls = lookup_symtab(symtab, *nm);
2539 gmx_fio_do_int(fio, ls);
2543 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2548 for (j = 0; (j < nstr); j++)
2550 do_symstr(fio, &(nm[j]), bRead, symtab);
2554 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2555 t_symtab *symtab, int file_version)
2559 for (j = 0; (j < n); j++)
2561 do_symstr(fio, &(ri[j].name), bRead, symtab);
2562 if (file_version >= 63)
2564 gmx_fio_do_int(fio, ri[j].nr);
2565 gmx_fio_do_uchar(fio, ri[j].ic);
2575 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2577 gmx_groups_t *groups)
2581 gmx_fio_do_int(fio, atoms->nr);
2582 gmx_fio_do_int(fio, atoms->nres);
2583 if (file_version < 57)
2585 gmx_fio_do_int(fio, groups->ngrpname);
2586 for (i = 0; i < egcNR; i++)
2588 groups->ngrpnr[i] = atoms->nr;
2589 snew(groups->grpnr[i], groups->ngrpnr[i]);
2594 snew(atoms->atom, atoms->nr);
2595 snew(atoms->atomname, atoms->nr);
2596 snew(atoms->atomtype, atoms->nr);
2597 snew(atoms->atomtypeB, atoms->nr);
2598 snew(atoms->resinfo, atoms->nres);
2599 if (file_version < 57)
2601 snew(groups->grpname, groups->ngrpname);
2603 atoms->pdbinfo = NULL;
2605 for (i = 0; (i < atoms->nr); i++)
2607 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2609 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2610 if (bRead && (file_version <= 20))
2612 for (i = 0; i < atoms->nr; i++)
2614 atoms->atomtype[i] = put_symtab(symtab, "?");
2615 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2620 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2621 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2623 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2625 if (file_version < 57)
2627 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2629 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2633 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2634 gmx_bool bRead, t_symtab *symtab,
2639 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2640 gmx_fio_do_int(fio, groups->ngrpname);
2643 snew(groups->grpname, groups->ngrpname);
2645 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2646 for (g = 0; g < egcNR; g++)
2648 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2649 if (groups->ngrpnr[g] == 0)
2653 groups->grpnr[g] = NULL;
2660 snew(groups->grpnr[g], groups->ngrpnr[g]);
2662 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2667 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2672 if (file_version > 25)
2674 gmx_fio_do_int(fio, atomtypes->nr);
2678 snew(atomtypes->radius, j);
2679 snew(atomtypes->vol, j);
2680 snew(atomtypes->surftens, j);
2681 snew(atomtypes->atomnumber, j);
2682 snew(atomtypes->gb_radius, j);
2683 snew(atomtypes->S_hct, j);
2685 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2686 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2687 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2688 if (file_version >= 40)
2690 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2692 if (file_version >= 60)
2694 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2695 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2700 /* File versions prior to 26 cannot do GBSA,
2701 * so they dont use this structure
2704 atomtypes->radius = NULL;
2705 atomtypes->vol = NULL;
2706 atomtypes->surftens = NULL;
2707 atomtypes->atomnumber = NULL;
2708 atomtypes->gb_radius = NULL;
2709 atomtypes->S_hct = NULL;
2713 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2719 gmx_fio_do_int(fio, symtab->nr);
2723 snew(symtab->symbuf, 1);
2724 symbuf = symtab->symbuf;
2725 symbuf->bufsize = nr;
2726 snew(symbuf->buf, nr);
2727 for (i = 0; (i < nr); i++)
2729 gmx_fio_do_string(fio, buf);
2730 symbuf->buf[i] = gmx_strdup(buf);
2735 symbuf = symtab->symbuf;
2736 while (symbuf != NULL)
2738 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2740 gmx_fio_do_string(fio, symbuf->buf[i]);
2743 symbuf = symbuf->next;
2747 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2752 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2754 int i, j, ngrid, gs, nelem;
2756 gmx_fio_do_int(fio, cmap_grid->ngrid);
2757 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2759 ngrid = cmap_grid->ngrid;
2760 gs = cmap_grid->grid_spacing;
2765 snew(cmap_grid->cmapdata, ngrid);
2767 for (i = 0; i < cmap_grid->ngrid; i++)
2769 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2773 for (i = 0; i < cmap_grid->ngrid; i++)
2775 for (j = 0; j < nelem; j++)
2777 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2778 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2779 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2780 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2786 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2788 int m, a, a0, a1, r;
2792 /* We always assign a new chain number, but save the chain id characters
2793 * for larger molecules.
2795 #define CHAIN_MIN_ATOMS 15
2799 for (m = 0; m < mols->nr; m++)
2801 a0 = mols->index[m];
2802 a1 = mols->index[m+1];
2803 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2812 for (a = a0; a < a1; a++)
2814 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2815 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2820 /* Blank out the chain id if there was only one chain */
2823 for (r = 0; r < atoms->nres; r++)
2825 atoms->resinfo[r].chainid = ' ';
2830 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2831 t_symtab *symtab, int file_version,
2832 gmx_groups_t *groups)
2836 if (file_version >= 57)
2838 do_symstr(fio, &(molt->name), bRead, symtab);
2841 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2843 if (bRead && gmx_debug_at)
2845 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2848 if (file_version >= 57)
2850 do_ilists(fio, molt->ilist, bRead, file_version);
2852 do_block(fio, &molt->cgs, bRead, file_version);
2853 if (bRead && gmx_debug_at)
2855 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2859 /* This used to be in the atoms struct */
2860 do_blocka(fio, &molt->excls, bRead, file_version);
2863 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2867 gmx_fio_do_int(fio, molb->type);
2868 gmx_fio_do_int(fio, molb->nmol);
2869 gmx_fio_do_int(fio, molb->natoms_mol);
2870 /* Position restraint coordinates */
2871 gmx_fio_do_int(fio, molb->nposres_xA);
2872 if (molb->nposres_xA > 0)
2876 snew(molb->posres_xA, molb->nposres_xA);
2878 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2880 gmx_fio_do_int(fio, molb->nposres_xB);
2881 if (molb->nposres_xB > 0)
2885 snew(molb->posres_xB, molb->nposres_xB);
2887 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2892 static t_block mtop_mols(gmx_mtop_t *mtop)
2898 for (mb = 0; mb < mtop->nmolblock; mb++)
2900 mols.nr += mtop->molblock[mb].nmol;
2902 mols.nalloc_index = mols.nr + 1;
2903 snew(mols.index, mols.nalloc_index);
2908 for (mb = 0; mb < mtop->nmolblock; mb++)
2910 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2912 a += mtop->molblock[mb].natoms_mol;
2921 static void add_posres_molblock(gmx_mtop_t *mtop)
2926 gmx_molblock_t *molb;
2929 /* posres reference positions are stored in ip->posres (if present) and
2930 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2931 posres.pos0A are identical to fbposres.pos0. */
2932 il = &mtop->moltype[0].ilist[F_POSRES];
2933 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2934 if (il->nr == 0 && ilfb->nr == 0)
2940 for (i = 0; i < il->nr; i += 2)
2942 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2943 am = max(am, il->iatoms[i+1]);
2944 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2945 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2946 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2951 /* This loop is required if we have only flat-bottomed posres:
2953 - bFE == FALSE (no B-state for flat-bottomed posres) */
2956 for (i = 0; i < ilfb->nr; i += 2)
2958 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2959 am = max(am, ilfb->iatoms[i+1]);
2962 /* Make the posres coordinate block end at a molecule end */
2964 while (am >= mtop->mols.index[mol+1])
2968 molb = &mtop->molblock[0];
2969 molb->nposres_xA = mtop->mols.index[mol+1];
2970 snew(molb->posres_xA, molb->nposres_xA);
2973 molb->nposres_xB = molb->nposres_xA;
2974 snew(molb->posres_xB, molb->nposres_xB);
2978 molb->nposres_xB = 0;
2980 for (i = 0; i < il->nr; i += 2)
2982 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2983 a = il->iatoms[i+1];
2984 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2985 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2986 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2989 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2990 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2991 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2996 /* If only flat-bottomed posres are present, take reference pos from them.
2997 Here: bFE == FALSE */
2998 for (i = 0; i < ilfb->nr; i += 2)
3000 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
3001 a = ilfb->iatoms[i+1];
3002 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
3003 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
3004 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
3009 static void set_disres_npair(gmx_mtop_t *mtop)
3012 gmx_mtop_ilistloop_t iloop;
3013 t_ilist *ilist, *il;
3017 ip = mtop->ffparams.iparams;
3019 iloop = gmx_mtop_ilistloop_init(mtop);
3020 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
3022 il = &ilist[F_DISRES];
3028 for (i = 0; i < il->nr; i += 3)
3031 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
3033 ip[a[i]].disres.npair = npair;
3041 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
3051 do_symtab(fio, &(mtop->symtab), bRead);
3054 pr_symtab(debug, 0, "symtab", &mtop->symtab);
3057 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
3059 if (file_version >= 57)
3061 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
3063 gmx_fio_do_int(fio, mtop->nmoltype);
3071 snew(mtop->moltype, mtop->nmoltype);
3072 if (file_version < 57)
3074 mtop->moltype[0].name = mtop->name;
3077 for (mt = 0; mt < mtop->nmoltype; mt++)
3079 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
3083 if (file_version >= 57)
3085 gmx_fio_do_int(fio, mtop->nmolblock);
3089 mtop->nmolblock = 1;
3093 snew(mtop->molblock, mtop->nmolblock);
3095 if (file_version >= 57)
3097 for (mb = 0; mb < mtop->nmolblock; mb++)
3099 do_molblock(fio, &mtop->molblock[mb], bRead);
3101 gmx_fio_do_int(fio, mtop->natoms);
3105 mtop->molblock[0].type = 0;
3106 mtop->molblock[0].nmol = 1;
3107 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
3108 mtop->molblock[0].nposres_xA = 0;
3109 mtop->molblock[0].nposres_xB = 0;
3112 if (file_version >= tpxv_IntermolecularBondeds)
3114 gmx_fio_do_gmx_bool(fio, mtop->bIntermolecularInteractions);
3115 if (mtop->bIntermolecularInteractions)
3119 snew(mtop->intermolecular_ilist, F_NRE);
3121 do_ilists(fio, mtop->intermolecular_ilist, bRead, file_version);
3126 mtop->bIntermolecularInteractions = FALSE;
3129 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
3132 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
3135 if (file_version < 57)
3137 /* Debug statements are inside do_idef */
3138 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
3139 mtop->natoms = mtop->moltype[0].atoms.nr;
3142 if (file_version >= 65)
3144 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
3148 mtop->ffparams.cmap_grid.ngrid = 0;
3149 mtop->ffparams.cmap_grid.grid_spacing = 0;
3150 mtop->ffparams.cmap_grid.cmapdata = NULL;
3153 if (file_version >= 57)
3155 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
3158 if (file_version < 57)
3160 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
3161 if (bRead && gmx_debug_at)
3163 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
3165 do_block(fio, &mtop->mols, bRead, file_version);
3166 /* Add the posres coordinates to the molblock */
3167 add_posres_molblock(mtop);
3171 if (file_version >= 57)
3173 done_block(&mtop->mols);
3174 mtop->mols = mtop_mols(mtop);
3178 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
3182 if (file_version < 51)
3184 /* Here used to be the shake blocks */
3185 do_blocka(fio, &dumb, bRead, file_version);
3198 close_symtab(&(mtop->symtab));
3202 /* If TopOnlyOK is TRUE then we can read even future versions
3203 * of tpx files, provided the file_generation hasn't changed.
3204 * If it is FALSE, we need the inputrecord too, and bail out
3205 * if the file is newer than the program.
3207 * The version and generation if the topology (see top of this file)
3208 * are returned in the two last arguments.
3210 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
3212 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
3213 gmx_bool TopOnlyOK, int *file_version,
3214 int *file_generation)
3217 char file_tag[STRLEN];
3224 gmx_fio_checktype(fio);
3225 gmx_fio_setdebug(fio, bDebugMode());
3227 /* XDR binary topology file */
3228 precision = sizeof(real);
3231 gmx_fio_do_string(fio, buf);
3232 if (strncmp(buf, "VERSION", 7))
3234 gmx_fatal(FARGS, "Can not read file %s,\n"
3235 " this file is from a GROMACS version which is older than 2.0\n"
3236 " Make a new one with grompp or use a gro or pdb file, if possible",
3237 gmx_fio_getname(fio));
3239 gmx_fio_do_int(fio, precision);
3240 bDouble = (precision == sizeof(double));
3241 if ((precision != sizeof(float)) && !bDouble)
3243 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3244 "instead of %d or %d",
3245 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3247 gmx_fio_setprecision(fio, bDouble);
3248 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3249 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3253 gmx_fio_write_string(fio, GromacsVersion());
3254 bDouble = (precision == sizeof(double));
3255 gmx_fio_setprecision(fio, bDouble);
3256 gmx_fio_do_int(fio, precision);
3258 sprintf(file_tag, "%s", tpx_tag);
3259 fgen = tpx_generation;
3262 /* Check versions! */
3263 gmx_fio_do_int(fio, fver);
3265 /* This is for backward compatibility with development versions 77-79
3266 * where the tag was, mistakenly, placed before the generation,
3267 * which would cause a segv instead of a proper error message
3268 * when reading the topology only from tpx with <77 code.
3270 if (fver >= 77 && fver <= 79)
3272 gmx_fio_do_string(fio, file_tag);
3277 gmx_fio_do_int(fio, fgen);
3286 gmx_fio_do_string(fio, file_tag);
3292 /* Versions before 77 don't have the tag, set it to release */
3293 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3296 if (strcmp(file_tag, tpx_tag) != 0)
3298 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3301 /* We only support reading tpx files with the same tag as the code
3302 * or tpx files with the release tag and with lower version number.
3304 if (strcmp(file_tag, TPX_TAG_RELEASE) != 0 && fver < tpx_version)
3306 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3307 gmx_fio_getname(fio), fver, file_tag,
3308 tpx_version, tpx_tag);
3313 if (file_version != NULL)
3315 *file_version = fver;
3317 if (file_generation != NULL)
3319 *file_generation = fgen;
3323 if ((fver <= tpx_incompatible_version) ||
3324 ((fver > tpx_version) && !TopOnlyOK) ||
3325 (fgen > tpx_generation) ||
3326 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3328 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3329 gmx_fio_getname(fio), fver, tpx_version);
3332 gmx_fio_do_int(fio, tpx->natoms);
3335 gmx_fio_do_int(fio, tpx->ngtc);
3343 gmx_fio_do_int(fio, idum);
3344 gmx_fio_do_real(fio, rdum);
3346 /*a better decision will eventually (5.0 or later) need to be made
3347 on how to treat the alchemical state of the system, which can now
3348 vary through a simulation, and cannot be completely described
3349 though a single lambda variable, or even a single state
3350 index. Eventually, should probably be a vector. MRS*/
3353 gmx_fio_do_int(fio, tpx->fep_state);
3355 gmx_fio_do_real(fio, tpx->lambda);
3356 gmx_fio_do_int(fio, tpx->bIr);
3357 gmx_fio_do_int(fio, tpx->bTop);
3358 gmx_fio_do_int(fio, tpx->bX);
3359 gmx_fio_do_int(fio, tpx->bV);
3360 gmx_fio_do_int(fio, tpx->bF);
3361 gmx_fio_do_int(fio, tpx->bBox);
3363 if ((fgen > tpx_generation))
3365 /* This can only happen if TopOnlyOK=TRUE */
3370 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3371 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3372 gmx_bool bXVallocated)
3378 int file_version, file_generation;
3382 gmx_bool bPeriodicMols;
3386 tpx.natoms = state->natoms;
3387 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3388 tpx.fep_state = state->fep_state;
3389 tpx.lambda = state->lambda[efptFEP];
3390 tpx.bIr = (ir != NULL);
3391 tpx.bTop = (mtop != NULL);
3392 tpx.bX = (state->x != NULL);
3393 tpx.bV = (state->v != NULL);
3394 tpx.bF = (f != NULL);
3398 TopOnlyOK = (ir == NULL);
3400 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3405 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3406 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3411 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3412 state->natoms = tpx.natoms;
3413 state->nalloc = tpx.natoms;
3419 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3423 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3425 do_test(fio, tpx.bBox, state->box);
3428 gmx_fio_ndo_rvec(fio, state->box, DIM);
3429 if (file_version >= 51)
3431 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3435 /* We initialize box_rel after reading the inputrec */
3436 clear_mat(state->box_rel);
3438 if (file_version >= 28)
3440 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3441 if (file_version < 56)
3444 gmx_fio_ndo_rvec(fio, mdum, DIM);
3449 if (state->ngtc > 0 && file_version >= 28)
3452 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3453 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3454 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3455 snew(dumv, state->ngtc);
3456 if (file_version < 69)
3458 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3460 /* These used to be the Berendsen tcoupl_lambda's */
3461 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3465 /* Prior to tpx version 26, the inputrec was here.
3466 * I moved it to enable partial forward-compatibility
3467 * for analysis/viewer programs.
3469 if (file_version < 26)
3471 do_test(fio, tpx.bIr, ir);
3476 do_inputrec(fio, ir, bRead, file_version,
3477 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3480 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3485 do_inputrec(fio, &dum_ir, bRead, file_version,
3486 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3489 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3491 done_inputrec(&dum_ir);
3497 do_test(fio, tpx.bTop, mtop);
3502 do_mtop(fio, mtop, bRead, file_version);
3506 do_mtop(fio, &dum_top, bRead, file_version);
3507 done_mtop(&dum_top, TRUE);
3510 do_test(fio, tpx.bX, state->x);
3515 state->flags |= (1<<estX);
3517 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3520 do_test(fio, tpx.bV, state->v);
3525 state->flags |= (1<<estV);
3527 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3530 do_test(fio, tpx.bF, f);
3533 gmx_fio_ndo_rvec(fio, f, state->natoms);
3536 /* Starting with tpx version 26, we have the inputrec
3537 * at the end of the file, so we can ignore it
3538 * if the file is never than the software (but still the
3539 * same generation - see comments at the top of this file.
3544 bPeriodicMols = FALSE;
3545 if (file_version >= 26)
3547 do_test(fio, tpx.bIr, ir);
3550 if (file_version >= 53)
3552 /* Removed the pbc info from do_inputrec, since we always want it */
3556 bPeriodicMols = ir->bPeriodicMols;
3558 gmx_fio_do_int(fio, ePBC);
3559 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3561 if (file_generation <= tpx_generation && ir)
3563 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3566 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3568 if (file_version < 51)
3570 set_box_rel(ir, state);
3572 if (file_version < 53)
3575 bPeriodicMols = ir->bPeriodicMols;
3578 if (bRead && ir && file_version >= 53)
3580 /* We need to do this after do_inputrec, since that initializes ir */
3582 ir->bPeriodicMols = bPeriodicMols;
3591 if (state->ngtc == 0)
3593 /* Reading old version without tcoupl state data: set it */
3594 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3596 if (tpx.bTop && mtop)
3598 if (file_version < 57)
3600 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3602 ir->eDisre = edrSimple;
3606 ir->eDisre = edrNone;
3609 set_disres_npair(mtop);
3613 if (tpx.bTop && mtop)
3615 gmx_mtop_finalize(mtop);
3618 if (file_version >= 57)
3622 env = getenv("GMX_NOCHARGEGROUPS");
3625 sscanf(env, "%d", &ienv);
3626 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3631 "Will make single atomic charge groups in non-solvent%s\n",
3632 ienv > 1 ? " and solvent" : "");
3633 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3635 fprintf(stderr, "\n");
3643 /************************************************************
3645 * The following routines are the exported ones
3647 ************************************************************/
3649 t_fileio *open_tpx(const char *fn, const char *mode)
3651 return gmx_fio_open(fn, mode);
3654 void close_tpx(t_fileio *fio)
3659 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3660 int *file_version, int *file_generation)
3664 fio = open_tpx(fn, "r");
3665 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3669 void write_tpx_state(const char *fn,
3670 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3674 fio = open_tpx(fn, "w");
3675 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3679 void read_tpx_state(const char *fn,
3680 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3684 fio = open_tpx(fn, "r");
3685 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3689 int read_tpx(const char *fn,
3690 t_inputrec *ir, matrix box, int *natoms,
3691 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3699 fio = open_tpx(fn, "r");
3700 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3702 *natoms = state.natoms;
3705 copy_mat(state.box, box);
3714 int read_tpx_top(const char *fn,
3715 t_inputrec *ir, matrix box, int *natoms,
3716 rvec *x, rvec *v, rvec *f, t_topology *top)
3722 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3724 *top = gmx_mtop_t_to_t_topology(&mtop);
3729 gmx_bool fn2bTPX(const char *file)
3731 return (efTPR == fn2ftp(file));
3734 static void done_gmx_groups_t(gmx_groups_t *g)
3738 for (i = 0; (i < egcNR); i++)
3740 if (NULL != g->grps[i].nm_ind)
3742 sfree(g->grps[i].nm_ind);
3743 g->grps[i].nm_ind = NULL;
3745 if (NULL != g->grpnr[i])
3751 /* The contents of this array is in symtab, don't free it here */
3755 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3756 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3759 int natoms, i, version, generation;
3760 gmx_bool bTop, bXNULL = FALSE;
3762 t_topology *topconv;
3765 bTop = fn2bTPX(infile);
3769 read_tpxheader(infile, &header, TRUE, &version, &generation);
3772 snew(*x, header.natoms);
3776 snew(*v, header.natoms);
3779 *ePBC = read_tpx(infile, NULL, box, &natoms,
3780 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3781 *top = gmx_mtop_t_to_t_topology(mtop);
3782 /* In this case we need to throw away the group data too */
3783 done_gmx_groups_t(&mtop->groups);
3785 strcpy(title, *top->name);
3786 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3790 get_stx_coordnum(infile, &natoms);
3791 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3802 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3810 aps = gmx_atomprop_init();
3811 for (i = 0; (i < natoms); i++)
3813 if (!gmx_atomprop_query(aps, epropMass,
3814 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3815 *top->atoms.atomname[i],
3816 &(top->atoms.atom[i].m)))
3820 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3821 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3822 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3823 *top->atoms.atomname[i]);
3827 gmx_atomprop_destroy(aps);
3829 top->idef.ntypes = -1;