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! */
46 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/utility/cstringutil.h"
48 #include "gmx_fatal.h"
61 #include "mtop_util.h"
63 #define TPX_TAG_RELEASE "release"
65 /*! \brief Tag string for the file format written to run input files
66 * written by this version of the code.
68 * Change this if you want to change the run input format in a feature
69 * branch. This ensures that there will not be different run input
70 * formats around which cannot be distinguished, while not causing
71 * problems rebasing the feature branch onto upstream changes. When
72 * merging with mainstream GROMACS, set this tag string back to
73 * TPX_TAG_RELEASE, and instead add an element to tpxv and set
74 * tpx_version to that.
76 static const char *tpx_tag = TPX_TAG_RELEASE;
78 /*! \brief Enum of values that describe the contents of a tpr file
79 * whose format matches a version number
81 * The enum helps the code be more self-documenting and ensure merges
82 * do not silently resolve when two patches make the same bump. When
83 * adding new functionality, add a new element to the end of this
84 * enumeration, change the definition of tpx_version, and write code
85 * below that does the right thing according to the value of
88 tpxv_ComputationalElectrophysiology = 96, /**< support for ion/water position swaps (computational electrophysiology) */
89 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to gmx_int64_t */
90 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
91 tpxv_InteractiveMolecularDynamics, /**< interactive molecular dynamics (IMD) */
92 tpxv_RemoveObsoleteParameters1 /**< remove optimize_fft, dihre_fc, nstcheckpoint */
95 /*! \brief Version number of the file format written to run input
96 * files by this version of the code.
98 * The tpx_version number should be increased whenever the file format
99 * in the main development branch changes, generally to the highest
100 * value present in tpxv. Backward compatibility for reading old run
101 * input files is maintained by checking this version number against
102 * that of the file and then using the correct code path.
104 * When developing a feature branch that needs to change the run input
105 * file format, change tpx_tag instead. */
106 static const int tpx_version = tpxv_RemoveObsoleteParameters1;
109 /* This number should only be increased when you edit the TOPOLOGY section
110 * or the HEADER of the tpx format.
111 * This way we can maintain forward compatibility too for all analysis tools
112 * and/or external programs that only need to know the atom/residue names,
113 * charges, and bond connectivity.
115 * It first appeared in tpx version 26, when I also moved the inputrecord
116 * to the end of the tpx file, so we can just skip it if we only
119 * In particular, it must be increased when adding new elements to
120 * ftupd, so that old code can read new .tpr files.
122 static const int tpx_generation = 26;
124 /* This number should be the most recent backwards incompatible version
125 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
127 static const int tpx_incompatible_version = 9;
131 /* Struct used to maintain tpx compatibility when function types are added */
133 int fvnr; /* file version number in which the function type first appeared */
134 int ftype; /* function type */
138 * The entries should be ordered in:
139 * 1. ascending file version number
140 * 2. ascending function type number
142 /*static const t_ftupd ftupd[] = {
143 { 20, F_CUBICBONDS },
147 { 22, F_DISRESVIOL },
153 { 26, F_DIHRESVIOL },
154 { 30, F_CROSS_BOND_BONDS },
155 { 30, F_CROSS_BOND_ANGLES },
156 { 30, F_UREY_BRADLEY },
157 { 30, F_POLARIZATION },
161 * The entries should be ordered in:
162 * 1. ascending function type number
163 * 2. ascending file version number
165 /* question; what is the purpose of the commented code above? */
166 static const t_ftupd ftupd[] = {
167 { 20, F_CUBICBONDS },
172 { 43, F_TABBONDSNC },
173 { 70, F_RESTRBONDS },
174 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
175 { 76, F_LINEAR_ANGLES },
176 { 30, F_CROSS_BOND_BONDS },
177 { 30, F_CROSS_BOND_ANGLES },
178 { 30, F_UREY_BRADLEY },
179 { 34, F_QUARTIC_ANGLES },
181 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
182 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
191 { 72, F_NPSOLVATION },
193 { 41, F_LJC_PAIRS_NB },
196 { 32, F_COUL_RECIP },
199 { 30, F_POLARIZATION },
202 { 22, F_DISRESVIOL },
206 { 26, F_DIHRESVIOL },
211 { 46, F_ECONSERVED },
212 { 69, F_VTEMP_NOLONGERUSED},
214 { 54, F_DVDL_CONSTR },
215 { 76, F_ANHARM_POL },
218 { 79, F_DVDL_BONDED, },
219 { 79, F_DVDL_RESTRAINT },
220 { 79, F_DVDL_TEMPERATURE },
222 #define NFTUPD asize(ftupd)
224 /* Needed for backward compatibility */
227 static void _do_section(t_fileio *fio, int key, gmx_bool bRead, const char *src,
233 if (gmx_fio_getftp(fio) == efTPA)
237 gmx_fio_write_string(fio, itemstr[key]);
238 bDbg = gmx_fio_getdebug(fio);
239 gmx_fio_setdebug(fio, FALSE);
240 gmx_fio_write_string(fio, comment_str[key]);
241 gmx_fio_setdebug(fio, bDbg);
245 if (gmx_fio_getdebug(fio))
247 fprintf(stderr, "Looking for section %s (%s, %d)",
248 itemstr[key], src, line);
253 gmx_fio_do_string(fio, buf);
255 while ((gmx_strcasecmp(buf, itemstr[key]) != 0));
257 if (gmx_strcasecmp(buf, itemstr[key]) != 0)
259 gmx_fatal(FARGS, "\nCould not find section heading %s", itemstr[key]);
261 else if (gmx_fio_getdebug(fio))
263 fprintf(stderr, " and found it\n");
269 #define do_section(fio, key, bRead) _do_section(fio, key, bRead, __FILE__, __LINE__)
271 /**************************************************************
273 * Now the higer level routines that do io of the structures and arrays
275 **************************************************************/
276 static void do_pullgrp_tpx_pre95(t_fileio *fio,
285 gmx_fio_do_int(fio, pgrp->nat);
288 snew(pgrp->ind, pgrp->nat);
290 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
291 gmx_fio_do_int(fio, pgrp->nweight);
294 snew(pgrp->weight, pgrp->nweight);
296 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
297 gmx_fio_do_int(fio, pgrp->pbcatom);
298 gmx_fio_do_rvec(fio, pcrd->vec);
299 clear_rvec(pcrd->origin);
300 gmx_fio_do_rvec(fio, tmp);
302 gmx_fio_do_real(fio, pcrd->rate);
303 gmx_fio_do_real(fio, pcrd->k);
304 if (file_version >= 56)
306 gmx_fio_do_real(fio, pcrd->kB);
314 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
318 gmx_fio_do_int(fio, pgrp->nat);
321 snew(pgrp->ind, pgrp->nat);
323 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
324 gmx_fio_do_int(fio, pgrp->nweight);
327 snew(pgrp->weight, pgrp->nweight);
329 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
330 gmx_fio_do_int(fio, pgrp->pbcatom);
333 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd)
337 gmx_fio_do_int(fio, pcrd->group[0]);
338 gmx_fio_do_int(fio, pcrd->group[1]);
339 gmx_fio_do_rvec(fio, pcrd->origin);
340 gmx_fio_do_rvec(fio, pcrd->vec);
341 gmx_fio_do_real(fio, pcrd->init);
342 gmx_fio_do_real(fio, pcrd->rate);
343 gmx_fio_do_real(fio, pcrd->k);
344 gmx_fio_do_real(fio, pcrd->kB);
347 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
349 /* i is used in the ndo_double macro*/
353 int n_lambda = fepvals->n_lambda;
355 /* reset the lambda calculation window */
356 fepvals->lambda_start_n = 0;
357 fepvals->lambda_stop_n = n_lambda;
358 if (file_version >= 79)
364 snew(expand->init_lambda_weights, n_lambda);
366 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
367 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
370 gmx_fio_do_int(fio, expand->nstexpanded);
371 gmx_fio_do_int(fio, expand->elmcmove);
372 gmx_fio_do_int(fio, expand->elamstats);
373 gmx_fio_do_int(fio, expand->lmc_repeats);
374 gmx_fio_do_int(fio, expand->gibbsdeltalam);
375 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
376 gmx_fio_do_int(fio, expand->lmc_seed);
377 gmx_fio_do_real(fio, expand->mc_temp);
378 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
379 gmx_fio_do_int(fio, expand->nstTij);
380 gmx_fio_do_int(fio, expand->minvarmin);
381 gmx_fio_do_int(fio, expand->c_range);
382 gmx_fio_do_real(fio, expand->wl_scale);
383 gmx_fio_do_real(fio, expand->wl_ratio);
384 gmx_fio_do_real(fio, expand->init_wl_delta);
385 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
386 gmx_fio_do_int(fio, expand->elmceq);
387 gmx_fio_do_int(fio, expand->equil_steps);
388 gmx_fio_do_int(fio, expand->equil_samples);
389 gmx_fio_do_int(fio, expand->equil_n_at_lam);
390 gmx_fio_do_real(fio, expand->equil_wl_delta);
391 gmx_fio_do_real(fio, expand->equil_ratio);
395 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
398 if (file_version >= 79)
400 gmx_fio_do_int(fio, simtemp->eSimTempScale);
401 gmx_fio_do_real(fio, simtemp->simtemp_high);
402 gmx_fio_do_real(fio, simtemp->simtemp_low);
407 snew(simtemp->temperatures, n_lambda);
409 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
414 static void do_imd(t_fileio *fio, t_IMD *imd, gmx_bool bRead)
416 gmx_fio_do_int(fio, imd->nat);
419 snew(imd->ind, imd->nat);
421 gmx_fio_ndo_int(fio, imd->ind, imd->nat);
424 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
426 /* i is defined in the ndo_double macro; use g to iterate. */
431 /* free energy values */
433 if (file_version >= 79)
435 gmx_fio_do_int(fio, fepvals->init_fep_state);
436 gmx_fio_do_double(fio, fepvals->init_lambda);
437 gmx_fio_do_double(fio, fepvals->delta_lambda);
439 else if (file_version >= 59)
441 gmx_fio_do_double(fio, fepvals->init_lambda);
442 gmx_fio_do_double(fio, fepvals->delta_lambda);
446 gmx_fio_do_real(fio, rdum);
447 fepvals->init_lambda = rdum;
448 gmx_fio_do_real(fio, rdum);
449 fepvals->delta_lambda = rdum;
451 if (file_version >= 79)
453 gmx_fio_do_int(fio, fepvals->n_lambda);
456 snew(fepvals->all_lambda, efptNR);
458 for (g = 0; g < efptNR; g++)
460 if (fepvals->n_lambda > 0)
464 snew(fepvals->all_lambda[g], fepvals->n_lambda);
466 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
467 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
469 else if (fepvals->init_lambda >= 0)
471 fepvals->separate_dvdl[efptFEP] = TRUE;
475 else if (file_version >= 64)
477 gmx_fio_do_int(fio, fepvals->n_lambda);
482 snew(fepvals->all_lambda, efptNR);
483 /* still allocate the all_lambda array's contents. */
484 for (g = 0; g < efptNR; g++)
486 if (fepvals->n_lambda > 0)
488 snew(fepvals->all_lambda[g], fepvals->n_lambda);
492 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
494 if (fepvals->init_lambda >= 0)
498 fepvals->separate_dvdl[efptFEP] = TRUE;
502 /* copy the contents of the efptFEP lambda component to all
503 the other components */
504 for (g = 0; g < efptNR; g++)
506 for (h = 0; h < fepvals->n_lambda; h++)
510 fepvals->all_lambda[g][h] =
511 fepvals->all_lambda[efptFEP][h];
520 fepvals->n_lambda = 0;
521 fepvals->all_lambda = NULL;
522 if (fepvals->init_lambda >= 0)
524 fepvals->separate_dvdl[efptFEP] = TRUE;
527 if (file_version >= 13)
529 gmx_fio_do_real(fio, fepvals->sc_alpha);
533 fepvals->sc_alpha = 0;
535 if (file_version >= 38)
537 gmx_fio_do_int(fio, fepvals->sc_power);
541 fepvals->sc_power = 2;
543 if (file_version >= 79)
545 gmx_fio_do_real(fio, fepvals->sc_r_power);
549 fepvals->sc_r_power = 6.0;
551 if (file_version >= 15)
553 gmx_fio_do_real(fio, fepvals->sc_sigma);
557 fepvals->sc_sigma = 0.3;
561 if (file_version >= 71)
563 fepvals->sc_sigma_min = fepvals->sc_sigma;
567 fepvals->sc_sigma_min = 0;
570 if (file_version >= 79)
572 gmx_fio_do_int(fio, fepvals->bScCoul);
576 fepvals->bScCoul = TRUE;
578 if (file_version >= 64)
580 gmx_fio_do_int(fio, fepvals->nstdhdl);
584 fepvals->nstdhdl = 1;
587 if (file_version >= 73)
589 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
590 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
594 fepvals->separate_dhdl_file = esepdhdlfileYES;
595 fepvals->dhdl_derivatives = edhdlderivativesYES;
597 if (file_version >= 71)
599 gmx_fio_do_int(fio, fepvals->dh_hist_size);
600 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
604 fepvals->dh_hist_size = 0;
605 fepvals->dh_hist_spacing = 0.1;
607 if (file_version >= 79)
609 gmx_fio_do_int(fio, fepvals->bPrintEnergy);
613 fepvals->bPrintEnergy = FALSE;
616 /* handle lambda_neighbors */
617 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
619 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
620 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
621 (fepvals->init_lambda < 0) )
623 fepvals->lambda_start_n = (fepvals->init_fep_state -
624 fepvals->lambda_neighbors);
625 fepvals->lambda_stop_n = (fepvals->init_fep_state +
626 fepvals->lambda_neighbors + 1);
627 if (fepvals->lambda_start_n < 0)
629 fepvals->lambda_start_n = 0;;
631 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
633 fepvals->lambda_stop_n = fepvals->n_lambda;
638 fepvals->lambda_start_n = 0;
639 fepvals->lambda_stop_n = fepvals->n_lambda;
644 fepvals->lambda_start_n = 0;
645 fepvals->lambda_stop_n = fepvals->n_lambda;
649 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
653 if (file_version >= 95)
655 gmx_fio_do_int(fio, pull->ngroup);
657 gmx_fio_do_int(fio, pull->ncoord);
658 if (file_version < 95)
660 pull->ngroup = pull->ncoord + 1;
662 gmx_fio_do_int(fio, pull->eGeom);
663 gmx_fio_do_ivec(fio, pull->dim);
664 gmx_fio_do_real(fio, pull->cyl_r1);
665 gmx_fio_do_real(fio, pull->cyl_r0);
666 gmx_fio_do_real(fio, pull->constr_tol);
667 if (file_version >= 95)
669 gmx_fio_do_int(fio, pull->bPrintRef);
671 gmx_fio_do_int(fio, pull->nstxout);
672 gmx_fio_do_int(fio, pull->nstfout);
675 snew(pull->group, pull->ngroup);
676 snew(pull->coord, pull->ncoord);
678 if (file_version < 95)
680 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
681 if (pull->eGeom == epullgDIRPBC)
683 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
685 if (pull->eGeom > epullgDIRPBC)
690 for (g = 0; g < pull->ngroup; g++)
692 /* We read and ignore a pull coordinate for group 0 */
693 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[max(g-1, 0)],
694 bRead, file_version);
697 pull->coord[g-1].group[0] = 0;
698 pull->coord[g-1].group[1] = g;
702 pull->bPrintRef = (pull->group[0].nat > 0);
706 for (g = 0; g < pull->ngroup; g++)
708 do_pull_group(fio, &pull->group[g], bRead);
710 for (g = 0; g < pull->ncoord; g++)
712 do_pull_coord(fio, &pull->coord[g]);
718 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
722 gmx_fio_do_int(fio, rotg->eType);
723 gmx_fio_do_int(fio, rotg->bMassW);
724 gmx_fio_do_int(fio, rotg->nat);
727 snew(rotg->ind, rotg->nat);
729 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
732 snew(rotg->x_ref, rotg->nat);
734 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
735 gmx_fio_do_rvec(fio, rotg->vec);
736 gmx_fio_do_rvec(fio, rotg->pivot);
737 gmx_fio_do_real(fio, rotg->rate);
738 gmx_fio_do_real(fio, rotg->k);
739 gmx_fio_do_real(fio, rotg->slab_dist);
740 gmx_fio_do_real(fio, rotg->min_gaussian);
741 gmx_fio_do_real(fio, rotg->eps);
742 gmx_fio_do_int(fio, rotg->eFittype);
743 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
744 gmx_fio_do_real(fio, rotg->PotAngle_step);
747 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
751 gmx_fio_do_int(fio, rot->ngrp);
752 gmx_fio_do_int(fio, rot->nstrout);
753 gmx_fio_do_int(fio, rot->nstsout);
756 snew(rot->grp, rot->ngrp);
758 for (g = 0; g < rot->ngrp; g++)
760 do_rotgrp(fio, &rot->grp[g], bRead);
765 static void do_swapcoords(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead)
770 gmx_fio_do_int(fio, swap->nat);
771 gmx_fio_do_int(fio, swap->nat_sol);
772 for (j = 0; j < 2; j++)
774 gmx_fio_do_int(fio, swap->nat_split[j]);
775 gmx_fio_do_int(fio, swap->massw_split[j]);
777 gmx_fio_do_int(fio, swap->nstswap);
778 gmx_fio_do_int(fio, swap->nAverage);
779 gmx_fio_do_real(fio, swap->threshold);
780 gmx_fio_do_real(fio, swap->cyl0r);
781 gmx_fio_do_real(fio, swap->cyl0u);
782 gmx_fio_do_real(fio, swap->cyl0l);
783 gmx_fio_do_real(fio, swap->cyl1r);
784 gmx_fio_do_real(fio, swap->cyl1u);
785 gmx_fio_do_real(fio, swap->cyl1l);
789 snew(swap->ind, swap->nat);
790 snew(swap->ind_sol, swap->nat_sol);
791 for (j = 0; j < 2; j++)
793 snew(swap->ind_split[j], swap->nat_split[j]);
797 gmx_fio_ndo_int(fio, swap->ind, swap->nat);
798 gmx_fio_ndo_int(fio, swap->ind_sol, swap->nat_sol);
799 for (j = 0; j < 2; j++)
801 gmx_fio_ndo_int(fio, swap->ind_split[j], swap->nat_split[j]);
804 for (j = 0; j < eCompNR; j++)
806 gmx_fio_do_int(fio, swap->nanions[j]);
807 gmx_fio_do_int(fio, swap->ncations[j]);
813 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
814 int file_version, real *fudgeQQ)
816 int i, j, k, *tmp, idum = 0;
819 gmx_bool bSimAnn, bdum = 0;
820 real zerotemptime, finish_t, init_temp, finish_temp;
822 if (file_version != tpx_version)
824 /* Give a warning about features that are not accessible */
825 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
826 file_version, tpx_version);
834 if (file_version == 0)
839 /* Basic inputrec stuff */
840 gmx_fio_do_int(fio, ir->eI);
841 if (file_version >= 62)
843 gmx_fio_do_int64(fio, ir->nsteps);
847 gmx_fio_do_int(fio, idum);
851 if (file_version > 25)
853 if (file_version >= 62)
855 gmx_fio_do_int64(fio, ir->init_step);
859 gmx_fio_do_int(fio, idum);
860 ir->init_step = idum;
868 if (file_version >= 58)
870 gmx_fio_do_int(fio, ir->simulation_part);
874 ir->simulation_part = 1;
877 if (file_version >= 67)
879 gmx_fio_do_int(fio, ir->nstcalcenergy);
883 ir->nstcalcenergy = 1;
885 if (file_version < 53)
887 /* The pbc info has been moved out of do_inputrec,
888 * since we always want it, also without reading the inputrec.
890 gmx_fio_do_int(fio, ir->ePBC);
891 if ((file_version <= 15) && (ir->ePBC == 2))
895 if (file_version >= 45)
897 gmx_fio_do_int(fio, ir->bPeriodicMols);
904 ir->bPeriodicMols = TRUE;
908 ir->bPeriodicMols = FALSE;
912 if (file_version >= 81)
914 gmx_fio_do_int(fio, ir->cutoff_scheme);
915 if (file_version < 94)
917 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
922 ir->cutoff_scheme = ecutsGROUP;
924 gmx_fio_do_int(fio, ir->ns_type);
925 gmx_fio_do_int(fio, ir->nstlist);
926 gmx_fio_do_int(fio, ir->ndelta);
927 if (file_version < 41)
929 gmx_fio_do_int(fio, idum);
930 gmx_fio_do_int(fio, idum);
932 if (file_version >= 45)
934 gmx_fio_do_real(fio, ir->rtpi);
940 gmx_fio_do_int(fio, ir->nstcomm);
941 if (file_version > 34)
943 gmx_fio_do_int(fio, ir->comm_mode);
945 else if (ir->nstcomm < 0)
947 ir->comm_mode = ecmANGULAR;
951 ir->comm_mode = ecmLINEAR;
953 ir->nstcomm = abs(ir->nstcomm);
955 /* ignore nstcheckpoint */
956 if (file_version > 25 && file_version < tpxv_RemoveObsoleteParameters1)
958 gmx_fio_do_int(fio, idum);
961 gmx_fio_do_int(fio, ir->nstcgsteep);
963 if (file_version >= 30)
965 gmx_fio_do_int(fio, ir->nbfgscorr);
972 gmx_fio_do_int(fio, ir->nstlog);
973 gmx_fio_do_int(fio, ir->nstxout);
974 gmx_fio_do_int(fio, ir->nstvout);
975 gmx_fio_do_int(fio, ir->nstfout);
976 gmx_fio_do_int(fio, ir->nstenergy);
977 gmx_fio_do_int(fio, ir->nstxout_compressed);
978 if (file_version >= 59)
980 gmx_fio_do_double(fio, ir->init_t);
981 gmx_fio_do_double(fio, ir->delta_t);
985 gmx_fio_do_real(fio, rdum);
987 gmx_fio_do_real(fio, rdum);
990 gmx_fio_do_real(fio, ir->x_compression_precision);
991 if (file_version < 19)
993 gmx_fio_do_int(fio, idum);
994 gmx_fio_do_int(fio, idum);
996 if (file_version < 18)
998 gmx_fio_do_int(fio, idum);
1000 if (file_version >= 81)
1002 gmx_fio_do_real(fio, ir->verletbuf_tol);
1006 ir->verletbuf_tol = 0;
1008 gmx_fio_do_real(fio, ir->rlist);
1009 if (file_version >= 67)
1011 gmx_fio_do_real(fio, ir->rlistlong);
1013 if (file_version >= 82 && file_version != 90)
1015 gmx_fio_do_int(fio, ir->nstcalclr);
1019 /* Calculate at NS steps */
1020 ir->nstcalclr = ir->nstlist;
1022 gmx_fio_do_int(fio, ir->coulombtype);
1023 if (file_version < 32 && ir->coulombtype == eelRF)
1025 ir->coulombtype = eelRF_NEC;
1027 if (file_version >= 81)
1029 gmx_fio_do_int(fio, ir->coulomb_modifier);
1033 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1035 gmx_fio_do_real(fio, ir->rcoulomb_switch);
1036 gmx_fio_do_real(fio, ir->rcoulomb);
1037 gmx_fio_do_int(fio, ir->vdwtype);
1038 if (file_version >= 81)
1040 gmx_fio_do_int(fio, ir->vdw_modifier);
1044 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1046 gmx_fio_do_real(fio, ir->rvdw_switch);
1047 gmx_fio_do_real(fio, ir->rvdw);
1048 if (file_version < 67)
1050 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1052 gmx_fio_do_int(fio, ir->eDispCorr);
1053 gmx_fio_do_real(fio, ir->epsilon_r);
1054 if (file_version >= 37)
1056 gmx_fio_do_real(fio, ir->epsilon_rf);
1060 if (EEL_RF(ir->coulombtype))
1062 ir->epsilon_rf = ir->epsilon_r;
1063 ir->epsilon_r = 1.0;
1067 ir->epsilon_rf = 1.0;
1070 if (file_version >= 29)
1072 gmx_fio_do_real(fio, ir->tabext);
1079 if (file_version > 25)
1081 gmx_fio_do_int(fio, ir->gb_algorithm);
1082 gmx_fio_do_int(fio, ir->nstgbradii);
1083 gmx_fio_do_real(fio, ir->rgbradii);
1084 gmx_fio_do_real(fio, ir->gb_saltconc);
1085 gmx_fio_do_int(fio, ir->implicit_solvent);
1089 ir->gb_algorithm = egbSTILL;
1092 ir->gb_saltconc = 0;
1093 ir->implicit_solvent = eisNO;
1095 if (file_version >= 55)
1097 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1098 gmx_fio_do_real(fio, ir->gb_obc_alpha);
1099 gmx_fio_do_real(fio, ir->gb_obc_beta);
1100 gmx_fio_do_real(fio, ir->gb_obc_gamma);
1101 if (file_version >= 60)
1103 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1104 gmx_fio_do_int(fio, ir->sa_algorithm);
1108 ir->gb_dielectric_offset = 0.009;
1109 ir->sa_algorithm = esaAPPROX;
1111 gmx_fio_do_real(fio, ir->sa_surface_tension);
1113 /* Override sa_surface_tension if it is not changed in the mpd-file */
1114 if (ir->sa_surface_tension < 0)
1116 if (ir->gb_algorithm == egbSTILL)
1118 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1120 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1122 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1129 /* Better use sensible values than insane (0.0) ones... */
1130 ir->gb_epsilon_solvent = 80;
1131 ir->gb_obc_alpha = 1.0;
1132 ir->gb_obc_beta = 0.8;
1133 ir->gb_obc_gamma = 4.85;
1134 ir->sa_surface_tension = 2.092;
1138 if (file_version >= 81)
1140 gmx_fio_do_real(fio, ir->fourier_spacing);
1144 ir->fourier_spacing = 0.0;
1146 gmx_fio_do_int(fio, ir->nkx);
1147 gmx_fio_do_int(fio, ir->nky);
1148 gmx_fio_do_int(fio, ir->nkz);
1149 gmx_fio_do_int(fio, ir->pme_order);
1150 gmx_fio_do_real(fio, ir->ewald_rtol);
1152 if (file_version >= 93)
1154 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1158 ir->ewald_rtol_lj = ir->ewald_rtol;
1161 if (file_version >= 24)
1163 gmx_fio_do_int(fio, ir->ewald_geometry);
1167 ir->ewald_geometry = eewg3D;
1170 if (file_version <= 17)
1172 ir->epsilon_surface = 0;
1173 if (file_version == 17)
1175 gmx_fio_do_int(fio, idum);
1180 gmx_fio_do_real(fio, ir->epsilon_surface);
1183 /* ignore bOptFFT */
1184 if (file_version < tpxv_RemoveObsoleteParameters1)
1186 gmx_fio_do_gmx_bool(fio, bdum);
1189 if (file_version >= 93)
1191 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1193 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1194 gmx_fio_do_int(fio, ir->etc);
1195 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1196 * but the values 0 and 1 still mean no and
1197 * berendsen temperature coupling, respectively.
1199 if (file_version >= 79)
1201 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1203 if (file_version >= 71)
1205 gmx_fio_do_int(fio, ir->nsttcouple);
1209 ir->nsttcouple = ir->nstcalcenergy;
1211 if (file_version <= 15)
1213 gmx_fio_do_int(fio, idum);
1215 if (file_version <= 17)
1217 gmx_fio_do_int(fio, ir->epct);
1218 if (file_version <= 15)
1222 ir->epct = epctSURFACETENSION;
1224 gmx_fio_do_int(fio, idum);
1227 /* we have removed the NO alternative at the beginning */
1231 ir->epct = epctISOTROPIC;
1235 ir->epc = epcBERENDSEN;
1240 gmx_fio_do_int(fio, ir->epc);
1241 gmx_fio_do_int(fio, ir->epct);
1243 if (file_version >= 71)
1245 gmx_fio_do_int(fio, ir->nstpcouple);
1249 ir->nstpcouple = ir->nstcalcenergy;
1251 gmx_fio_do_real(fio, ir->tau_p);
1252 if (file_version <= 15)
1254 gmx_fio_do_rvec(fio, vdum);
1255 clear_mat(ir->ref_p);
1256 for (i = 0; i < DIM; i++)
1258 ir->ref_p[i][i] = vdum[i];
1263 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1264 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1265 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1267 if (file_version <= 15)
1269 gmx_fio_do_rvec(fio, vdum);
1270 clear_mat(ir->compress);
1271 for (i = 0; i < DIM; i++)
1273 ir->compress[i][i] = vdum[i];
1278 gmx_fio_do_rvec(fio, ir->compress[XX]);
1279 gmx_fio_do_rvec(fio, ir->compress[YY]);
1280 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1282 if (file_version >= 47)
1284 gmx_fio_do_int(fio, ir->refcoord_scaling);
1285 gmx_fio_do_rvec(fio, ir->posres_com);
1286 gmx_fio_do_rvec(fio, ir->posres_comB);
1290 ir->refcoord_scaling = erscNO;
1291 clear_rvec(ir->posres_com);
1292 clear_rvec(ir->posres_comB);
1294 if ((file_version > 25) && (file_version < 79))
1296 gmx_fio_do_int(fio, ir->andersen_seed);
1300 ir->andersen_seed = 0;
1302 if (file_version < 26)
1304 gmx_fio_do_gmx_bool(fio, bSimAnn);
1305 gmx_fio_do_real(fio, zerotemptime);
1308 if (file_version < 37)
1310 gmx_fio_do_real(fio, rdum);
1313 gmx_fio_do_real(fio, ir->shake_tol);
1314 if (file_version < 54)
1316 gmx_fio_do_real(fio, *fudgeQQ);
1319 gmx_fio_do_int(fio, ir->efep);
1320 if (file_version <= 14 && ir->efep != efepNO)
1324 do_fepvals(fio, ir->fepvals, bRead, file_version);
1326 if (file_version >= 79)
1328 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1331 ir->bSimTemp = TRUE;
1336 ir->bSimTemp = FALSE;
1340 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1343 if (file_version >= 79)
1345 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1348 ir->bExpanded = TRUE;
1352 ir->bExpanded = FALSE;
1357 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1359 if (file_version >= 57)
1361 gmx_fio_do_int(fio, ir->eDisre);
1363 gmx_fio_do_int(fio, ir->eDisreWeighting);
1364 if (file_version < 22)
1366 if (ir->eDisreWeighting == 0)
1368 ir->eDisreWeighting = edrwEqual;
1372 ir->eDisreWeighting = edrwConservative;
1375 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1376 gmx_fio_do_real(fio, ir->dr_fc);
1377 gmx_fio_do_real(fio, ir->dr_tau);
1378 gmx_fio_do_int(fio, ir->nstdisreout);
1379 if (file_version >= 22)
1381 gmx_fio_do_real(fio, ir->orires_fc);
1382 gmx_fio_do_real(fio, ir->orires_tau);
1383 gmx_fio_do_int(fio, ir->nstorireout);
1389 ir->nstorireout = 0;
1392 /* ignore dihre_fc */
1393 if (file_version >= 26 && file_version < 79)
1395 gmx_fio_do_real(fio, rdum);
1396 if (file_version < 56)
1398 gmx_fio_do_real(fio, rdum);
1399 gmx_fio_do_int(fio, idum);
1403 gmx_fio_do_real(fio, ir->em_stepsize);
1404 gmx_fio_do_real(fio, ir->em_tol);
1405 if (file_version >= 22)
1407 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1411 ir->bShakeSOR = TRUE;
1413 if (file_version >= 11)
1415 gmx_fio_do_int(fio, ir->niter);
1420 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1423 if (file_version >= 21)
1425 gmx_fio_do_real(fio, ir->fc_stepsize);
1429 ir->fc_stepsize = 0;
1431 gmx_fio_do_int(fio, ir->eConstrAlg);
1432 gmx_fio_do_int(fio, ir->nProjOrder);
1433 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1434 if (file_version <= 14)
1436 gmx_fio_do_int(fio, idum);
1438 if (file_version >= 26)
1440 gmx_fio_do_int(fio, ir->nLincsIter);
1445 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1448 if (file_version < 33)
1450 gmx_fio_do_real(fio, bd_temp);
1452 gmx_fio_do_real(fio, ir->bd_fric);
1453 if (file_version >= tpxv_Use64BitRandomSeed)
1455 gmx_fio_do_int64(fio, ir->ld_seed);
1459 gmx_fio_do_int(fio, idum);
1462 if (file_version >= 33)
1464 for (i = 0; i < DIM; i++)
1466 gmx_fio_do_rvec(fio, ir->deform[i]);
1471 for (i = 0; i < DIM; i++)
1473 clear_rvec(ir->deform[i]);
1476 if (file_version >= 14)
1478 gmx_fio_do_real(fio, ir->cos_accel);
1484 gmx_fio_do_int(fio, ir->userint1);
1485 gmx_fio_do_int(fio, ir->userint2);
1486 gmx_fio_do_int(fio, ir->userint3);
1487 gmx_fio_do_int(fio, ir->userint4);
1488 gmx_fio_do_real(fio, ir->userreal1);
1489 gmx_fio_do_real(fio, ir->userreal2);
1490 gmx_fio_do_real(fio, ir->userreal3);
1491 gmx_fio_do_real(fio, ir->userreal4);
1494 if (file_version >= 77)
1496 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1501 snew(ir->adress, 1);
1503 gmx_fio_do_int(fio, ir->adress->type);
1504 gmx_fio_do_real(fio, ir->adress->const_wf);
1505 gmx_fio_do_real(fio, ir->adress->ex_width);
1506 gmx_fio_do_real(fio, ir->adress->hy_width);
1507 gmx_fio_do_int(fio, ir->adress->icor);
1508 gmx_fio_do_int(fio, ir->adress->site);
1509 gmx_fio_do_rvec(fio, ir->adress->refs);
1510 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1511 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1512 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1513 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1517 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1519 if (ir->adress->n_tf_grps > 0)
1521 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1525 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1527 if (ir->adress->n_energy_grps > 0)
1529 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1535 ir->bAdress = FALSE;
1539 if (file_version >= 48)
1541 gmx_fio_do_int(fio, ir->ePull);
1542 if (ir->ePull != epullNO)
1548 do_pull(fio, ir->pull, bRead, file_version);
1553 ir->ePull = epullNO;
1556 /* Enforced rotation */
1557 if (file_version >= 74)
1559 gmx_fio_do_int(fio, ir->bRot);
1560 if (ir->bRot == TRUE)
1566 do_rot(fio, ir->rot, bRead);
1574 /* Interactive molecular dynamics */
1575 if (file_version >= tpxv_InteractiveMolecularDynamics)
1577 gmx_fio_do_int(fio, ir->bIMD);
1578 if (TRUE == ir->bIMD)
1584 do_imd(fio, ir->imd, bRead);
1589 /* We don't support IMD sessions for old .tpr files */
1594 gmx_fio_do_int(fio, ir->opts.ngtc);
1595 if (file_version >= 69)
1597 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1601 ir->opts.nhchainlength = 1;
1603 gmx_fio_do_int(fio, ir->opts.ngacc);
1604 gmx_fio_do_int(fio, ir->opts.ngfrz);
1605 gmx_fio_do_int(fio, ir->opts.ngener);
1609 snew(ir->opts.nrdf, ir->opts.ngtc);
1610 snew(ir->opts.ref_t, ir->opts.ngtc);
1611 snew(ir->opts.annealing, ir->opts.ngtc);
1612 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1613 snew(ir->opts.anneal_time, ir->opts.ngtc);
1614 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1615 snew(ir->opts.tau_t, ir->opts.ngtc);
1616 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1617 snew(ir->opts.acc, ir->opts.ngacc);
1618 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1620 if (ir->opts.ngtc > 0)
1622 if (bRead && file_version < 13)
1624 snew(tmp, ir->opts.ngtc);
1625 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1626 for (i = 0; i < ir->opts.ngtc; i++)
1628 ir->opts.nrdf[i] = tmp[i];
1634 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1636 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1637 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1638 if (file_version < 33 && ir->eI == eiBD)
1640 for (i = 0; i < ir->opts.ngtc; i++)
1642 ir->opts.tau_t[i] = bd_temp;
1646 if (ir->opts.ngfrz > 0)
1648 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1650 if (ir->opts.ngacc > 0)
1652 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1654 if (file_version >= 12)
1656 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1657 ir->opts.ngener*ir->opts.ngener);
1660 if (bRead && file_version < 26)
1662 for (i = 0; i < ir->opts.ngtc; i++)
1666 ir->opts.annealing[i] = eannSINGLE;
1667 ir->opts.anneal_npoints[i] = 2;
1668 snew(ir->opts.anneal_time[i], 2);
1669 snew(ir->opts.anneal_temp[i], 2);
1670 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1671 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1672 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1673 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1674 ir->opts.anneal_time[i][0] = ir->init_t;
1675 ir->opts.anneal_time[i][1] = finish_t;
1676 ir->opts.anneal_temp[i][0] = init_temp;
1677 ir->opts.anneal_temp[i][1] = finish_temp;
1681 ir->opts.annealing[i] = eannNO;
1682 ir->opts.anneal_npoints[i] = 0;
1688 /* file version 26 or later */
1689 /* First read the lists with annealing and npoints for each group */
1690 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1691 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1692 for (j = 0; j < (ir->opts.ngtc); j++)
1694 k = ir->opts.anneal_npoints[j];
1697 snew(ir->opts.anneal_time[j], k);
1698 snew(ir->opts.anneal_temp[j], k);
1700 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1701 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1705 if (file_version >= 45)
1707 gmx_fio_do_int(fio, ir->nwall);
1708 gmx_fio_do_int(fio, ir->wall_type);
1709 if (file_version >= 50)
1711 gmx_fio_do_real(fio, ir->wall_r_linpot);
1715 ir->wall_r_linpot = -1;
1717 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1718 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1719 gmx_fio_do_real(fio, ir->wall_density[0]);
1720 gmx_fio_do_real(fio, ir->wall_density[1]);
1721 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1727 ir->wall_atomtype[0] = -1;
1728 ir->wall_atomtype[1] = -1;
1729 ir->wall_density[0] = 0;
1730 ir->wall_density[1] = 0;
1731 ir->wall_ewald_zfac = 3;
1733 /* Cosine stuff for electric fields */
1734 for (j = 0; (j < DIM); j++)
1736 gmx_fio_do_int(fio, ir->ex[j].n);
1737 gmx_fio_do_int(fio, ir->et[j].n);
1740 snew(ir->ex[j].a, ir->ex[j].n);
1741 snew(ir->ex[j].phi, ir->ex[j].n);
1742 snew(ir->et[j].a, ir->et[j].n);
1743 snew(ir->et[j].phi, ir->et[j].n);
1745 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1746 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1747 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1748 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1752 if (file_version >= tpxv_ComputationalElectrophysiology)
1754 gmx_fio_do_int(fio, ir->eSwapCoords);
1755 if (ir->eSwapCoords != eswapNO)
1761 do_swapcoords(fio, ir->swap, bRead);
1766 if (file_version >= 39)
1768 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1769 gmx_fio_do_int(fio, ir->QMMMscheme);
1770 gmx_fio_do_real(fio, ir->scalefactor);
1771 gmx_fio_do_int(fio, ir->opts.ngQM);
1774 snew(ir->opts.QMmethod, ir->opts.ngQM);
1775 snew(ir->opts.QMbasis, ir->opts.ngQM);
1776 snew(ir->opts.QMcharge, ir->opts.ngQM);
1777 snew(ir->opts.QMmult, ir->opts.ngQM);
1778 snew(ir->opts.bSH, ir->opts.ngQM);
1779 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1780 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1781 snew(ir->opts.SAon, ir->opts.ngQM);
1782 snew(ir->opts.SAoff, ir->opts.ngQM);
1783 snew(ir->opts.SAsteps, ir->opts.ngQM);
1784 snew(ir->opts.bOPT, ir->opts.ngQM);
1785 snew(ir->opts.bTS, ir->opts.ngQM);
1787 if (ir->opts.ngQM > 0)
1789 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1790 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1791 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1792 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1793 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1794 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1795 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1796 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1797 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1798 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1799 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1800 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1802 /* end of QMMM stuff */
1807 static void do_harm(t_fileio *fio, t_iparams *iparams)
1809 gmx_fio_do_real(fio, iparams->harmonic.rA);
1810 gmx_fio_do_real(fio, iparams->harmonic.krA);
1811 gmx_fio_do_real(fio, iparams->harmonic.rB);
1812 gmx_fio_do_real(fio, iparams->harmonic.krB);
1815 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1816 gmx_bool bRead, int file_version)
1823 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1833 do_harm(fio, iparams);
1834 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1836 /* Correct incorrect storage of parameters */
1837 iparams->pdihs.phiB = iparams->pdihs.phiA;
1838 iparams->pdihs.cpB = iparams->pdihs.cpA;
1842 gmx_fio_do_real(fio, iparams->harmonic.rA);
1843 gmx_fio_do_real(fio, iparams->harmonic.krA);
1845 case F_LINEAR_ANGLES:
1846 gmx_fio_do_real(fio, iparams->linangle.klinA);
1847 gmx_fio_do_real(fio, iparams->linangle.aA);
1848 gmx_fio_do_real(fio, iparams->linangle.klinB);
1849 gmx_fio_do_real(fio, iparams->linangle.aB);
1852 gmx_fio_do_real(fio, iparams->fene.bm);
1853 gmx_fio_do_real(fio, iparams->fene.kb);
1857 gmx_fio_do_real(fio, iparams->restraint.lowA);
1858 gmx_fio_do_real(fio, iparams->restraint.up1A);
1859 gmx_fio_do_real(fio, iparams->restraint.up2A);
1860 gmx_fio_do_real(fio, iparams->restraint.kA);
1861 gmx_fio_do_real(fio, iparams->restraint.lowB);
1862 gmx_fio_do_real(fio, iparams->restraint.up1B);
1863 gmx_fio_do_real(fio, iparams->restraint.up2B);
1864 gmx_fio_do_real(fio, iparams->restraint.kB);
1870 gmx_fio_do_real(fio, iparams->tab.kA);
1871 gmx_fio_do_int(fio, iparams->tab.table);
1872 gmx_fio_do_real(fio, iparams->tab.kB);
1874 case F_CROSS_BOND_BONDS:
1875 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1876 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1877 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1879 case F_CROSS_BOND_ANGLES:
1880 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1881 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1882 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1883 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1885 case F_UREY_BRADLEY:
1886 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1887 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1888 gmx_fio_do_real(fio, iparams->u_b.r13A);
1889 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1890 if (file_version >= 79)
1892 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1893 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1894 gmx_fio_do_real(fio, iparams->u_b.r13B);
1895 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1899 iparams->u_b.thetaB = iparams->u_b.thetaA;
1900 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1901 iparams->u_b.r13B = iparams->u_b.r13A;
1902 iparams->u_b.kUBB = iparams->u_b.kUBA;
1905 case F_QUARTIC_ANGLES:
1906 gmx_fio_do_real(fio, iparams->qangle.theta);
1907 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1910 gmx_fio_do_real(fio, iparams->bham.a);
1911 gmx_fio_do_real(fio, iparams->bham.b);
1912 gmx_fio_do_real(fio, iparams->bham.c);
1915 gmx_fio_do_real(fio, iparams->morse.b0A);
1916 gmx_fio_do_real(fio, iparams->morse.cbA);
1917 gmx_fio_do_real(fio, iparams->morse.betaA);
1918 if (file_version >= 79)
1920 gmx_fio_do_real(fio, iparams->morse.b0B);
1921 gmx_fio_do_real(fio, iparams->morse.cbB);
1922 gmx_fio_do_real(fio, iparams->morse.betaB);
1926 iparams->morse.b0B = iparams->morse.b0A;
1927 iparams->morse.cbB = iparams->morse.cbA;
1928 iparams->morse.betaB = iparams->morse.betaA;
1932 gmx_fio_do_real(fio, iparams->cubic.b0);
1933 gmx_fio_do_real(fio, iparams->cubic.kb);
1934 gmx_fio_do_real(fio, iparams->cubic.kcub);
1938 case F_POLARIZATION:
1939 gmx_fio_do_real(fio, iparams->polarize.alpha);
1942 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1943 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1944 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1947 if (file_version < 31)
1949 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1951 gmx_fio_do_real(fio, iparams->wpol.al_x);
1952 gmx_fio_do_real(fio, iparams->wpol.al_y);
1953 gmx_fio_do_real(fio, iparams->wpol.al_z);
1954 gmx_fio_do_real(fio, iparams->wpol.rOH);
1955 gmx_fio_do_real(fio, iparams->wpol.rHH);
1956 gmx_fio_do_real(fio, iparams->wpol.rOD);
1959 gmx_fio_do_real(fio, iparams->thole.a);
1960 gmx_fio_do_real(fio, iparams->thole.alpha1);
1961 gmx_fio_do_real(fio, iparams->thole.alpha2);
1962 gmx_fio_do_real(fio, iparams->thole.rfac);
1965 gmx_fio_do_real(fio, iparams->lj.c6);
1966 gmx_fio_do_real(fio, iparams->lj.c12);
1969 gmx_fio_do_real(fio, iparams->lj14.c6A);
1970 gmx_fio_do_real(fio, iparams->lj14.c12A);
1971 gmx_fio_do_real(fio, iparams->lj14.c6B);
1972 gmx_fio_do_real(fio, iparams->lj14.c12B);
1975 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1976 gmx_fio_do_real(fio, iparams->ljc14.qi);
1977 gmx_fio_do_real(fio, iparams->ljc14.qj);
1978 gmx_fio_do_real(fio, iparams->ljc14.c6);
1979 gmx_fio_do_real(fio, iparams->ljc14.c12);
1981 case F_LJC_PAIRS_NB:
1982 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1983 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1984 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1985 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1991 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1992 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1993 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1995 /* Read the incorrectly stored multiplicity */
1996 gmx_fio_do_real(fio, iparams->harmonic.rB);
1997 gmx_fio_do_real(fio, iparams->harmonic.krB);
1998 iparams->pdihs.phiB = iparams->pdihs.phiA;
1999 iparams->pdihs.cpB = iparams->pdihs.cpA;
2003 gmx_fio_do_real(fio, iparams->pdihs.phiB);
2004 gmx_fio_do_real(fio, iparams->pdihs.cpB);
2005 gmx_fio_do_int(fio, iparams->pdihs.mult);
2009 gmx_fio_do_real(fio, iparams->pdihs.phiA);
2010 gmx_fio_do_real(fio, iparams->pdihs.cpA);
2013 gmx_fio_do_int(fio, iparams->disres.label);
2014 gmx_fio_do_int(fio, iparams->disres.type);
2015 gmx_fio_do_real(fio, iparams->disres.low);
2016 gmx_fio_do_real(fio, iparams->disres.up1);
2017 gmx_fio_do_real(fio, iparams->disres.up2);
2018 gmx_fio_do_real(fio, iparams->disres.kfac);
2021 gmx_fio_do_int(fio, iparams->orires.ex);
2022 gmx_fio_do_int(fio, iparams->orires.label);
2023 gmx_fio_do_int(fio, iparams->orires.power);
2024 gmx_fio_do_real(fio, iparams->orires.c);
2025 gmx_fio_do_real(fio, iparams->orires.obs);
2026 gmx_fio_do_real(fio, iparams->orires.kfac);
2029 if (file_version < 82)
2031 gmx_fio_do_int(fio, idum);
2032 gmx_fio_do_int(fio, idum);
2034 gmx_fio_do_real(fio, iparams->dihres.phiA);
2035 gmx_fio_do_real(fio, iparams->dihres.dphiA);
2036 gmx_fio_do_real(fio, iparams->dihres.kfacA);
2037 if (file_version >= 82)
2039 gmx_fio_do_real(fio, iparams->dihres.phiB);
2040 gmx_fio_do_real(fio, iparams->dihres.dphiB);
2041 gmx_fio_do_real(fio, iparams->dihres.kfacB);
2045 iparams->dihres.phiB = iparams->dihres.phiA;
2046 iparams->dihres.dphiB = iparams->dihres.dphiA;
2047 iparams->dihres.kfacB = iparams->dihres.kfacA;
2051 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
2052 gmx_fio_do_rvec(fio, iparams->posres.fcA);
2053 if (bRead && file_version < 27)
2055 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
2056 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
2060 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
2061 gmx_fio_do_rvec(fio, iparams->posres.fcB);
2065 gmx_fio_do_int(fio, iparams->fbposres.geom);
2066 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
2067 gmx_fio_do_real(fio, iparams->fbposres.r);
2068 gmx_fio_do_real(fio, iparams->fbposres.k);
2071 gmx_fio_ndo_real(fio, iparams->cbtdihs.cbtcA, NR_CBTDIHS);
2074 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2075 if (file_version >= 25)
2077 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2081 /* Fourier dihedrals are internally represented
2082 * as Ryckaert-Bellemans since those are faster to compute.
2084 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2085 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2089 gmx_fio_do_real(fio, iparams->constr.dA);
2090 gmx_fio_do_real(fio, iparams->constr.dB);
2093 gmx_fio_do_real(fio, iparams->settle.doh);
2094 gmx_fio_do_real(fio, iparams->settle.dhh);
2097 gmx_fio_do_real(fio, iparams->vsite.a);
2102 gmx_fio_do_real(fio, iparams->vsite.a);
2103 gmx_fio_do_real(fio, iparams->vsite.b);
2108 gmx_fio_do_real(fio, iparams->vsite.a);
2109 gmx_fio_do_real(fio, iparams->vsite.b);
2110 gmx_fio_do_real(fio, iparams->vsite.c);
2113 gmx_fio_do_int(fio, iparams->vsiten.n);
2114 gmx_fio_do_real(fio, iparams->vsiten.a);
2119 /* We got rid of some parameters in version 68 */
2120 if (bRead && file_version < 68)
2122 gmx_fio_do_real(fio, rdum);
2123 gmx_fio_do_real(fio, rdum);
2124 gmx_fio_do_real(fio, rdum);
2125 gmx_fio_do_real(fio, rdum);
2127 gmx_fio_do_real(fio, iparams->gb.sar);
2128 gmx_fio_do_real(fio, iparams->gb.st);
2129 gmx_fio_do_real(fio, iparams->gb.pi);
2130 gmx_fio_do_real(fio, iparams->gb.gbr);
2131 gmx_fio_do_real(fio, iparams->gb.bmlt);
2134 gmx_fio_do_int(fio, iparams->cmap.cmapA);
2135 gmx_fio_do_int(fio, iparams->cmap.cmapB);
2138 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
2139 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
2143 gmx_fio_unset_comment(fio);
2147 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
2154 gmx_fio_set_comment(fio, interaction_function[ftype].name);
2156 if (file_version < 44)
2158 for (i = 0; i < MAXNODES; i++)
2160 gmx_fio_do_int(fio, idum);
2163 gmx_fio_do_int(fio, ilist->nr);
2166 snew(ilist->iatoms, ilist->nr);
2168 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2171 gmx_fio_unset_comment(fio);
2175 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2176 gmx_bool bRead, int file_version)
2181 gmx_fio_do_int(fio, ffparams->atnr);
2182 if (file_version < 57)
2184 gmx_fio_do_int(fio, idum);
2186 gmx_fio_do_int(fio, ffparams->ntypes);
2189 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2190 ffparams->atnr, ffparams->ntypes);
2194 snew(ffparams->functype, ffparams->ntypes);
2195 snew(ffparams->iparams, ffparams->ntypes);
2197 /* Read/write all the function types */
2198 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2201 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2204 if (file_version >= 66)
2206 gmx_fio_do_double(fio, ffparams->reppow);
2210 ffparams->reppow = 12.0;
2213 if (file_version >= 57)
2215 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2218 /* Check whether all these function types are supported by the code.
2219 * In practice the code is backwards compatible, which means that the
2220 * numbering may have to be altered from old numbering to new numbering
2222 for (i = 0; (i < ffparams->ntypes); i++)
2226 /* Loop over file versions */
2227 for (k = 0; (k < NFTUPD); k++)
2229 /* Compare the read file_version to the update table */
2230 if ((file_version < ftupd[k].fvnr) &&
2231 (ffparams->functype[i] >= ftupd[k].ftype))
2233 ffparams->functype[i] += 1;
2236 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2237 i, ffparams->functype[i],
2238 interaction_function[ftupd[k].ftype].longname);
2245 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2249 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2254 static void add_settle_atoms(t_ilist *ilist)
2258 /* Settle used to only store the first atom: add the other two */
2259 srenew(ilist->iatoms, 2*ilist->nr);
2260 for (i = ilist->nr/2-1; i >= 0; i--)
2262 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2263 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2264 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2265 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2267 ilist->nr = 2*ilist->nr;
2270 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2273 int i, j, renum[F_NRE];
2277 for (j = 0; (j < F_NRE); j++)
2282 for (k = 0; k < NFTUPD; k++)
2284 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2293 ilist[j].iatoms = NULL;
2297 do_ilist(fio, &ilist[j], bRead, file_version, j);
2298 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2300 add_settle_atoms(&ilist[j]);
2304 if (bRead && gmx_debug_at)
2305 pr_ilist(debug,0,interaction_function[j].longname,
2306 functype,&ilist[j],TRUE);
2311 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2312 gmx_bool bRead, int file_version)
2314 do_ffparams(fio, ffparams, bRead, file_version);
2316 if (file_version >= 54)
2318 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2321 do_ilists(fio, molt->ilist, bRead, file_version);
2324 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2326 int i, idum, dum_nra, *dum_a;
2328 if (file_version < 44)
2330 for (i = 0; i < MAXNODES; i++)
2332 gmx_fio_do_int(fio, idum);
2335 gmx_fio_do_int(fio, block->nr);
2336 if (file_version < 51)
2338 gmx_fio_do_int(fio, dum_nra);
2342 if ((block->nalloc_index > 0) && (NULL != block->index))
2344 sfree(block->index);
2346 block->nalloc_index = block->nr+1;
2347 snew(block->index, block->nalloc_index);
2349 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2351 if (file_version < 51 && dum_nra > 0)
2353 snew(dum_a, dum_nra);
2354 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2359 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2364 if (file_version < 44)
2366 for (i = 0; i < MAXNODES; i++)
2368 gmx_fio_do_int(fio, idum);
2371 gmx_fio_do_int(fio, block->nr);
2372 gmx_fio_do_int(fio, block->nra);
2375 block->nalloc_index = block->nr+1;
2376 snew(block->index, block->nalloc_index);
2377 block->nalloc_a = block->nra;
2378 snew(block->a, block->nalloc_a);
2380 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2381 gmx_fio_ndo_int(fio, block->a, block->nra);
2384 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2385 int file_version, gmx_groups_t *groups, int atnr)
2389 gmx_fio_do_real(fio, atom->m);
2390 gmx_fio_do_real(fio, atom->q);
2391 gmx_fio_do_real(fio, atom->mB);
2392 gmx_fio_do_real(fio, atom->qB);
2393 gmx_fio_do_ushort(fio, atom->type);
2394 gmx_fio_do_ushort(fio, atom->typeB);
2395 gmx_fio_do_int(fio, atom->ptype);
2396 gmx_fio_do_int(fio, atom->resind);
2397 if (file_version >= 52)
2399 gmx_fio_do_int(fio, atom->atomnumber);
2403 atom->atomnumber = NOTSET;
2405 if (file_version < 23)
2409 else if (file_version < 39)
2418 if (file_version < 57)
2420 unsigned char uchar[egcNR];
2421 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2422 for (i = myngrp; (i < ngrp); i++)
2426 /* Copy the old data format to the groups struct */
2427 for (i = 0; i < ngrp; i++)
2429 groups->grpnr[i][atnr] = uchar[i];
2434 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2439 if (file_version < 23)
2443 else if (file_version < 39)
2452 for (j = 0; (j < ngrp); j++)
2456 gmx_fio_do_int(fio, grps[j].nr);
2459 snew(grps[j].nm_ind, grps[j].nr);
2461 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2466 snew(grps[j].nm_ind, grps[j].nr);
2471 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2477 gmx_fio_do_int(fio, ls);
2478 *nm = get_symtab_handle(symtab, ls);
2482 ls = lookup_symtab(symtab, *nm);
2483 gmx_fio_do_int(fio, ls);
2487 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2492 for (j = 0; (j < nstr); j++)
2494 do_symstr(fio, &(nm[j]), bRead, symtab);
2498 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2499 t_symtab *symtab, int file_version)
2503 for (j = 0; (j < n); j++)
2505 do_symstr(fio, &(ri[j].name), bRead, symtab);
2506 if (file_version >= 63)
2508 gmx_fio_do_int(fio, ri[j].nr);
2509 gmx_fio_do_uchar(fio, ri[j].ic);
2519 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2521 gmx_groups_t *groups)
2525 gmx_fio_do_int(fio, atoms->nr);
2526 gmx_fio_do_int(fio, atoms->nres);
2527 if (file_version < 57)
2529 gmx_fio_do_int(fio, groups->ngrpname);
2530 for (i = 0; i < egcNR; i++)
2532 groups->ngrpnr[i] = atoms->nr;
2533 snew(groups->grpnr[i], groups->ngrpnr[i]);
2538 snew(atoms->atom, atoms->nr);
2539 snew(atoms->atomname, atoms->nr);
2540 snew(atoms->atomtype, atoms->nr);
2541 snew(atoms->atomtypeB, atoms->nr);
2542 snew(atoms->resinfo, atoms->nres);
2543 if (file_version < 57)
2545 snew(groups->grpname, groups->ngrpname);
2547 atoms->pdbinfo = NULL;
2549 for (i = 0; (i < atoms->nr); i++)
2551 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2553 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2554 if (bRead && (file_version <= 20))
2556 for (i = 0; i < atoms->nr; i++)
2558 atoms->atomtype[i] = put_symtab(symtab, "?");
2559 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2564 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2565 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2567 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2569 if (file_version < 57)
2571 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2573 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2577 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2578 gmx_bool bRead, t_symtab *symtab,
2583 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2584 gmx_fio_do_int(fio, groups->ngrpname);
2587 snew(groups->grpname, groups->ngrpname);
2589 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2590 for (g = 0; g < egcNR; g++)
2592 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2593 if (groups->ngrpnr[g] == 0)
2597 groups->grpnr[g] = NULL;
2604 snew(groups->grpnr[g], groups->ngrpnr[g]);
2606 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2611 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2616 if (file_version > 25)
2618 gmx_fio_do_int(fio, atomtypes->nr);
2622 snew(atomtypes->radius, j);
2623 snew(atomtypes->vol, j);
2624 snew(atomtypes->surftens, j);
2625 snew(atomtypes->atomnumber, j);
2626 snew(atomtypes->gb_radius, j);
2627 snew(atomtypes->S_hct, j);
2629 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2630 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2631 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2632 if (file_version >= 40)
2634 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2636 if (file_version >= 60)
2638 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2639 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2644 /* File versions prior to 26 cannot do GBSA,
2645 * so they dont use this structure
2648 atomtypes->radius = NULL;
2649 atomtypes->vol = NULL;
2650 atomtypes->surftens = NULL;
2651 atomtypes->atomnumber = NULL;
2652 atomtypes->gb_radius = NULL;
2653 atomtypes->S_hct = NULL;
2657 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2663 gmx_fio_do_int(fio, symtab->nr);
2667 snew(symtab->symbuf, 1);
2668 symbuf = symtab->symbuf;
2669 symbuf->bufsize = nr;
2670 snew(symbuf->buf, nr);
2671 for (i = 0; (i < nr); i++)
2673 gmx_fio_do_string(fio, buf);
2674 symbuf->buf[i] = strdup(buf);
2679 symbuf = symtab->symbuf;
2680 while (symbuf != NULL)
2682 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2684 gmx_fio_do_string(fio, symbuf->buf[i]);
2687 symbuf = symbuf->next;
2691 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2696 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2698 int i, j, ngrid, gs, nelem;
2700 gmx_fio_do_int(fio, cmap_grid->ngrid);
2701 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2703 ngrid = cmap_grid->ngrid;
2704 gs = cmap_grid->grid_spacing;
2709 snew(cmap_grid->cmapdata, ngrid);
2711 for (i = 0; i < cmap_grid->ngrid; i++)
2713 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2717 for (i = 0; i < cmap_grid->ngrid; i++)
2719 for (j = 0; j < nelem; j++)
2721 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2722 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2723 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2724 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2730 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2732 int m, a, a0, a1, r;
2736 /* We always assign a new chain number, but save the chain id characters
2737 * for larger molecules.
2739 #define CHAIN_MIN_ATOMS 15
2743 for (m = 0; m < mols->nr; m++)
2745 a0 = mols->index[m];
2746 a1 = mols->index[m+1];
2747 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2756 for (a = a0; a < a1; a++)
2758 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2759 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2764 /* Blank out the chain id if there was only one chain */
2767 for (r = 0; r < atoms->nres; r++)
2769 atoms->resinfo[r].chainid = ' ';
2774 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2775 t_symtab *symtab, int file_version,
2776 gmx_groups_t *groups)
2780 if (file_version >= 57)
2782 do_symstr(fio, &(molt->name), bRead, symtab);
2785 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2787 if (bRead && gmx_debug_at)
2789 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2792 if (file_version >= 57)
2794 do_ilists(fio, molt->ilist, bRead, file_version);
2796 do_block(fio, &molt->cgs, bRead, file_version);
2797 if (bRead && gmx_debug_at)
2799 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2803 /* This used to be in the atoms struct */
2804 do_blocka(fio, &molt->excls, bRead, file_version);
2807 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2811 gmx_fio_do_int(fio, molb->type);
2812 gmx_fio_do_int(fio, molb->nmol);
2813 gmx_fio_do_int(fio, molb->natoms_mol);
2814 /* Position restraint coordinates */
2815 gmx_fio_do_int(fio, molb->nposres_xA);
2816 if (molb->nposres_xA > 0)
2820 snew(molb->posres_xA, molb->nposres_xA);
2822 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2824 gmx_fio_do_int(fio, molb->nposres_xB);
2825 if (molb->nposres_xB > 0)
2829 snew(molb->posres_xB, molb->nposres_xB);
2831 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2836 static t_block mtop_mols(gmx_mtop_t *mtop)
2842 for (mb = 0; mb < mtop->nmolblock; mb++)
2844 mols.nr += mtop->molblock[mb].nmol;
2846 mols.nalloc_index = mols.nr + 1;
2847 snew(mols.index, mols.nalloc_index);
2852 for (mb = 0; mb < mtop->nmolblock; mb++)
2854 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2856 a += mtop->molblock[mb].natoms_mol;
2865 static void add_posres_molblock(gmx_mtop_t *mtop)
2870 gmx_molblock_t *molb;
2873 /* posres reference positions are stored in ip->posres (if present) and
2874 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2875 posres.pos0A are identical to fbposres.pos0. */
2876 il = &mtop->moltype[0].ilist[F_POSRES];
2877 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2878 if (il->nr == 0 && ilfb->nr == 0)
2884 for (i = 0; i < il->nr; i += 2)
2886 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2887 am = max(am, il->iatoms[i+1]);
2888 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2889 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2890 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2895 /* This loop is required if we have only flat-bottomed posres:
2897 - bFE == FALSE (no B-state for flat-bottomed posres) */
2900 for (i = 0; i < ilfb->nr; i += 2)
2902 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2903 am = max(am, ilfb->iatoms[i+1]);
2906 /* Make the posres coordinate block end at a molecule end */
2908 while (am >= mtop->mols.index[mol+1])
2912 molb = &mtop->molblock[0];
2913 molb->nposres_xA = mtop->mols.index[mol+1];
2914 snew(molb->posres_xA, molb->nposres_xA);
2917 molb->nposres_xB = molb->nposres_xA;
2918 snew(molb->posres_xB, molb->nposres_xB);
2922 molb->nposres_xB = 0;
2924 for (i = 0; i < il->nr; i += 2)
2926 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2927 a = il->iatoms[i+1];
2928 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2929 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2930 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2933 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2934 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2935 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2940 /* If only flat-bottomed posres are present, take reference pos from them.
2941 Here: bFE == FALSE */
2942 for (i = 0; i < ilfb->nr; i += 2)
2944 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2945 a = ilfb->iatoms[i+1];
2946 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2947 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2948 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2953 static void set_disres_npair(gmx_mtop_t *mtop)
2960 ip = mtop->ffparams.iparams;
2962 for (mt = 0; mt < mtop->nmoltype; mt++)
2964 il = &mtop->moltype[mt].ilist[F_DISRES];
2969 for (i = 0; i < il->nr; i += 3)
2972 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2974 ip[a[i]].disres.npair = npair;
2982 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2992 do_symtab(fio, &(mtop->symtab), bRead);
2995 pr_symtab(debug, 0, "symtab", &mtop->symtab);
2998 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
3000 if (file_version >= 57)
3002 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
3004 gmx_fio_do_int(fio, mtop->nmoltype);
3012 snew(mtop->moltype, mtop->nmoltype);
3013 if (file_version < 57)
3015 mtop->moltype[0].name = mtop->name;
3018 for (mt = 0; mt < mtop->nmoltype; mt++)
3020 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
3024 if (file_version >= 57)
3026 gmx_fio_do_int(fio, mtop->nmolblock);
3030 mtop->nmolblock = 1;
3034 snew(mtop->molblock, mtop->nmolblock);
3036 if (file_version >= 57)
3038 for (mb = 0; mb < mtop->nmolblock; mb++)
3040 do_molblock(fio, &mtop->molblock[mb], bRead);
3042 gmx_fio_do_int(fio, mtop->natoms);
3046 mtop->molblock[0].type = 0;
3047 mtop->molblock[0].nmol = 1;
3048 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
3049 mtop->molblock[0].nposres_xA = 0;
3050 mtop->molblock[0].nposres_xB = 0;
3053 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
3056 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
3059 if (file_version < 57)
3061 /* Debug statements are inside do_idef */
3062 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
3063 mtop->natoms = mtop->moltype[0].atoms.nr;
3066 if (file_version >= 65)
3068 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
3072 mtop->ffparams.cmap_grid.ngrid = 0;
3073 mtop->ffparams.cmap_grid.grid_spacing = 0;
3074 mtop->ffparams.cmap_grid.cmapdata = NULL;
3077 if (file_version >= 57)
3079 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
3082 if (file_version < 57)
3084 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
3085 if (bRead && gmx_debug_at)
3087 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
3089 do_block(fio, &mtop->mols, bRead, file_version);
3090 /* Add the posres coordinates to the molblock */
3091 add_posres_molblock(mtop);
3095 if (file_version >= 57)
3097 done_block(&mtop->mols);
3098 mtop->mols = mtop_mols(mtop);
3102 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
3106 if (file_version < 51)
3108 /* Here used to be the shake blocks */
3109 do_blocka(fio, &dumb, bRead, file_version);
3122 close_symtab(&(mtop->symtab));
3126 /* If TopOnlyOK is TRUE then we can read even future versions
3127 * of tpx files, provided the file_generation hasn't changed.
3128 * If it is FALSE, we need the inputrecord too, and bail out
3129 * if the file is newer than the program.
3131 * The version and generation if the topology (see top of this file)
3132 * are returned in the two last arguments.
3134 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
3136 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
3137 gmx_bool TopOnlyOK, int *file_version,
3138 int *file_generation)
3141 char file_tag[STRLEN];
3148 gmx_fio_checktype(fio);
3149 gmx_fio_setdebug(fio, bDebugMode());
3151 /* NEW! XDR tpb file */
3152 precision = sizeof(real);
3155 gmx_fio_do_string(fio, buf);
3156 if (strncmp(buf, "VERSION", 7))
3158 gmx_fatal(FARGS, "Can not read file %s,\n"
3159 " this file is from a Gromacs version which is older than 2.0\n"
3160 " Make a new one with grompp or use a gro or pdb file, if possible",
3161 gmx_fio_getname(fio));
3163 gmx_fio_do_int(fio, precision);
3164 bDouble = (precision == sizeof(double));
3165 if ((precision != sizeof(float)) && !bDouble)
3167 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3168 "instead of %d or %d",
3169 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3171 gmx_fio_setprecision(fio, bDouble);
3172 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3173 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3177 gmx_fio_write_string(fio, GromacsVersion());
3178 bDouble = (precision == sizeof(double));
3179 gmx_fio_setprecision(fio, bDouble);
3180 gmx_fio_do_int(fio, precision);
3182 sprintf(file_tag, "%s", tpx_tag);
3183 fgen = tpx_generation;
3186 /* Check versions! */
3187 gmx_fio_do_int(fio, fver);
3189 /* This is for backward compatibility with development versions 77-79
3190 * where the tag was, mistakenly, placed before the generation,
3191 * which would cause a segv instead of a proper error message
3192 * when reading the topology only from tpx with <77 code.
3194 if (fver >= 77 && fver <= 79)
3196 gmx_fio_do_string(fio, file_tag);
3201 gmx_fio_do_int(fio, fgen);
3210 gmx_fio_do_string(fio, file_tag);
3216 /* Versions before 77 don't have the tag, set it to release */
3217 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3220 if (strcmp(file_tag, tpx_tag) != 0)
3222 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3225 /* We only support reading tpx files with the same tag as the code
3226 * or tpx files with the release tag and with lower version number.
3228 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
3230 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3231 gmx_fio_getname(fio), fver, file_tag,
3232 tpx_version, tpx_tag);
3237 if (file_version != NULL)
3239 *file_version = fver;
3241 if (file_generation != NULL)
3243 *file_generation = fgen;
3247 if ((fver <= tpx_incompatible_version) ||
3248 ((fver > tpx_version) && !TopOnlyOK) ||
3249 (fgen > tpx_generation) ||
3250 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3252 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3253 gmx_fio_getname(fio), fver, tpx_version);
3256 do_section(fio, eitemHEADER, bRead);
3257 gmx_fio_do_int(fio, tpx->natoms);
3260 gmx_fio_do_int(fio, tpx->ngtc);
3268 gmx_fio_do_int(fio, idum);
3269 gmx_fio_do_real(fio, rdum);
3271 /*a better decision will eventually (5.0 or later) need to be made
3272 on how to treat the alchemical state of the system, which can now
3273 vary through a simulation, and cannot be completely described
3274 though a single lambda variable, or even a single state
3275 index. Eventually, should probably be a vector. MRS*/
3278 gmx_fio_do_int(fio, tpx->fep_state);
3280 gmx_fio_do_real(fio, tpx->lambda);
3281 gmx_fio_do_int(fio, tpx->bIr);
3282 gmx_fio_do_int(fio, tpx->bTop);
3283 gmx_fio_do_int(fio, tpx->bX);
3284 gmx_fio_do_int(fio, tpx->bV);
3285 gmx_fio_do_int(fio, tpx->bF);
3286 gmx_fio_do_int(fio, tpx->bBox);
3288 if ((fgen > tpx_generation))
3290 /* This can only happen if TopOnlyOK=TRUE */
3295 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3296 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3297 gmx_bool bXVallocated)
3303 int file_version, file_generation;
3307 gmx_bool bPeriodicMols;
3311 tpx.natoms = state->natoms;
3312 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3313 tpx.fep_state = state->fep_state;
3314 tpx.lambda = state->lambda[efptFEP];
3315 tpx.bIr = (ir != NULL);
3316 tpx.bTop = (mtop != NULL);
3317 tpx.bX = (state->x != NULL);
3318 tpx.bV = (state->v != NULL);
3319 tpx.bF = (f != NULL);
3323 TopOnlyOK = (ir == NULL);
3325 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3330 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3331 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3336 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3337 state->natoms = tpx.natoms;
3338 state->nalloc = tpx.natoms;
3344 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3348 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3350 do_test(fio, tpx.bBox, state->box);
3351 do_section(fio, eitemBOX, bRead);
3354 gmx_fio_ndo_rvec(fio, state->box, DIM);
3355 if (file_version >= 51)
3357 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3361 /* We initialize box_rel after reading the inputrec */
3362 clear_mat(state->box_rel);
3364 if (file_version >= 28)
3366 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3367 if (file_version < 56)
3370 gmx_fio_ndo_rvec(fio, mdum, DIM);
3375 if (state->ngtc > 0 && file_version >= 28)
3378 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3379 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3380 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3381 snew(dumv, state->ngtc);
3382 if (file_version < 69)
3384 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3386 /* These used to be the Berendsen tcoupl_lambda's */
3387 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3391 /* Prior to tpx version 26, the inputrec was here.
3392 * I moved it to enable partial forward-compatibility
3393 * for analysis/viewer programs.
3395 if (file_version < 26)
3397 do_test(fio, tpx.bIr, ir);
3398 do_section(fio, eitemIR, bRead);
3403 do_inputrec(fio, ir, bRead, file_version,
3404 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3407 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3412 do_inputrec(fio, &dum_ir, bRead, file_version,
3413 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3416 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3418 done_inputrec(&dum_ir);
3424 do_test(fio, tpx.bTop, mtop);
3425 do_section(fio, eitemTOP, bRead);
3430 do_mtop(fio, mtop, bRead, file_version);
3434 do_mtop(fio, &dum_top, bRead, file_version);
3435 done_mtop(&dum_top, TRUE);
3438 do_test(fio, tpx.bX, state->x);
3439 do_section(fio, eitemX, bRead);
3444 state->flags |= (1<<estX);
3446 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3449 do_test(fio, tpx.bV, state->v);
3450 do_section(fio, eitemV, bRead);
3455 state->flags |= (1<<estV);
3457 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3460 do_test(fio, tpx.bF, f);
3461 do_section(fio, eitemF, bRead);
3464 gmx_fio_ndo_rvec(fio, f, state->natoms);
3467 /* Starting with tpx version 26, we have the inputrec
3468 * at the end of the file, so we can ignore it
3469 * if the file is never than the software (but still the
3470 * same generation - see comments at the top of this file.
3475 bPeriodicMols = FALSE;
3476 if (file_version >= 26)
3478 do_test(fio, tpx.bIr, ir);
3479 do_section(fio, eitemIR, bRead);
3482 if (file_version >= 53)
3484 /* Removed the pbc info from do_inputrec, since we always want it */
3488 bPeriodicMols = ir->bPeriodicMols;
3490 gmx_fio_do_int(fio, ePBC);
3491 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3493 if (file_generation <= tpx_generation && ir)
3495 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3498 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3500 if (file_version < 51)
3502 set_box_rel(ir, state);
3504 if (file_version < 53)
3507 bPeriodicMols = ir->bPeriodicMols;
3510 if (bRead && ir && file_version >= 53)
3512 /* We need to do this after do_inputrec, since that initializes ir */
3514 ir->bPeriodicMols = bPeriodicMols;
3523 if (state->ngtc == 0)
3525 /* Reading old version without tcoupl state data: set it */
3526 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3528 if (tpx.bTop && mtop)
3530 if (file_version < 57)
3532 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3534 ir->eDisre = edrSimple;
3538 ir->eDisre = edrNone;
3541 set_disres_npair(mtop);
3545 if (tpx.bTop && mtop)
3547 gmx_mtop_finalize(mtop);
3550 if (file_version >= 57)
3554 env = getenv("GMX_NOCHARGEGROUPS");
3557 sscanf(env, "%d", &ienv);
3558 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3563 "Will make single atomic charge groups in non-solvent%s\n",
3564 ienv > 1 ? " and solvent" : "");
3565 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3567 fprintf(stderr, "\n");
3575 /************************************************************
3577 * The following routines are the exported ones
3579 ************************************************************/
3581 t_fileio *open_tpx(const char *fn, const char *mode)
3583 return gmx_fio_open(fn, mode);
3586 void close_tpx(t_fileio *fio)
3591 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3592 int *file_version, int *file_generation)
3596 fio = open_tpx(fn, "r");
3597 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3601 void write_tpx_state(const char *fn,
3602 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3606 fio = open_tpx(fn, "w");
3607 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3611 void read_tpx_state(const char *fn,
3612 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3616 fio = open_tpx(fn, "r");
3617 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3621 int read_tpx(const char *fn,
3622 t_inputrec *ir, matrix box, int *natoms,
3623 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3631 fio = open_tpx(fn, "r");
3632 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3634 *natoms = state.natoms;
3637 copy_mat(state.box, box);
3646 int read_tpx_top(const char *fn,
3647 t_inputrec *ir, matrix box, int *natoms,
3648 rvec *x, rvec *v, rvec *f, t_topology *top)
3654 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3656 *top = gmx_mtop_t_to_t_topology(&mtop);
3661 gmx_bool fn2bTPX(const char *file)
3663 switch (fn2ftp(file))
3674 static void done_gmx_groups_t(gmx_groups_t *g)
3678 for (i = 0; (i < egcNR); i++)
3680 if (NULL != g->grps[i].nm_ind)
3682 sfree(g->grps[i].nm_ind);
3683 g->grps[i].nm_ind = NULL;
3685 if (NULL != g->grpnr[i])
3691 /* The contents of this array is in symtab, don't free it here */
3695 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3696 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3699 int natoms, i, version, generation;
3700 gmx_bool bTop, bXNULL = FALSE;
3702 t_topology *topconv;
3705 bTop = fn2bTPX(infile);
3709 read_tpxheader(infile, &header, TRUE, &version, &generation);
3712 snew(*x, header.natoms);
3716 snew(*v, header.natoms);
3719 *ePBC = read_tpx(infile, NULL, box, &natoms,
3720 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3721 *top = gmx_mtop_t_to_t_topology(mtop);
3722 /* In this case we need to throw away the group data too */
3723 done_gmx_groups_t(&mtop->groups);
3725 strcpy(title, *top->name);
3726 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3730 get_stx_coordnum(infile, &natoms);
3731 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3742 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3750 aps = gmx_atomprop_init();
3751 for (i = 0; (i < natoms); i++)
3753 if (!gmx_atomprop_query(aps, epropMass,
3754 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3755 *top->atoms.atomname[i],
3756 &(top->atoms.atom[i].m)))
3760 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3761 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3762 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3763 *top->atoms.atomname[i]);
3767 gmx_atomprop_destroy(aps);
3769 top->idef.ntypes = -1;