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, 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.
41 /* This file is completely threadsafe - keep it that way! */
48 #include "gromacs/utility/futil.h"
56 #include "mtop_util.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/topology/block.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/smalloc.h"
66 #define TPX_TAG_RELEASE "release"
68 /*! \brief Tag string for the file format written to run input files
69 * written by this version of the code.
71 * Change this if you want to change the run input format in a feature
72 * branch. This ensures that there will not be different run input
73 * formats around which cannot be distinguished, while not causing
74 * problems rebasing the feature branch onto upstream changes. When
75 * merging with mainstream GROMACS, set this tag string back to
76 * TPX_TAG_RELEASE, and instead add an element to tpxv and set
77 * tpx_version to that.
79 static const char *tpx_tag = TPX_TAG_RELEASE;
81 /*! \brief Enum of values that describe the contents of a tpr file
82 * whose format matches a version number
84 * The enum helps the code be more self-documenting and ensure merges
85 * do not silently resolve when two patches make the same bump. When
86 * adding new functionality, add a new element to the end of this
87 * enumeration, change the definition of tpx_version, and write code
88 * below that does the right thing according to the value of
91 tpxv_ComputationalElectrophysiology = 96, /**< support for ion/water position swaps (computational electrophysiology) */
92 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to gmx_int64_t */
93 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
94 tpxv_InteractiveMolecularDynamics /**< interactive molecular dynamics (IMD) */
97 /*! \brief Version number of the file format written to run input
98 * files by this version of the code.
100 * The tpx_version number should be increased whenever the file format
101 * in the main development branch changes, generally to the highest
102 * value present in tpxv. Backward compatibility for reading old run
103 * input files is maintained by checking this version number against
104 * that of the file and then using the correct code path.
106 * When developing a feature branch that needs to change the run input
107 * file format, change tpx_tag instead. */
108 static const int tpx_version = tpxv_InteractiveMolecularDynamics;
111 /* This number should only be increased when you edit the TOPOLOGY section
112 * or the HEADER of the tpx format.
113 * This way we can maintain forward compatibility too for all analysis tools
114 * and/or external programs that only need to know the atom/residue names,
115 * charges, and bond connectivity.
117 * It first appeared in tpx version 26, when I also moved the inputrecord
118 * to the end of the tpx file, so we can just skip it if we only
121 * In particular, it must be increased when adding new elements to
122 * ftupd, so that old code can read new .tpr files.
124 static const int tpx_generation = 26;
126 /* This number should be the most recent backwards incompatible version
127 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
129 static const int tpx_incompatible_version = 9;
133 /* Struct used to maintain tpx compatibility when function types are added */
135 int fvnr; /* file version number in which the function type first appeared */
136 int ftype; /* function type */
140 * The entries should be ordered in:
141 * 1. ascending file version number
142 * 2. ascending function type number
144 /*static const t_ftupd ftupd[] = {
145 { 20, F_CUBICBONDS },
149 { 22, F_DISRESVIOL },
155 { 26, F_DIHRESVIOL },
156 { 30, F_CROSS_BOND_BONDS },
157 { 30, F_CROSS_BOND_ANGLES },
158 { 30, F_UREY_BRADLEY },
159 { 30, F_POLARIZATION },
163 * The entries should be ordered in:
164 * 1. ascending function type number
165 * 2. ascending file version number
167 /* question; what is the purpose of the commented code above? */
168 static const t_ftupd ftupd[] = {
169 { 20, F_CUBICBONDS },
174 { 43, F_TABBONDSNC },
175 { 70, F_RESTRBONDS },
176 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
177 { 76, F_LINEAR_ANGLES },
178 { 30, F_CROSS_BOND_BONDS },
179 { 30, F_CROSS_BOND_ANGLES },
180 { 30, F_UREY_BRADLEY },
181 { 34, F_QUARTIC_ANGLES },
183 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
184 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
193 { 72, F_NPSOLVATION },
195 { 41, F_LJC_PAIRS_NB },
198 { 32, F_COUL_RECIP },
201 { 30, F_POLARIZATION },
204 { 22, F_DISRESVIOL },
208 { 26, F_DIHRESVIOL },
213 { 46, F_ECONSERVED },
214 { 69, F_VTEMP_NOLONGERUSED},
216 { 54, F_DVDL_CONSTR },
217 { 76, F_ANHARM_POL },
220 { 79, F_DVDL_BONDED, },
221 { 79, F_DVDL_RESTRAINT },
222 { 79, F_DVDL_TEMPERATURE },
224 #define NFTUPD asize(ftupd)
226 /* Needed for backward compatibility */
229 static void _do_section(t_fileio *fio, int key, gmx_bool bRead, const char *src,
235 if (gmx_fio_getftp(fio) == efTPA)
239 gmx_fio_write_string(fio, itemstr[key]);
240 bDbg = gmx_fio_getdebug(fio);
241 gmx_fio_setdebug(fio, FALSE);
242 gmx_fio_write_string(fio, comment_str[key]);
243 gmx_fio_setdebug(fio, bDbg);
247 if (gmx_fio_getdebug(fio))
249 fprintf(stderr, "Looking for section %s (%s, %d)",
250 itemstr[key], src, line);
255 gmx_fio_do_string(fio, buf);
257 while ((gmx_strcasecmp(buf, itemstr[key]) != 0));
259 if (gmx_strcasecmp(buf, itemstr[key]) != 0)
261 gmx_fatal(FARGS, "\nCould not find section heading %s", itemstr[key]);
263 else if (gmx_fio_getdebug(fio))
265 fprintf(stderr, " and found it\n");
271 #define do_section(fio, key, bRead) _do_section(fio, key, bRead, __FILE__, __LINE__)
273 /**************************************************************
275 * Now the higer level routines that do io of the structures and arrays
277 **************************************************************/
278 static void do_pullgrp_tpx_pre95(t_fileio *fio,
287 gmx_fio_do_int(fio, pgrp->nat);
290 snew(pgrp->ind, pgrp->nat);
292 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
293 gmx_fio_do_int(fio, pgrp->nweight);
296 snew(pgrp->weight, pgrp->nweight);
298 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
299 gmx_fio_do_int(fio, pgrp->pbcatom);
300 gmx_fio_do_rvec(fio, pcrd->vec);
301 clear_rvec(pcrd->origin);
302 gmx_fio_do_rvec(fio, tmp);
304 gmx_fio_do_real(fio, pcrd->rate);
305 gmx_fio_do_real(fio, pcrd->k);
306 if (file_version >= 56)
308 gmx_fio_do_real(fio, pcrd->kB);
316 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
320 gmx_fio_do_int(fio, pgrp->nat);
323 snew(pgrp->ind, pgrp->nat);
325 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
326 gmx_fio_do_int(fio, pgrp->nweight);
329 snew(pgrp->weight, pgrp->nweight);
331 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
332 gmx_fio_do_int(fio, pgrp->pbcatom);
335 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd)
339 gmx_fio_do_int(fio, pcrd->group[0]);
340 gmx_fio_do_int(fio, pcrd->group[1]);
341 gmx_fio_do_rvec(fio, pcrd->origin);
342 gmx_fio_do_rvec(fio, pcrd->vec);
343 gmx_fio_do_real(fio, pcrd->init);
344 gmx_fio_do_real(fio, pcrd->rate);
345 gmx_fio_do_real(fio, pcrd->k);
346 gmx_fio_do_real(fio, pcrd->kB);
349 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
351 /* i is used in the ndo_double macro*/
355 int n_lambda = fepvals->n_lambda;
357 /* reset the lambda calculation window */
358 fepvals->lambda_start_n = 0;
359 fepvals->lambda_stop_n = n_lambda;
360 if (file_version >= 79)
366 snew(expand->init_lambda_weights, n_lambda);
368 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
369 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
372 gmx_fio_do_int(fio, expand->nstexpanded);
373 gmx_fio_do_int(fio, expand->elmcmove);
374 gmx_fio_do_int(fio, expand->elamstats);
375 gmx_fio_do_int(fio, expand->lmc_repeats);
376 gmx_fio_do_int(fio, expand->gibbsdeltalam);
377 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
378 gmx_fio_do_int(fio, expand->lmc_seed);
379 gmx_fio_do_real(fio, expand->mc_temp);
380 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
381 gmx_fio_do_int(fio, expand->nstTij);
382 gmx_fio_do_int(fio, expand->minvarmin);
383 gmx_fio_do_int(fio, expand->c_range);
384 gmx_fio_do_real(fio, expand->wl_scale);
385 gmx_fio_do_real(fio, expand->wl_ratio);
386 gmx_fio_do_real(fio, expand->init_wl_delta);
387 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
388 gmx_fio_do_int(fio, expand->elmceq);
389 gmx_fio_do_int(fio, expand->equil_steps);
390 gmx_fio_do_int(fio, expand->equil_samples);
391 gmx_fio_do_int(fio, expand->equil_n_at_lam);
392 gmx_fio_do_real(fio, expand->equil_wl_delta);
393 gmx_fio_do_real(fio, expand->equil_ratio);
397 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
400 if (file_version >= 79)
402 gmx_fio_do_int(fio, simtemp->eSimTempScale);
403 gmx_fio_do_real(fio, simtemp->simtemp_high);
404 gmx_fio_do_real(fio, simtemp->simtemp_low);
409 snew(simtemp->temperatures, n_lambda);
411 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
416 static void do_imd(t_fileio *fio, t_IMD *imd, gmx_bool bRead)
418 gmx_fio_do_int(fio, imd->nat);
421 snew(imd->ind, imd->nat);
423 gmx_fio_ndo_int(fio, imd->ind, imd->nat);
426 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
428 /* i is defined in the ndo_double macro; use g to iterate. */
433 /* free energy values */
435 if (file_version >= 79)
437 gmx_fio_do_int(fio, fepvals->init_fep_state);
438 gmx_fio_do_double(fio, fepvals->init_lambda);
439 gmx_fio_do_double(fio, fepvals->delta_lambda);
441 else if (file_version >= 59)
443 gmx_fio_do_double(fio, fepvals->init_lambda);
444 gmx_fio_do_double(fio, fepvals->delta_lambda);
448 gmx_fio_do_real(fio, rdum);
449 fepvals->init_lambda = rdum;
450 gmx_fio_do_real(fio, rdum);
451 fepvals->delta_lambda = rdum;
453 if (file_version >= 79)
455 gmx_fio_do_int(fio, fepvals->n_lambda);
458 snew(fepvals->all_lambda, efptNR);
460 for (g = 0; g < efptNR; g++)
462 if (fepvals->n_lambda > 0)
466 snew(fepvals->all_lambda[g], fepvals->n_lambda);
468 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
469 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
471 else if (fepvals->init_lambda >= 0)
473 fepvals->separate_dvdl[efptFEP] = TRUE;
477 else if (file_version >= 64)
479 gmx_fio_do_int(fio, fepvals->n_lambda);
484 snew(fepvals->all_lambda, efptNR);
485 /* still allocate the all_lambda array's contents. */
486 for (g = 0; g < efptNR; g++)
488 if (fepvals->n_lambda > 0)
490 snew(fepvals->all_lambda[g], fepvals->n_lambda);
494 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
496 if (fepvals->init_lambda >= 0)
500 fepvals->separate_dvdl[efptFEP] = TRUE;
504 /* copy the contents of the efptFEP lambda component to all
505 the other components */
506 for (g = 0; g < efptNR; g++)
508 for (h = 0; h < fepvals->n_lambda; h++)
512 fepvals->all_lambda[g][h] =
513 fepvals->all_lambda[efptFEP][h];
522 fepvals->n_lambda = 0;
523 fepvals->all_lambda = NULL;
524 if (fepvals->init_lambda >= 0)
526 fepvals->separate_dvdl[efptFEP] = TRUE;
529 if (file_version >= 13)
531 gmx_fio_do_real(fio, fepvals->sc_alpha);
535 fepvals->sc_alpha = 0;
537 if (file_version >= 38)
539 gmx_fio_do_int(fio, fepvals->sc_power);
543 fepvals->sc_power = 2;
545 if (file_version >= 79)
547 gmx_fio_do_real(fio, fepvals->sc_r_power);
551 fepvals->sc_r_power = 6.0;
553 if (file_version >= 15)
555 gmx_fio_do_real(fio, fepvals->sc_sigma);
559 fepvals->sc_sigma = 0.3;
563 if (file_version >= 71)
565 fepvals->sc_sigma_min = fepvals->sc_sigma;
569 fepvals->sc_sigma_min = 0;
572 if (file_version >= 79)
574 gmx_fio_do_int(fio, fepvals->bScCoul);
578 fepvals->bScCoul = TRUE;
580 if (file_version >= 64)
582 gmx_fio_do_int(fio, fepvals->nstdhdl);
586 fepvals->nstdhdl = 1;
589 if (file_version >= 73)
591 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
592 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
596 fepvals->separate_dhdl_file = esepdhdlfileYES;
597 fepvals->dhdl_derivatives = edhdlderivativesYES;
599 if (file_version >= 71)
601 gmx_fio_do_int(fio, fepvals->dh_hist_size);
602 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
606 fepvals->dh_hist_size = 0;
607 fepvals->dh_hist_spacing = 0.1;
609 if (file_version >= 79)
611 gmx_fio_do_int(fio, fepvals->bPrintEnergy);
615 fepvals->bPrintEnergy = FALSE;
618 /* handle lambda_neighbors */
619 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
621 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
622 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
623 (fepvals->init_lambda < 0) )
625 fepvals->lambda_start_n = (fepvals->init_fep_state -
626 fepvals->lambda_neighbors);
627 fepvals->lambda_stop_n = (fepvals->init_fep_state +
628 fepvals->lambda_neighbors + 1);
629 if (fepvals->lambda_start_n < 0)
631 fepvals->lambda_start_n = 0;;
633 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
635 fepvals->lambda_stop_n = fepvals->n_lambda;
640 fepvals->lambda_start_n = 0;
641 fepvals->lambda_stop_n = fepvals->n_lambda;
646 fepvals->lambda_start_n = 0;
647 fepvals->lambda_stop_n = fepvals->n_lambda;
651 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
655 if (file_version >= 95)
657 gmx_fio_do_int(fio, pull->ngroup);
659 gmx_fio_do_int(fio, pull->ncoord);
660 if (file_version < 95)
662 pull->ngroup = pull->ncoord + 1;
664 gmx_fio_do_int(fio, pull->eGeom);
665 gmx_fio_do_ivec(fio, pull->dim);
666 gmx_fio_do_real(fio, pull->cyl_r1);
667 gmx_fio_do_real(fio, pull->cyl_r0);
668 gmx_fio_do_real(fio, pull->constr_tol);
669 if (file_version >= 95)
671 gmx_fio_do_int(fio, pull->bPrintRef);
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 (pull->eGeom == epullgDIRPBC)
685 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
687 if (pull->eGeom > 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[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->bPrintRef = (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]);
720 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
724 gmx_fio_do_int(fio, rotg->eType);
725 gmx_fio_do_int(fio, rotg->bMassW);
726 gmx_fio_do_int(fio, rotg->nat);
729 snew(rotg->ind, rotg->nat);
731 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
734 snew(rotg->x_ref, rotg->nat);
736 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
737 gmx_fio_do_rvec(fio, rotg->vec);
738 gmx_fio_do_rvec(fio, rotg->pivot);
739 gmx_fio_do_real(fio, rotg->rate);
740 gmx_fio_do_real(fio, rotg->k);
741 gmx_fio_do_real(fio, rotg->slab_dist);
742 gmx_fio_do_real(fio, rotg->min_gaussian);
743 gmx_fio_do_real(fio, rotg->eps);
744 gmx_fio_do_int(fio, rotg->eFittype);
745 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
746 gmx_fio_do_real(fio, rotg->PotAngle_step);
749 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
753 gmx_fio_do_int(fio, rot->ngrp);
754 gmx_fio_do_int(fio, rot->nstrout);
755 gmx_fio_do_int(fio, rot->nstsout);
758 snew(rot->grp, rot->ngrp);
760 for (g = 0; g < rot->ngrp; g++)
762 do_rotgrp(fio, &rot->grp[g], bRead);
767 static void do_swapcoords(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead)
772 gmx_fio_do_int(fio, swap->nat);
773 gmx_fio_do_int(fio, swap->nat_sol);
774 for (j = 0; j < 2; j++)
776 gmx_fio_do_int(fio, swap->nat_split[j]);
777 gmx_fio_do_int(fio, swap->massw_split[j]);
779 gmx_fio_do_int(fio, swap->nstswap);
780 gmx_fio_do_int(fio, swap->nAverage);
781 gmx_fio_do_real(fio, swap->threshold);
782 gmx_fio_do_real(fio, swap->cyl0r);
783 gmx_fio_do_real(fio, swap->cyl0u);
784 gmx_fio_do_real(fio, swap->cyl0l);
785 gmx_fio_do_real(fio, swap->cyl1r);
786 gmx_fio_do_real(fio, swap->cyl1u);
787 gmx_fio_do_real(fio, swap->cyl1l);
791 snew(swap->ind, swap->nat);
792 snew(swap->ind_sol, swap->nat_sol);
793 for (j = 0; j < 2; j++)
795 snew(swap->ind_split[j], swap->nat_split[j]);
799 gmx_fio_ndo_int(fio, swap->ind, swap->nat);
800 gmx_fio_ndo_int(fio, swap->ind_sol, swap->nat_sol);
801 for (j = 0; j < 2; j++)
803 gmx_fio_ndo_int(fio, swap->ind_split[j], swap->nat_split[j]);
806 for (j = 0; j < eCompNR; j++)
808 gmx_fio_do_int(fio, swap->nanions[j]);
809 gmx_fio_do_int(fio, swap->ncations[j]);
815 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
816 int file_version, real *fudgeQQ)
818 int i, j, k, *tmp, idum = 0;
822 real zerotemptime, finish_t, init_temp, finish_temp;
824 if (file_version != tpx_version)
826 /* Give a warning about features that are not accessible */
827 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
828 file_version, tpx_version);
836 if (file_version == 0)
841 /* Basic inputrec stuff */
842 gmx_fio_do_int(fio, ir->eI);
843 if (file_version >= 62)
845 gmx_fio_do_int64(fio, ir->nsteps);
849 gmx_fio_do_int(fio, idum);
852 if (file_version > 25)
854 if (file_version >= 62)
856 gmx_fio_do_int64(fio, ir->init_step);
860 gmx_fio_do_int(fio, idum);
861 ir->init_step = idum;
869 if (file_version >= 58)
871 gmx_fio_do_int(fio, ir->simulation_part);
875 ir->simulation_part = 1;
878 if (file_version >= 67)
880 gmx_fio_do_int(fio, ir->nstcalcenergy);
884 ir->nstcalcenergy = 1;
886 if (file_version < 53)
888 /* The pbc info has been moved out of do_inputrec,
889 * since we always want it, also without reading the inputrec.
891 gmx_fio_do_int(fio, ir->ePBC);
892 if ((file_version <= 15) && (ir->ePBC == 2))
896 if (file_version >= 45)
898 gmx_fio_do_int(fio, ir->bPeriodicMols);
905 ir->bPeriodicMols = TRUE;
909 ir->bPeriodicMols = FALSE;
913 if (file_version >= 81)
915 gmx_fio_do_int(fio, ir->cutoff_scheme);
916 if (file_version < 94)
918 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
923 ir->cutoff_scheme = ecutsGROUP;
925 gmx_fio_do_int(fio, ir->ns_type);
926 gmx_fio_do_int(fio, ir->nstlist);
927 gmx_fio_do_int(fio, ir->ndelta);
928 if (file_version < 41)
930 gmx_fio_do_int(fio, idum);
931 gmx_fio_do_int(fio, idum);
933 if (file_version >= 45)
935 gmx_fio_do_real(fio, ir->rtpi);
941 gmx_fio_do_int(fio, ir->nstcomm);
942 if (file_version > 34)
944 gmx_fio_do_int(fio, ir->comm_mode);
946 else if (ir->nstcomm < 0)
948 ir->comm_mode = ecmANGULAR;
952 ir->comm_mode = ecmLINEAR;
954 ir->nstcomm = abs(ir->nstcomm);
956 if (file_version > 25)
958 gmx_fio_do_int(fio, ir->nstcheckpoint);
962 ir->nstcheckpoint = 0;
965 gmx_fio_do_int(fio, ir->nstcgsteep);
967 if (file_version >= 30)
969 gmx_fio_do_int(fio, ir->nbfgscorr);
976 gmx_fio_do_int(fio, ir->nstlog);
977 gmx_fio_do_int(fio, ir->nstxout);
978 gmx_fio_do_int(fio, ir->nstvout);
979 gmx_fio_do_int(fio, ir->nstfout);
980 gmx_fio_do_int(fio, ir->nstenergy);
981 gmx_fio_do_int(fio, ir->nstxout_compressed);
982 if (file_version >= 59)
984 gmx_fio_do_double(fio, ir->init_t);
985 gmx_fio_do_double(fio, ir->delta_t);
989 gmx_fio_do_real(fio, rdum);
991 gmx_fio_do_real(fio, rdum);
994 gmx_fio_do_real(fio, ir->x_compression_precision);
995 if (file_version < 19)
997 gmx_fio_do_int(fio, idum);
998 gmx_fio_do_int(fio, idum);
1000 if (file_version < 18)
1002 gmx_fio_do_int(fio, idum);
1004 if (file_version >= 81)
1006 gmx_fio_do_real(fio, ir->verletbuf_tol);
1010 ir->verletbuf_tol = 0;
1012 gmx_fio_do_real(fio, ir->rlist);
1013 if (file_version >= 67)
1015 gmx_fio_do_real(fio, ir->rlistlong);
1017 if (file_version >= 82 && file_version != 90)
1019 gmx_fio_do_int(fio, ir->nstcalclr);
1023 /* Calculate at NS steps */
1024 ir->nstcalclr = ir->nstlist;
1026 gmx_fio_do_int(fio, ir->coulombtype);
1027 if (file_version < 32 && ir->coulombtype == eelRF)
1029 ir->coulombtype = eelRF_NEC;
1031 if (file_version >= 81)
1033 gmx_fio_do_int(fio, ir->coulomb_modifier);
1037 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1039 gmx_fio_do_real(fio, ir->rcoulomb_switch);
1040 gmx_fio_do_real(fio, ir->rcoulomb);
1041 gmx_fio_do_int(fio, ir->vdwtype);
1042 if (file_version >= 81)
1044 gmx_fio_do_int(fio, ir->vdw_modifier);
1048 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1050 gmx_fio_do_real(fio, ir->rvdw_switch);
1051 gmx_fio_do_real(fio, ir->rvdw);
1052 if (file_version < 67)
1054 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1056 gmx_fio_do_int(fio, ir->eDispCorr);
1057 gmx_fio_do_real(fio, ir->epsilon_r);
1058 if (file_version >= 37)
1060 gmx_fio_do_real(fio, ir->epsilon_rf);
1064 if (EEL_RF(ir->coulombtype))
1066 ir->epsilon_rf = ir->epsilon_r;
1067 ir->epsilon_r = 1.0;
1071 ir->epsilon_rf = 1.0;
1074 if (file_version >= 29)
1076 gmx_fio_do_real(fio, ir->tabext);
1083 if (file_version > 25)
1085 gmx_fio_do_int(fio, ir->gb_algorithm);
1086 gmx_fio_do_int(fio, ir->nstgbradii);
1087 gmx_fio_do_real(fio, ir->rgbradii);
1088 gmx_fio_do_real(fio, ir->gb_saltconc);
1089 gmx_fio_do_int(fio, ir->implicit_solvent);
1093 ir->gb_algorithm = egbSTILL;
1096 ir->gb_saltconc = 0;
1097 ir->implicit_solvent = eisNO;
1099 if (file_version >= 55)
1101 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1102 gmx_fio_do_real(fio, ir->gb_obc_alpha);
1103 gmx_fio_do_real(fio, ir->gb_obc_beta);
1104 gmx_fio_do_real(fio, ir->gb_obc_gamma);
1105 if (file_version >= 60)
1107 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1108 gmx_fio_do_int(fio, ir->sa_algorithm);
1112 ir->gb_dielectric_offset = 0.009;
1113 ir->sa_algorithm = esaAPPROX;
1115 gmx_fio_do_real(fio, ir->sa_surface_tension);
1117 /* Override sa_surface_tension if it is not changed in the mpd-file */
1118 if (ir->sa_surface_tension < 0)
1120 if (ir->gb_algorithm == egbSTILL)
1122 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1124 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1126 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1133 /* Better use sensible values than insane (0.0) ones... */
1134 ir->gb_epsilon_solvent = 80;
1135 ir->gb_obc_alpha = 1.0;
1136 ir->gb_obc_beta = 0.8;
1137 ir->gb_obc_gamma = 4.85;
1138 ir->sa_surface_tension = 2.092;
1142 if (file_version >= 81)
1144 gmx_fio_do_real(fio, ir->fourier_spacing);
1148 ir->fourier_spacing = 0.0;
1150 gmx_fio_do_int(fio, ir->nkx);
1151 gmx_fio_do_int(fio, ir->nky);
1152 gmx_fio_do_int(fio, ir->nkz);
1153 gmx_fio_do_int(fio, ir->pme_order);
1154 gmx_fio_do_real(fio, ir->ewald_rtol);
1156 if (file_version >= 93)
1158 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1162 ir->ewald_rtol_lj = ir->ewald_rtol;
1165 if (file_version >= 24)
1167 gmx_fio_do_int(fio, ir->ewald_geometry);
1171 ir->ewald_geometry = eewg3D;
1174 if (file_version <= 17)
1176 ir->epsilon_surface = 0;
1177 if (file_version == 17)
1179 gmx_fio_do_int(fio, idum);
1184 gmx_fio_do_real(fio, ir->epsilon_surface);
1187 gmx_fio_do_gmx_bool(fio, ir->bOptFFT);
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;
1391 if (file_version >= 26 && file_version < 79)
1393 gmx_fio_do_real(fio, ir->dihre_fc);
1394 if (file_version < 56)
1396 gmx_fio_do_real(fio, rdum);
1397 gmx_fio_do_int(fio, idum);
1405 gmx_fio_do_real(fio, ir->em_stepsize);
1406 gmx_fio_do_real(fio, ir->em_tol);
1407 if (file_version >= 22)
1409 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1413 ir->bShakeSOR = TRUE;
1415 if (file_version >= 11)
1417 gmx_fio_do_int(fio, ir->niter);
1422 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1425 if (file_version >= 21)
1427 gmx_fio_do_real(fio, ir->fc_stepsize);
1431 ir->fc_stepsize = 0;
1433 gmx_fio_do_int(fio, ir->eConstrAlg);
1434 gmx_fio_do_int(fio, ir->nProjOrder);
1435 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1436 if (file_version <= 14)
1438 gmx_fio_do_int(fio, idum);
1440 if (file_version >= 26)
1442 gmx_fio_do_int(fio, ir->nLincsIter);
1447 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1450 if (file_version < 33)
1452 gmx_fio_do_real(fio, bd_temp);
1454 gmx_fio_do_real(fio, ir->bd_fric);
1455 if (file_version >= tpxv_Use64BitRandomSeed)
1457 gmx_fio_do_int64(fio, ir->ld_seed);
1461 gmx_fio_do_int(fio, idum);
1464 if (file_version >= 33)
1466 for (i = 0; i < DIM; i++)
1468 gmx_fio_do_rvec(fio, ir->deform[i]);
1473 for (i = 0; i < DIM; i++)
1475 clear_rvec(ir->deform[i]);
1478 if (file_version >= 14)
1480 gmx_fio_do_real(fio, ir->cos_accel);
1486 gmx_fio_do_int(fio, ir->userint1);
1487 gmx_fio_do_int(fio, ir->userint2);
1488 gmx_fio_do_int(fio, ir->userint3);
1489 gmx_fio_do_int(fio, ir->userint4);
1490 gmx_fio_do_real(fio, ir->userreal1);
1491 gmx_fio_do_real(fio, ir->userreal2);
1492 gmx_fio_do_real(fio, ir->userreal3);
1493 gmx_fio_do_real(fio, ir->userreal4);
1496 if (file_version >= 77)
1498 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1503 snew(ir->adress, 1);
1505 gmx_fio_do_int(fio, ir->adress->type);
1506 gmx_fio_do_real(fio, ir->adress->const_wf);
1507 gmx_fio_do_real(fio, ir->adress->ex_width);
1508 gmx_fio_do_real(fio, ir->adress->hy_width);
1509 gmx_fio_do_int(fio, ir->adress->icor);
1510 gmx_fio_do_int(fio, ir->adress->site);
1511 gmx_fio_do_rvec(fio, ir->adress->refs);
1512 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1513 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1514 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1515 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1519 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1521 if (ir->adress->n_tf_grps > 0)
1523 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1527 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1529 if (ir->adress->n_energy_grps > 0)
1531 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1537 ir->bAdress = FALSE;
1541 if (file_version >= 48)
1543 gmx_fio_do_int(fio, ir->ePull);
1544 if (ir->ePull != epullNO)
1550 do_pull(fio, ir->pull, bRead, file_version);
1555 ir->ePull = epullNO;
1558 /* Enforced rotation */
1559 if (file_version >= 74)
1561 gmx_fio_do_int(fio, ir->bRot);
1562 if (ir->bRot == TRUE)
1568 do_rot(fio, ir->rot, bRead);
1576 /* Interactive molecular dynamics */
1577 if (file_version >= tpxv_InteractiveMolecularDynamics)
1579 gmx_fio_do_int(fio, ir->bIMD);
1580 if (TRUE == ir->bIMD)
1586 do_imd(fio, ir->imd, bRead);
1591 /* We don't support IMD sessions for old .tpr files */
1596 gmx_fio_do_int(fio, ir->opts.ngtc);
1597 if (file_version >= 69)
1599 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1603 ir->opts.nhchainlength = 1;
1605 gmx_fio_do_int(fio, ir->opts.ngacc);
1606 gmx_fio_do_int(fio, ir->opts.ngfrz);
1607 gmx_fio_do_int(fio, ir->opts.ngener);
1611 snew(ir->opts.nrdf, ir->opts.ngtc);
1612 snew(ir->opts.ref_t, ir->opts.ngtc);
1613 snew(ir->opts.annealing, ir->opts.ngtc);
1614 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1615 snew(ir->opts.anneal_time, ir->opts.ngtc);
1616 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1617 snew(ir->opts.tau_t, ir->opts.ngtc);
1618 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1619 snew(ir->opts.acc, ir->opts.ngacc);
1620 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1622 if (ir->opts.ngtc > 0)
1624 if (bRead && file_version < 13)
1626 snew(tmp, ir->opts.ngtc);
1627 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1628 for (i = 0; i < ir->opts.ngtc; i++)
1630 ir->opts.nrdf[i] = tmp[i];
1636 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1638 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1639 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1640 if (file_version < 33 && ir->eI == eiBD)
1642 for (i = 0; i < ir->opts.ngtc; i++)
1644 ir->opts.tau_t[i] = bd_temp;
1648 if (ir->opts.ngfrz > 0)
1650 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1652 if (ir->opts.ngacc > 0)
1654 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1656 if (file_version >= 12)
1658 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1659 ir->opts.ngener*ir->opts.ngener);
1662 if (bRead && file_version < 26)
1664 for (i = 0; i < ir->opts.ngtc; i++)
1668 ir->opts.annealing[i] = eannSINGLE;
1669 ir->opts.anneal_npoints[i] = 2;
1670 snew(ir->opts.anneal_time[i], 2);
1671 snew(ir->opts.anneal_temp[i], 2);
1672 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1673 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1674 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1675 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1676 ir->opts.anneal_time[i][0] = ir->init_t;
1677 ir->opts.anneal_time[i][1] = finish_t;
1678 ir->opts.anneal_temp[i][0] = init_temp;
1679 ir->opts.anneal_temp[i][1] = finish_temp;
1683 ir->opts.annealing[i] = eannNO;
1684 ir->opts.anneal_npoints[i] = 0;
1690 /* file version 26 or later */
1691 /* First read the lists with annealing and npoints for each group */
1692 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1693 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1694 for (j = 0; j < (ir->opts.ngtc); j++)
1696 k = ir->opts.anneal_npoints[j];
1699 snew(ir->opts.anneal_time[j], k);
1700 snew(ir->opts.anneal_temp[j], k);
1702 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1703 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1707 if (file_version >= 45)
1709 gmx_fio_do_int(fio, ir->nwall);
1710 gmx_fio_do_int(fio, ir->wall_type);
1711 if (file_version >= 50)
1713 gmx_fio_do_real(fio, ir->wall_r_linpot);
1717 ir->wall_r_linpot = -1;
1719 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1720 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1721 gmx_fio_do_real(fio, ir->wall_density[0]);
1722 gmx_fio_do_real(fio, ir->wall_density[1]);
1723 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1729 ir->wall_atomtype[0] = -1;
1730 ir->wall_atomtype[1] = -1;
1731 ir->wall_density[0] = 0;
1732 ir->wall_density[1] = 0;
1733 ir->wall_ewald_zfac = 3;
1735 /* Cosine stuff for electric fields */
1736 for (j = 0; (j < DIM); j++)
1738 gmx_fio_do_int(fio, ir->ex[j].n);
1739 gmx_fio_do_int(fio, ir->et[j].n);
1742 snew(ir->ex[j].a, ir->ex[j].n);
1743 snew(ir->ex[j].phi, ir->ex[j].n);
1744 snew(ir->et[j].a, ir->et[j].n);
1745 snew(ir->et[j].phi, ir->et[j].n);
1747 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1748 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1749 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1750 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1754 if (file_version >= tpxv_ComputationalElectrophysiology)
1756 gmx_fio_do_int(fio, ir->eSwapCoords);
1757 if (ir->eSwapCoords != eswapNO)
1763 do_swapcoords(fio, ir->swap, bRead);
1768 if (file_version >= 39)
1770 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1771 gmx_fio_do_int(fio, ir->QMMMscheme);
1772 gmx_fio_do_real(fio, ir->scalefactor);
1773 gmx_fio_do_int(fio, ir->opts.ngQM);
1776 snew(ir->opts.QMmethod, ir->opts.ngQM);
1777 snew(ir->opts.QMbasis, ir->opts.ngQM);
1778 snew(ir->opts.QMcharge, ir->opts.ngQM);
1779 snew(ir->opts.QMmult, ir->opts.ngQM);
1780 snew(ir->opts.bSH, ir->opts.ngQM);
1781 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1782 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1783 snew(ir->opts.SAon, ir->opts.ngQM);
1784 snew(ir->opts.SAoff, ir->opts.ngQM);
1785 snew(ir->opts.SAsteps, ir->opts.ngQM);
1786 snew(ir->opts.bOPT, ir->opts.ngQM);
1787 snew(ir->opts.bTS, ir->opts.ngQM);
1789 if (ir->opts.ngQM > 0)
1791 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1792 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1793 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1794 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1795 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1796 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1797 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1798 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1799 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1800 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1801 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1802 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1804 /* end of QMMM stuff */
1809 static void do_harm(t_fileio *fio, t_iparams *iparams)
1811 gmx_fio_do_real(fio, iparams->harmonic.rA);
1812 gmx_fio_do_real(fio, iparams->harmonic.krA);
1813 gmx_fio_do_real(fio, iparams->harmonic.rB);
1814 gmx_fio_do_real(fio, iparams->harmonic.krB);
1817 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1818 gmx_bool bRead, int file_version)
1825 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1835 do_harm(fio, iparams);
1836 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1838 /* Correct incorrect storage of parameters */
1839 iparams->pdihs.phiB = iparams->pdihs.phiA;
1840 iparams->pdihs.cpB = iparams->pdihs.cpA;
1844 gmx_fio_do_real(fio, iparams->harmonic.rA);
1845 gmx_fio_do_real(fio, iparams->harmonic.krA);
1847 case F_LINEAR_ANGLES:
1848 gmx_fio_do_real(fio, iparams->linangle.klinA);
1849 gmx_fio_do_real(fio, iparams->linangle.aA);
1850 gmx_fio_do_real(fio, iparams->linangle.klinB);
1851 gmx_fio_do_real(fio, iparams->linangle.aB);
1854 gmx_fio_do_real(fio, iparams->fene.bm);
1855 gmx_fio_do_real(fio, iparams->fene.kb);
1859 gmx_fio_do_real(fio, iparams->restraint.lowA);
1860 gmx_fio_do_real(fio, iparams->restraint.up1A);
1861 gmx_fio_do_real(fio, iparams->restraint.up2A);
1862 gmx_fio_do_real(fio, iparams->restraint.kA);
1863 gmx_fio_do_real(fio, iparams->restraint.lowB);
1864 gmx_fio_do_real(fio, iparams->restraint.up1B);
1865 gmx_fio_do_real(fio, iparams->restraint.up2B);
1866 gmx_fio_do_real(fio, iparams->restraint.kB);
1872 gmx_fio_do_real(fio, iparams->tab.kA);
1873 gmx_fio_do_int(fio, iparams->tab.table);
1874 gmx_fio_do_real(fio, iparams->tab.kB);
1876 case F_CROSS_BOND_BONDS:
1877 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1878 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1879 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1881 case F_CROSS_BOND_ANGLES:
1882 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1883 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1884 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1885 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1887 case F_UREY_BRADLEY:
1888 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1889 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1890 gmx_fio_do_real(fio, iparams->u_b.r13A);
1891 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1892 if (file_version >= 79)
1894 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1895 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1896 gmx_fio_do_real(fio, iparams->u_b.r13B);
1897 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1901 iparams->u_b.thetaB = iparams->u_b.thetaA;
1902 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1903 iparams->u_b.r13B = iparams->u_b.r13A;
1904 iparams->u_b.kUBB = iparams->u_b.kUBA;
1907 case F_QUARTIC_ANGLES:
1908 gmx_fio_do_real(fio, iparams->qangle.theta);
1909 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1912 gmx_fio_do_real(fio, iparams->bham.a);
1913 gmx_fio_do_real(fio, iparams->bham.b);
1914 gmx_fio_do_real(fio, iparams->bham.c);
1917 gmx_fio_do_real(fio, iparams->morse.b0A);
1918 gmx_fio_do_real(fio, iparams->morse.cbA);
1919 gmx_fio_do_real(fio, iparams->morse.betaA);
1920 if (file_version >= 79)
1922 gmx_fio_do_real(fio, iparams->morse.b0B);
1923 gmx_fio_do_real(fio, iparams->morse.cbB);
1924 gmx_fio_do_real(fio, iparams->morse.betaB);
1928 iparams->morse.b0B = iparams->morse.b0A;
1929 iparams->morse.cbB = iparams->morse.cbA;
1930 iparams->morse.betaB = iparams->morse.betaA;
1934 gmx_fio_do_real(fio, iparams->cubic.b0);
1935 gmx_fio_do_real(fio, iparams->cubic.kb);
1936 gmx_fio_do_real(fio, iparams->cubic.kcub);
1940 case F_POLARIZATION:
1941 gmx_fio_do_real(fio, iparams->polarize.alpha);
1944 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1945 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1946 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1949 if (file_version < 31)
1951 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1953 gmx_fio_do_real(fio, iparams->wpol.al_x);
1954 gmx_fio_do_real(fio, iparams->wpol.al_y);
1955 gmx_fio_do_real(fio, iparams->wpol.al_z);
1956 gmx_fio_do_real(fio, iparams->wpol.rOH);
1957 gmx_fio_do_real(fio, iparams->wpol.rHH);
1958 gmx_fio_do_real(fio, iparams->wpol.rOD);
1961 gmx_fio_do_real(fio, iparams->thole.a);
1962 gmx_fio_do_real(fio, iparams->thole.alpha1);
1963 gmx_fio_do_real(fio, iparams->thole.alpha2);
1964 gmx_fio_do_real(fio, iparams->thole.rfac);
1967 gmx_fio_do_real(fio, iparams->lj.c6);
1968 gmx_fio_do_real(fio, iparams->lj.c12);
1971 gmx_fio_do_real(fio, iparams->lj14.c6A);
1972 gmx_fio_do_real(fio, iparams->lj14.c12A);
1973 gmx_fio_do_real(fio, iparams->lj14.c6B);
1974 gmx_fio_do_real(fio, iparams->lj14.c12B);
1977 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1978 gmx_fio_do_real(fio, iparams->ljc14.qi);
1979 gmx_fio_do_real(fio, iparams->ljc14.qj);
1980 gmx_fio_do_real(fio, iparams->ljc14.c6);
1981 gmx_fio_do_real(fio, iparams->ljc14.c12);
1983 case F_LJC_PAIRS_NB:
1984 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1985 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1986 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1987 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1993 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1994 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1995 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1997 /* Read the incorrectly stored multiplicity */
1998 gmx_fio_do_real(fio, iparams->harmonic.rB);
1999 gmx_fio_do_real(fio, iparams->harmonic.krB);
2000 iparams->pdihs.phiB = iparams->pdihs.phiA;
2001 iparams->pdihs.cpB = iparams->pdihs.cpA;
2005 gmx_fio_do_real(fio, iparams->pdihs.phiB);
2006 gmx_fio_do_real(fio, iparams->pdihs.cpB);
2007 gmx_fio_do_int(fio, iparams->pdihs.mult);
2011 gmx_fio_do_real(fio, iparams->pdihs.phiA);
2012 gmx_fio_do_real(fio, iparams->pdihs.cpA);
2015 gmx_fio_do_int(fio, iparams->disres.label);
2016 gmx_fio_do_int(fio, iparams->disres.type);
2017 gmx_fio_do_real(fio, iparams->disres.low);
2018 gmx_fio_do_real(fio, iparams->disres.up1);
2019 gmx_fio_do_real(fio, iparams->disres.up2);
2020 gmx_fio_do_real(fio, iparams->disres.kfac);
2023 gmx_fio_do_int(fio, iparams->orires.ex);
2024 gmx_fio_do_int(fio, iparams->orires.label);
2025 gmx_fio_do_int(fio, iparams->orires.power);
2026 gmx_fio_do_real(fio, iparams->orires.c);
2027 gmx_fio_do_real(fio, iparams->orires.obs);
2028 gmx_fio_do_real(fio, iparams->orires.kfac);
2031 if (file_version < 82)
2033 gmx_fio_do_int(fio, idum);
2034 gmx_fio_do_int(fio, idum);
2036 gmx_fio_do_real(fio, iparams->dihres.phiA);
2037 gmx_fio_do_real(fio, iparams->dihres.dphiA);
2038 gmx_fio_do_real(fio, iparams->dihres.kfacA);
2039 if (file_version >= 82)
2041 gmx_fio_do_real(fio, iparams->dihres.phiB);
2042 gmx_fio_do_real(fio, iparams->dihres.dphiB);
2043 gmx_fio_do_real(fio, iparams->dihres.kfacB);
2047 iparams->dihres.phiB = iparams->dihres.phiA;
2048 iparams->dihres.dphiB = iparams->dihres.dphiA;
2049 iparams->dihres.kfacB = iparams->dihres.kfacA;
2053 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
2054 gmx_fio_do_rvec(fio, iparams->posres.fcA);
2055 if (bRead && file_version < 27)
2057 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
2058 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
2062 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
2063 gmx_fio_do_rvec(fio, iparams->posres.fcB);
2067 gmx_fio_do_int(fio, iparams->fbposres.geom);
2068 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
2069 gmx_fio_do_real(fio, iparams->fbposres.r);
2070 gmx_fio_do_real(fio, iparams->fbposres.k);
2073 gmx_fio_ndo_real(fio, iparams->cbtdihs.cbtcA, NR_CBTDIHS);
2076 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2077 if (file_version >= 25)
2079 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2083 /* Fourier dihedrals are internally represented
2084 * as Ryckaert-Bellemans since those are faster to compute.
2086 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2087 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2091 gmx_fio_do_real(fio, iparams->constr.dA);
2092 gmx_fio_do_real(fio, iparams->constr.dB);
2095 gmx_fio_do_real(fio, iparams->settle.doh);
2096 gmx_fio_do_real(fio, iparams->settle.dhh);
2099 gmx_fio_do_real(fio, iparams->vsite.a);
2104 gmx_fio_do_real(fio, iparams->vsite.a);
2105 gmx_fio_do_real(fio, iparams->vsite.b);
2110 gmx_fio_do_real(fio, iparams->vsite.a);
2111 gmx_fio_do_real(fio, iparams->vsite.b);
2112 gmx_fio_do_real(fio, iparams->vsite.c);
2115 gmx_fio_do_int(fio, iparams->vsiten.n);
2116 gmx_fio_do_real(fio, iparams->vsiten.a);
2121 /* We got rid of some parameters in version 68 */
2122 if (bRead && file_version < 68)
2124 gmx_fio_do_real(fio, rdum);
2125 gmx_fio_do_real(fio, rdum);
2126 gmx_fio_do_real(fio, rdum);
2127 gmx_fio_do_real(fio, rdum);
2129 gmx_fio_do_real(fio, iparams->gb.sar);
2130 gmx_fio_do_real(fio, iparams->gb.st);
2131 gmx_fio_do_real(fio, iparams->gb.pi);
2132 gmx_fio_do_real(fio, iparams->gb.gbr);
2133 gmx_fio_do_real(fio, iparams->gb.bmlt);
2136 gmx_fio_do_int(fio, iparams->cmap.cmapA);
2137 gmx_fio_do_int(fio, iparams->cmap.cmapB);
2140 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
2141 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
2145 gmx_fio_unset_comment(fio);
2149 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
2156 gmx_fio_set_comment(fio, interaction_function[ftype].name);
2158 if (file_version < 44)
2160 for (i = 0; i < MAXNODES; i++)
2162 gmx_fio_do_int(fio, idum);
2165 gmx_fio_do_int(fio, ilist->nr);
2168 snew(ilist->iatoms, ilist->nr);
2170 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2173 gmx_fio_unset_comment(fio);
2177 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2178 gmx_bool bRead, int file_version)
2183 gmx_fio_do_int(fio, ffparams->atnr);
2184 if (file_version < 57)
2186 gmx_fio_do_int(fio, idum);
2188 gmx_fio_do_int(fio, ffparams->ntypes);
2191 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2192 ffparams->atnr, ffparams->ntypes);
2196 snew(ffparams->functype, ffparams->ntypes);
2197 snew(ffparams->iparams, ffparams->ntypes);
2199 /* Read/write all the function types */
2200 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2203 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2206 if (file_version >= 66)
2208 gmx_fio_do_double(fio, ffparams->reppow);
2212 ffparams->reppow = 12.0;
2215 if (file_version >= 57)
2217 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2220 /* Check whether all these function types are supported by the code.
2221 * In practice the code is backwards compatible, which means that the
2222 * numbering may have to be altered from old numbering to new numbering
2224 for (i = 0; (i < ffparams->ntypes); i++)
2228 /* Loop over file versions */
2229 for (k = 0; (k < NFTUPD); k++)
2231 /* Compare the read file_version to the update table */
2232 if ((file_version < ftupd[k].fvnr) &&
2233 (ffparams->functype[i] >= ftupd[k].ftype))
2235 ffparams->functype[i] += 1;
2238 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2239 i, ffparams->functype[i],
2240 interaction_function[ftupd[k].ftype].longname);
2247 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2251 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2256 static void add_settle_atoms(t_ilist *ilist)
2260 /* Settle used to only store the first atom: add the other two */
2261 srenew(ilist->iatoms, 2*ilist->nr);
2262 for (i = ilist->nr/2-1; i >= 0; i--)
2264 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2265 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2266 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2267 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2269 ilist->nr = 2*ilist->nr;
2272 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2275 int i, j, renum[F_NRE];
2279 for (j = 0; (j < F_NRE); j++)
2284 for (k = 0; k < NFTUPD; k++)
2286 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2295 ilist[j].iatoms = NULL;
2299 do_ilist(fio, &ilist[j], bRead, file_version, j);
2300 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2302 add_settle_atoms(&ilist[j]);
2306 if (bRead && gmx_debug_at)
2307 pr_ilist(debug,0,interaction_function[j].longname,
2308 functype,&ilist[j],TRUE);
2313 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2314 gmx_bool bRead, int file_version)
2316 do_ffparams(fio, ffparams, bRead, file_version);
2318 if (file_version >= 54)
2320 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2323 do_ilists(fio, molt->ilist, bRead, file_version);
2326 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2328 int i, idum, dum_nra, *dum_a;
2330 if (file_version < 44)
2332 for (i = 0; i < MAXNODES; i++)
2334 gmx_fio_do_int(fio, idum);
2337 gmx_fio_do_int(fio, block->nr);
2338 if (file_version < 51)
2340 gmx_fio_do_int(fio, dum_nra);
2344 if ((block->nalloc_index > 0) && (NULL != block->index))
2346 sfree(block->index);
2348 block->nalloc_index = block->nr+1;
2349 snew(block->index, block->nalloc_index);
2351 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2353 if (file_version < 51 && dum_nra > 0)
2355 snew(dum_a, dum_nra);
2356 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2361 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2366 if (file_version < 44)
2368 for (i = 0; i < MAXNODES; i++)
2370 gmx_fio_do_int(fio, idum);
2373 gmx_fio_do_int(fio, block->nr);
2374 gmx_fio_do_int(fio, block->nra);
2377 block->nalloc_index = block->nr+1;
2378 snew(block->index, block->nalloc_index);
2379 block->nalloc_a = block->nra;
2380 snew(block->a, block->nalloc_a);
2382 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2383 gmx_fio_ndo_int(fio, block->a, block->nra);
2386 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2387 int file_version, gmx_groups_t *groups, int atnr)
2391 gmx_fio_do_real(fio, atom->m);
2392 gmx_fio_do_real(fio, atom->q);
2393 gmx_fio_do_real(fio, atom->mB);
2394 gmx_fio_do_real(fio, atom->qB);
2395 gmx_fio_do_ushort(fio, atom->type);
2396 gmx_fio_do_ushort(fio, atom->typeB);
2397 gmx_fio_do_int(fio, atom->ptype);
2398 gmx_fio_do_int(fio, atom->resind);
2399 if (file_version >= 52)
2401 gmx_fio_do_int(fio, atom->atomnumber);
2405 atom->atomnumber = NOTSET;
2407 if (file_version < 23)
2411 else if (file_version < 39)
2420 if (file_version < 57)
2422 unsigned char uchar[egcNR];
2423 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2424 for (i = myngrp; (i < ngrp); i++)
2428 /* Copy the old data format to the groups struct */
2429 for (i = 0; i < ngrp; i++)
2431 groups->grpnr[i][atnr] = uchar[i];
2436 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2441 if (file_version < 23)
2445 else if (file_version < 39)
2454 for (j = 0; (j < ngrp); j++)
2458 gmx_fio_do_int(fio, grps[j].nr);
2461 snew(grps[j].nm_ind, grps[j].nr);
2463 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2468 snew(grps[j].nm_ind, grps[j].nr);
2473 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2479 gmx_fio_do_int(fio, ls);
2480 *nm = get_symtab_handle(symtab, ls);
2484 ls = lookup_symtab(symtab, *nm);
2485 gmx_fio_do_int(fio, ls);
2489 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2494 for (j = 0; (j < nstr); j++)
2496 do_symstr(fio, &(nm[j]), bRead, symtab);
2500 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2501 t_symtab *symtab, int file_version)
2505 for (j = 0; (j < n); j++)
2507 do_symstr(fio, &(ri[j].name), bRead, symtab);
2508 if (file_version >= 63)
2510 gmx_fio_do_int(fio, ri[j].nr);
2511 gmx_fio_do_uchar(fio, ri[j].ic);
2521 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2523 gmx_groups_t *groups)
2527 gmx_fio_do_int(fio, atoms->nr);
2528 gmx_fio_do_int(fio, atoms->nres);
2529 if (file_version < 57)
2531 gmx_fio_do_int(fio, groups->ngrpname);
2532 for (i = 0; i < egcNR; i++)
2534 groups->ngrpnr[i] = atoms->nr;
2535 snew(groups->grpnr[i], groups->ngrpnr[i]);
2540 snew(atoms->atom, atoms->nr);
2541 snew(atoms->atomname, atoms->nr);
2542 snew(atoms->atomtype, atoms->nr);
2543 snew(atoms->atomtypeB, atoms->nr);
2544 snew(atoms->resinfo, atoms->nres);
2545 if (file_version < 57)
2547 snew(groups->grpname, groups->ngrpname);
2549 atoms->pdbinfo = NULL;
2551 for (i = 0; (i < atoms->nr); i++)
2553 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2555 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2556 if (bRead && (file_version <= 20))
2558 for (i = 0; i < atoms->nr; i++)
2560 atoms->atomtype[i] = put_symtab(symtab, "?");
2561 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2566 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2567 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2569 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2571 if (file_version < 57)
2573 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2575 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2579 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2580 gmx_bool bRead, t_symtab *symtab,
2585 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2586 gmx_fio_do_int(fio, groups->ngrpname);
2589 snew(groups->grpname, groups->ngrpname);
2591 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2592 for (g = 0; g < egcNR; g++)
2594 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2595 if (groups->ngrpnr[g] == 0)
2599 groups->grpnr[g] = NULL;
2606 snew(groups->grpnr[g], groups->ngrpnr[g]);
2608 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2613 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2618 if (file_version > 25)
2620 gmx_fio_do_int(fio, atomtypes->nr);
2624 snew(atomtypes->radius, j);
2625 snew(atomtypes->vol, j);
2626 snew(atomtypes->surftens, j);
2627 snew(atomtypes->atomnumber, j);
2628 snew(atomtypes->gb_radius, j);
2629 snew(atomtypes->S_hct, j);
2631 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2632 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2633 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2634 if (file_version >= 40)
2636 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2638 if (file_version >= 60)
2640 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2641 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2646 /* File versions prior to 26 cannot do GBSA,
2647 * so they dont use this structure
2650 atomtypes->radius = NULL;
2651 atomtypes->vol = NULL;
2652 atomtypes->surftens = NULL;
2653 atomtypes->atomnumber = NULL;
2654 atomtypes->gb_radius = NULL;
2655 atomtypes->S_hct = NULL;
2659 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2665 gmx_fio_do_int(fio, symtab->nr);
2669 snew(symtab->symbuf, 1);
2670 symbuf = symtab->symbuf;
2671 symbuf->bufsize = nr;
2672 snew(symbuf->buf, nr);
2673 for (i = 0; (i < nr); i++)
2675 gmx_fio_do_string(fio, buf);
2676 symbuf->buf[i] = strdup(buf);
2681 symbuf = symtab->symbuf;
2682 while (symbuf != NULL)
2684 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2686 gmx_fio_do_string(fio, symbuf->buf[i]);
2689 symbuf = symbuf->next;
2693 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2698 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2700 int i, j, ngrid, gs, nelem;
2702 gmx_fio_do_int(fio, cmap_grid->ngrid);
2703 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2705 ngrid = cmap_grid->ngrid;
2706 gs = cmap_grid->grid_spacing;
2711 snew(cmap_grid->cmapdata, ngrid);
2713 for (i = 0; i < cmap_grid->ngrid; i++)
2715 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2719 for (i = 0; i < cmap_grid->ngrid; i++)
2721 for (j = 0; j < nelem; j++)
2723 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2724 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2725 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2726 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2732 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2734 int m, a, a0, a1, r;
2738 /* We always assign a new chain number, but save the chain id characters
2739 * for larger molecules.
2741 #define CHAIN_MIN_ATOMS 15
2745 for (m = 0; m < mols->nr; m++)
2747 a0 = mols->index[m];
2748 a1 = mols->index[m+1];
2749 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2758 for (a = a0; a < a1; a++)
2760 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2761 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2766 /* Blank out the chain id if there was only one chain */
2769 for (r = 0; r < atoms->nres; r++)
2771 atoms->resinfo[r].chainid = ' ';
2776 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2777 t_symtab *symtab, int file_version,
2778 gmx_groups_t *groups)
2782 if (file_version >= 57)
2784 do_symstr(fio, &(molt->name), bRead, symtab);
2787 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2789 if (bRead && gmx_debug_at)
2791 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2794 if (file_version >= 57)
2796 do_ilists(fio, molt->ilist, bRead, file_version);
2798 do_block(fio, &molt->cgs, bRead, file_version);
2799 if (bRead && gmx_debug_at)
2801 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2805 /* This used to be in the atoms struct */
2806 do_blocka(fio, &molt->excls, bRead, file_version);
2809 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2813 gmx_fio_do_int(fio, molb->type);
2814 gmx_fio_do_int(fio, molb->nmol);
2815 gmx_fio_do_int(fio, molb->natoms_mol);
2816 /* Position restraint coordinates */
2817 gmx_fio_do_int(fio, molb->nposres_xA);
2818 if (molb->nposres_xA > 0)
2822 snew(molb->posres_xA, molb->nposres_xA);
2824 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2826 gmx_fio_do_int(fio, molb->nposres_xB);
2827 if (molb->nposres_xB > 0)
2831 snew(molb->posres_xB, molb->nposres_xB);
2833 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2838 static t_block mtop_mols(gmx_mtop_t *mtop)
2844 for (mb = 0; mb < mtop->nmolblock; mb++)
2846 mols.nr += mtop->molblock[mb].nmol;
2848 mols.nalloc_index = mols.nr + 1;
2849 snew(mols.index, mols.nalloc_index);
2854 for (mb = 0; mb < mtop->nmolblock; mb++)
2856 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2858 a += mtop->molblock[mb].natoms_mol;
2867 static void add_posres_molblock(gmx_mtop_t *mtop)
2872 gmx_molblock_t *molb;
2875 /* posres reference positions are stored in ip->posres (if present) and
2876 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2877 posres.pos0A are identical to fbposres.pos0. */
2878 il = &mtop->moltype[0].ilist[F_POSRES];
2879 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2880 if (il->nr == 0 && ilfb->nr == 0)
2886 for (i = 0; i < il->nr; i += 2)
2888 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2889 am = max(am, il->iatoms[i+1]);
2890 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2891 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2892 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2897 /* This loop is required if we have only flat-bottomed posres:
2899 - bFE == FALSE (no B-state for flat-bottomed posres) */
2902 for (i = 0; i < ilfb->nr; i += 2)
2904 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2905 am = max(am, ilfb->iatoms[i+1]);
2908 /* Make the posres coordinate block end at a molecule end */
2910 while (am >= mtop->mols.index[mol+1])
2914 molb = &mtop->molblock[0];
2915 molb->nposres_xA = mtop->mols.index[mol+1];
2916 snew(molb->posres_xA, molb->nposres_xA);
2919 molb->nposres_xB = molb->nposres_xA;
2920 snew(molb->posres_xB, molb->nposres_xB);
2924 molb->nposres_xB = 0;
2926 for (i = 0; i < il->nr; i += 2)
2928 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2929 a = il->iatoms[i+1];
2930 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2931 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2932 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2935 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2936 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2937 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2942 /* If only flat-bottomed posres are present, take reference pos from them.
2943 Here: bFE == FALSE */
2944 for (i = 0; i < ilfb->nr; i += 2)
2946 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2947 a = ilfb->iatoms[i+1];
2948 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2949 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2950 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2955 static void set_disres_npair(gmx_mtop_t *mtop)
2962 ip = mtop->ffparams.iparams;
2964 for (mt = 0; mt < mtop->nmoltype; mt++)
2966 il = &mtop->moltype[mt].ilist[F_DISRES];
2971 for (i = 0; i < il->nr; i += 3)
2974 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2976 ip[a[i]].disres.npair = npair;
2984 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2994 do_symtab(fio, &(mtop->symtab), bRead);
2997 pr_symtab(debug, 0, "symtab", &mtop->symtab);
3000 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
3002 if (file_version >= 57)
3004 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
3006 gmx_fio_do_int(fio, mtop->nmoltype);
3014 snew(mtop->moltype, mtop->nmoltype);
3015 if (file_version < 57)
3017 mtop->moltype[0].name = mtop->name;
3020 for (mt = 0; mt < mtop->nmoltype; mt++)
3022 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
3026 if (file_version >= 57)
3028 gmx_fio_do_int(fio, mtop->nmolblock);
3032 mtop->nmolblock = 1;
3036 snew(mtop->molblock, mtop->nmolblock);
3038 if (file_version >= 57)
3040 for (mb = 0; mb < mtop->nmolblock; mb++)
3042 do_molblock(fio, &mtop->molblock[mb], bRead);
3044 gmx_fio_do_int(fio, mtop->natoms);
3048 mtop->molblock[0].type = 0;
3049 mtop->molblock[0].nmol = 1;
3050 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
3051 mtop->molblock[0].nposres_xA = 0;
3052 mtop->molblock[0].nposres_xB = 0;
3055 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
3058 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
3061 if (file_version < 57)
3063 /* Debug statements are inside do_idef */
3064 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
3065 mtop->natoms = mtop->moltype[0].atoms.nr;
3068 if (file_version >= 65)
3070 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
3074 mtop->ffparams.cmap_grid.ngrid = 0;
3075 mtop->ffparams.cmap_grid.grid_spacing = 0;
3076 mtop->ffparams.cmap_grid.cmapdata = NULL;
3079 if (file_version >= 57)
3081 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
3084 if (file_version < 57)
3086 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
3087 if (bRead && gmx_debug_at)
3089 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
3091 do_block(fio, &mtop->mols, bRead, file_version);
3092 /* Add the posres coordinates to the molblock */
3093 add_posres_molblock(mtop);
3097 if (file_version >= 57)
3099 done_block(&mtop->mols);
3100 mtop->mols = mtop_mols(mtop);
3104 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
3108 if (file_version < 51)
3110 /* Here used to be the shake blocks */
3111 do_blocka(fio, &dumb, bRead, file_version);
3124 close_symtab(&(mtop->symtab));
3128 /* If TopOnlyOK is TRUE then we can read even future versions
3129 * of tpx files, provided the file_generation hasn't changed.
3130 * If it is FALSE, we need the inputrecord too, and bail out
3131 * if the file is newer than the program.
3133 * The version and generation if the topology (see top of this file)
3134 * are returned in the two last arguments.
3136 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
3138 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
3139 gmx_bool TopOnlyOK, int *file_version,
3140 int *file_generation)
3143 char file_tag[STRLEN];
3150 gmx_fio_checktype(fio);
3151 gmx_fio_setdebug(fio, bDebugMode());
3153 /* NEW! XDR tpb file */
3154 precision = sizeof(real);
3157 gmx_fio_do_string(fio, buf);
3158 if (strncmp(buf, "VERSION", 7))
3160 gmx_fatal(FARGS, "Can not read file %s,\n"
3161 " this file is from a Gromacs version which is older than 2.0\n"
3162 " Make a new one with grompp or use a gro or pdb file, if possible",
3163 gmx_fio_getname(fio));
3165 gmx_fio_do_int(fio, precision);
3166 bDouble = (precision == sizeof(double));
3167 if ((precision != sizeof(float)) && !bDouble)
3169 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3170 "instead of %d or %d",
3171 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3173 gmx_fio_setprecision(fio, bDouble);
3174 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3175 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3179 gmx_fio_write_string(fio, GromacsVersion());
3180 bDouble = (precision == sizeof(double));
3181 gmx_fio_setprecision(fio, bDouble);
3182 gmx_fio_do_int(fio, precision);
3184 sprintf(file_tag, "%s", tpx_tag);
3185 fgen = tpx_generation;
3188 /* Check versions! */
3189 gmx_fio_do_int(fio, fver);
3191 /* This is for backward compatibility with development versions 77-79
3192 * where the tag was, mistakenly, placed before the generation,
3193 * which would cause a segv instead of a proper error message
3194 * when reading the topology only from tpx with <77 code.
3196 if (fver >= 77 && fver <= 79)
3198 gmx_fio_do_string(fio, file_tag);
3203 gmx_fio_do_int(fio, fgen);
3212 gmx_fio_do_string(fio, file_tag);
3218 /* Versions before 77 don't have the tag, set it to release */
3219 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3222 if (strcmp(file_tag, tpx_tag) != 0)
3224 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3227 /* We only support reading tpx files with the same tag as the code
3228 * or tpx files with the release tag and with lower version number.
3230 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
3232 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3233 gmx_fio_getname(fio), fver, file_tag,
3234 tpx_version, tpx_tag);
3239 if (file_version != NULL)
3241 *file_version = fver;
3243 if (file_generation != NULL)
3245 *file_generation = fgen;
3249 if ((fver <= tpx_incompatible_version) ||
3250 ((fver > tpx_version) && !TopOnlyOK) ||
3251 (fgen > tpx_generation) ||
3252 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3254 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3255 gmx_fio_getname(fio), fver, tpx_version);
3258 do_section(fio, eitemHEADER, bRead);
3259 gmx_fio_do_int(fio, tpx->natoms);
3262 gmx_fio_do_int(fio, tpx->ngtc);
3270 gmx_fio_do_int(fio, idum);
3271 gmx_fio_do_real(fio, rdum);
3273 /*a better decision will eventually (5.0 or later) need to be made
3274 on how to treat the alchemical state of the system, which can now
3275 vary through a simulation, and cannot be completely described
3276 though a single lambda variable, or even a single state
3277 index. Eventually, should probably be a vector. MRS*/
3280 gmx_fio_do_int(fio, tpx->fep_state);
3282 gmx_fio_do_real(fio, tpx->lambda);
3283 gmx_fio_do_int(fio, tpx->bIr);
3284 gmx_fio_do_int(fio, tpx->bTop);
3285 gmx_fio_do_int(fio, tpx->bX);
3286 gmx_fio_do_int(fio, tpx->bV);
3287 gmx_fio_do_int(fio, tpx->bF);
3288 gmx_fio_do_int(fio, tpx->bBox);
3290 if ((fgen > tpx_generation))
3292 /* This can only happen if TopOnlyOK=TRUE */
3297 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3298 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3299 gmx_bool bXVallocated)
3305 int file_version, file_generation;
3309 gmx_bool bPeriodicMols;
3313 tpx.natoms = state->natoms;
3314 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3315 tpx.fep_state = state->fep_state;
3316 tpx.lambda = state->lambda[efptFEP];
3317 tpx.bIr = (ir != NULL);
3318 tpx.bTop = (mtop != NULL);
3319 tpx.bX = (state->x != NULL);
3320 tpx.bV = (state->v != NULL);
3321 tpx.bF = (f != NULL);
3325 TopOnlyOK = (ir == NULL);
3327 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3332 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3333 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3338 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3339 state->natoms = tpx.natoms;
3340 state->nalloc = tpx.natoms;
3346 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3350 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3352 do_test(fio, tpx.bBox, state->box);
3353 do_section(fio, eitemBOX, bRead);
3356 gmx_fio_ndo_rvec(fio, state->box, DIM);
3357 if (file_version >= 51)
3359 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3363 /* We initialize box_rel after reading the inputrec */
3364 clear_mat(state->box_rel);
3366 if (file_version >= 28)
3368 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3369 if (file_version < 56)
3372 gmx_fio_ndo_rvec(fio, mdum, DIM);
3377 if (state->ngtc > 0 && file_version >= 28)
3380 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3381 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3382 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3383 snew(dumv, state->ngtc);
3384 if (file_version < 69)
3386 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3388 /* These used to be the Berendsen tcoupl_lambda's */
3389 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3393 /* Prior to tpx version 26, the inputrec was here.
3394 * I moved it to enable partial forward-compatibility
3395 * for analysis/viewer programs.
3397 if (file_version < 26)
3399 do_test(fio, tpx.bIr, ir);
3400 do_section(fio, eitemIR, bRead);
3405 do_inputrec(fio, ir, bRead, file_version,
3406 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3409 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3414 do_inputrec(fio, &dum_ir, bRead, file_version,
3415 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3418 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3420 done_inputrec(&dum_ir);
3426 do_test(fio, tpx.bTop, mtop);
3427 do_section(fio, eitemTOP, bRead);
3432 do_mtop(fio, mtop, bRead, file_version);
3436 do_mtop(fio, &dum_top, bRead, file_version);
3437 done_mtop(&dum_top, TRUE);
3440 do_test(fio, tpx.bX, state->x);
3441 do_section(fio, eitemX, bRead);
3446 state->flags |= (1<<estX);
3448 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3451 do_test(fio, tpx.bV, state->v);
3452 do_section(fio, eitemV, bRead);
3457 state->flags |= (1<<estV);
3459 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3462 do_test(fio, tpx.bF, f);
3463 do_section(fio, eitemF, bRead);
3466 gmx_fio_ndo_rvec(fio, f, state->natoms);
3469 /* Starting with tpx version 26, we have the inputrec
3470 * at the end of the file, so we can ignore it
3471 * if the file is never than the software (but still the
3472 * same generation - see comments at the top of this file.
3477 bPeriodicMols = FALSE;
3478 if (file_version >= 26)
3480 do_test(fio, tpx.bIr, ir);
3481 do_section(fio, eitemIR, bRead);
3484 if (file_version >= 53)
3486 /* Removed the pbc info from do_inputrec, since we always want it */
3490 bPeriodicMols = ir->bPeriodicMols;
3492 gmx_fio_do_int(fio, ePBC);
3493 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3495 if (file_generation <= tpx_generation && ir)
3497 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3500 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3502 if (file_version < 51)
3504 set_box_rel(ir, state);
3506 if (file_version < 53)
3509 bPeriodicMols = ir->bPeriodicMols;
3512 if (bRead && ir && file_version >= 53)
3514 /* We need to do this after do_inputrec, since that initializes ir */
3516 ir->bPeriodicMols = bPeriodicMols;
3525 if (state->ngtc == 0)
3527 /* Reading old version without tcoupl state data: set it */
3528 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3530 if (tpx.bTop && mtop)
3532 if (file_version < 57)
3534 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3536 ir->eDisre = edrSimple;
3540 ir->eDisre = edrNone;
3543 set_disres_npair(mtop);
3547 if (tpx.bTop && mtop)
3549 gmx_mtop_finalize(mtop);
3552 if (file_version >= 57)
3556 env = getenv("GMX_NOCHARGEGROUPS");
3559 sscanf(env, "%d", &ienv);
3560 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3565 "Will make single atomic charge groups in non-solvent%s\n",
3566 ienv > 1 ? " and solvent" : "");
3567 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3569 fprintf(stderr, "\n");
3577 /************************************************************
3579 * The following routines are the exported ones
3581 ************************************************************/
3583 t_fileio *open_tpx(const char *fn, const char *mode)
3585 return gmx_fio_open(fn, mode);
3588 void close_tpx(t_fileio *fio)
3593 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3594 int *file_version, int *file_generation)
3598 fio = open_tpx(fn, "r");
3599 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3603 void write_tpx_state(const char *fn,
3604 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3608 fio = open_tpx(fn, "w");
3609 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3613 void read_tpx_state(const char *fn,
3614 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3618 fio = open_tpx(fn, "r");
3619 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3623 int read_tpx(const char *fn,
3624 t_inputrec *ir, matrix box, int *natoms,
3625 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3633 fio = open_tpx(fn, "r");
3634 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3636 *natoms = state.natoms;
3639 copy_mat(state.box, box);
3648 int read_tpx_top(const char *fn,
3649 t_inputrec *ir, matrix box, int *natoms,
3650 rvec *x, rvec *v, rvec *f, t_topology *top)
3656 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3658 *top = gmx_mtop_t_to_t_topology(&mtop);
3663 gmx_bool fn2bTPX(const char *file)
3665 switch (fn2ftp(file))
3676 static void done_gmx_groups_t(gmx_groups_t *g)
3680 for (i = 0; (i < egcNR); i++)
3682 if (NULL != g->grps[i].nm_ind)
3684 sfree(g->grps[i].nm_ind);
3685 g->grps[i].nm_ind = NULL;
3687 if (NULL != g->grpnr[i])
3693 /* The contents of this array is in symtab, don't free it here */
3697 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3698 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3701 int natoms, i, version, generation;
3702 gmx_bool bTop, bXNULL = FALSE;
3704 t_topology *topconv;
3707 bTop = fn2bTPX(infile);
3711 read_tpxheader(infile, &header, TRUE, &version, &generation);
3714 snew(*x, header.natoms);
3718 snew(*v, header.natoms);
3721 *ePBC = read_tpx(infile, NULL, box, &natoms,
3722 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3723 *top = gmx_mtop_t_to_t_topology(mtop);
3724 /* In this case we need to throw away the group data too */
3725 done_gmx_groups_t(&mtop->groups);
3727 strcpy(title, *top->name);
3728 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3732 get_stx_coordnum(infile, &natoms);
3733 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3744 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3752 aps = gmx_atomprop_init();
3753 for (i = 0; (i < natoms); i++)
3755 if (!gmx_atomprop_query(aps, epropMass,
3756 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3757 *top->atoms.atomname[i],
3758 &(top->atoms.atom[i].m)))
3762 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3763 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3764 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3765 *top->atoms.atomname[i]);
3769 gmx_atomprop_destroy(aps);
3771 top->idef.ntypes = -1;