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 "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/futil.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) */
94 /*! \brief Version number of the file format written to run input
95 * files by this version of the code.
97 * The tpx_version number should be increased whenever the file format
98 * in the main development branch changes, generally to the highest
99 * value present in tpxv. Backward compatibility for reading old run
100 * input files is maintained by checking this version number against
101 * that of the file and then using the correct code path.
103 * When developing a feature branch that needs to change the run input
104 * file format, change tpx_tag instead. */
105 static const int tpx_version = tpxv_InteractiveMolecularDynamics;
108 /* This number should only be increased when you edit the TOPOLOGY section
109 * or the HEADER of the tpx format.
110 * This way we can maintain forward compatibility too for all analysis tools
111 * and/or external programs that only need to know the atom/residue names,
112 * charges, and bond connectivity.
114 * It first appeared in tpx version 26, when I also moved the inputrecord
115 * to the end of the tpx file, so we can just skip it if we only
118 * In particular, it must be increased when adding new elements to
119 * ftupd, so that old code can read new .tpr files.
121 static const int tpx_generation = 26;
123 /* This number should be the most recent backwards incompatible version
124 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
126 static const int tpx_incompatible_version = 9;
130 /* Struct used to maintain tpx compatibility when function types are added */
132 int fvnr; /* file version number in which the function type first appeared */
133 int ftype; /* function type */
137 * The entries should be ordered in:
138 * 1. ascending file version number
139 * 2. ascending function type number
141 /*static const t_ftupd ftupd[] = {
142 { 20, F_CUBICBONDS },
146 { 22, F_DISRESVIOL },
152 { 26, F_DIHRESVIOL },
153 { 30, F_CROSS_BOND_BONDS },
154 { 30, F_CROSS_BOND_ANGLES },
155 { 30, F_UREY_BRADLEY },
156 { 30, F_POLARIZATION },
160 * The entries should be ordered in:
161 * 1. ascending function type number
162 * 2. ascending file version number
164 /* question; what is the purpose of the commented code above? */
165 static const t_ftupd ftupd[] = {
166 { 20, F_CUBICBONDS },
171 { 43, F_TABBONDSNC },
172 { 70, F_RESTRBONDS },
173 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
174 { 76, F_LINEAR_ANGLES },
175 { 30, F_CROSS_BOND_BONDS },
176 { 30, F_CROSS_BOND_ANGLES },
177 { 30, F_UREY_BRADLEY },
178 { 34, F_QUARTIC_ANGLES },
180 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
181 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
190 { 72, F_NPSOLVATION },
192 { 41, F_LJC_PAIRS_NB },
195 { 32, F_COUL_RECIP },
198 { 30, F_POLARIZATION },
201 { 22, F_DISRESVIOL },
205 { 26, F_DIHRESVIOL },
210 { 46, F_ECONSERVED },
211 { 69, F_VTEMP_NOLONGERUSED},
213 { 54, F_DVDL_CONSTR },
214 { 76, F_ANHARM_POL },
217 { 79, F_DVDL_BONDED, },
218 { 79, F_DVDL_RESTRAINT },
219 { 79, F_DVDL_TEMPERATURE },
221 #define NFTUPD asize(ftupd)
223 /* Needed for backward compatibility */
226 static void _do_section(t_fileio *fio, int key, gmx_bool bRead, const char *src,
232 if (gmx_fio_getftp(fio) == efTPA)
236 gmx_fio_write_string(fio, itemstr[key]);
237 bDbg = gmx_fio_getdebug(fio);
238 gmx_fio_setdebug(fio, FALSE);
239 gmx_fio_write_string(fio, comment_str[key]);
240 gmx_fio_setdebug(fio, bDbg);
244 if (gmx_fio_getdebug(fio))
246 fprintf(stderr, "Looking for section %s (%s, %d)",
247 itemstr[key], src, line);
252 gmx_fio_do_string(fio, buf);
254 while ((gmx_strcasecmp(buf, itemstr[key]) != 0));
256 if (gmx_strcasecmp(buf, itemstr[key]) != 0)
258 gmx_fatal(FARGS, "\nCould not find section heading %s", itemstr[key]);
260 else if (gmx_fio_getdebug(fio))
262 fprintf(stderr, " and found it\n");
268 #define do_section(fio, key, bRead) _do_section(fio, key, bRead, __FILE__, __LINE__)
270 /**************************************************************
272 * Now the higer level routines that do io of the structures and arrays
274 **************************************************************/
275 static void do_pullgrp_tpx_pre95(t_fileio *fio,
284 gmx_fio_do_int(fio, pgrp->nat);
287 snew(pgrp->ind, pgrp->nat);
289 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
290 gmx_fio_do_int(fio, pgrp->nweight);
293 snew(pgrp->weight, pgrp->nweight);
295 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
296 gmx_fio_do_int(fio, pgrp->pbcatom);
297 gmx_fio_do_rvec(fio, pcrd->vec);
298 clear_rvec(pcrd->origin);
299 gmx_fio_do_rvec(fio, tmp);
301 gmx_fio_do_real(fio, pcrd->rate);
302 gmx_fio_do_real(fio, pcrd->k);
303 if (file_version >= 56)
305 gmx_fio_do_real(fio, pcrd->kB);
313 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
317 gmx_fio_do_int(fio, pgrp->nat);
320 snew(pgrp->ind, pgrp->nat);
322 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
323 gmx_fio_do_int(fio, pgrp->nweight);
326 snew(pgrp->weight, pgrp->nweight);
328 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
329 gmx_fio_do_int(fio, pgrp->pbcatom);
332 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd)
336 gmx_fio_do_int(fio, pcrd->group[0]);
337 gmx_fio_do_int(fio, pcrd->group[1]);
338 gmx_fio_do_rvec(fio, pcrd->origin);
339 gmx_fio_do_rvec(fio, pcrd->vec);
340 gmx_fio_do_real(fio, pcrd->init);
341 gmx_fio_do_real(fio, pcrd->rate);
342 gmx_fio_do_real(fio, pcrd->k);
343 gmx_fio_do_real(fio, pcrd->kB);
346 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
348 /* i is used in the ndo_double macro*/
352 int n_lambda = fepvals->n_lambda;
354 /* reset the lambda calculation window */
355 fepvals->lambda_start_n = 0;
356 fepvals->lambda_stop_n = n_lambda;
357 if (file_version >= 79)
363 snew(expand->init_lambda_weights, n_lambda);
365 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
366 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
369 gmx_fio_do_int(fio, expand->nstexpanded);
370 gmx_fio_do_int(fio, expand->elmcmove);
371 gmx_fio_do_int(fio, expand->elamstats);
372 gmx_fio_do_int(fio, expand->lmc_repeats);
373 gmx_fio_do_int(fio, expand->gibbsdeltalam);
374 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
375 gmx_fio_do_int(fio, expand->lmc_seed);
376 gmx_fio_do_real(fio, expand->mc_temp);
377 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
378 gmx_fio_do_int(fio, expand->nstTij);
379 gmx_fio_do_int(fio, expand->minvarmin);
380 gmx_fio_do_int(fio, expand->c_range);
381 gmx_fio_do_real(fio, expand->wl_scale);
382 gmx_fio_do_real(fio, expand->wl_ratio);
383 gmx_fio_do_real(fio, expand->init_wl_delta);
384 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
385 gmx_fio_do_int(fio, expand->elmceq);
386 gmx_fio_do_int(fio, expand->equil_steps);
387 gmx_fio_do_int(fio, expand->equil_samples);
388 gmx_fio_do_int(fio, expand->equil_n_at_lam);
389 gmx_fio_do_real(fio, expand->equil_wl_delta);
390 gmx_fio_do_real(fio, expand->equil_ratio);
394 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
397 if (file_version >= 79)
399 gmx_fio_do_int(fio, simtemp->eSimTempScale);
400 gmx_fio_do_real(fio, simtemp->simtemp_high);
401 gmx_fio_do_real(fio, simtemp->simtemp_low);
406 snew(simtemp->temperatures, n_lambda);
408 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
413 static void do_imd(t_fileio *fio, t_IMD *imd, gmx_bool bRead)
415 gmx_fio_do_int(fio, imd->nat);
418 snew(imd->ind, imd->nat);
420 gmx_fio_ndo_int(fio, imd->ind, imd->nat);
423 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
425 /* i is defined in the ndo_double macro; use g to iterate. */
430 /* free energy values */
432 if (file_version >= 79)
434 gmx_fio_do_int(fio, fepvals->init_fep_state);
435 gmx_fio_do_double(fio, fepvals->init_lambda);
436 gmx_fio_do_double(fio, fepvals->delta_lambda);
438 else if (file_version >= 59)
440 gmx_fio_do_double(fio, fepvals->init_lambda);
441 gmx_fio_do_double(fio, fepvals->delta_lambda);
445 gmx_fio_do_real(fio, rdum);
446 fepvals->init_lambda = rdum;
447 gmx_fio_do_real(fio, rdum);
448 fepvals->delta_lambda = rdum;
450 if (file_version >= 79)
452 gmx_fio_do_int(fio, fepvals->n_lambda);
455 snew(fepvals->all_lambda, efptNR);
457 for (g = 0; g < efptNR; g++)
459 if (fepvals->n_lambda > 0)
463 snew(fepvals->all_lambda[g], fepvals->n_lambda);
465 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
466 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
468 else if (fepvals->init_lambda >= 0)
470 fepvals->separate_dvdl[efptFEP] = TRUE;
474 else if (file_version >= 64)
476 gmx_fio_do_int(fio, fepvals->n_lambda);
481 snew(fepvals->all_lambda, efptNR);
482 /* still allocate the all_lambda array's contents. */
483 for (g = 0; g < efptNR; g++)
485 if (fepvals->n_lambda > 0)
487 snew(fepvals->all_lambda[g], fepvals->n_lambda);
491 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
493 if (fepvals->init_lambda >= 0)
497 fepvals->separate_dvdl[efptFEP] = TRUE;
501 /* copy the contents of the efptFEP lambda component to all
502 the other components */
503 for (g = 0; g < efptNR; g++)
505 for (h = 0; h < fepvals->n_lambda; h++)
509 fepvals->all_lambda[g][h] =
510 fepvals->all_lambda[efptFEP][h];
519 fepvals->n_lambda = 0;
520 fepvals->all_lambda = NULL;
521 if (fepvals->init_lambda >= 0)
523 fepvals->separate_dvdl[efptFEP] = TRUE;
526 if (file_version >= 13)
528 gmx_fio_do_real(fio, fepvals->sc_alpha);
532 fepvals->sc_alpha = 0;
534 if (file_version >= 38)
536 gmx_fio_do_int(fio, fepvals->sc_power);
540 fepvals->sc_power = 2;
542 if (file_version >= 79)
544 gmx_fio_do_real(fio, fepvals->sc_r_power);
548 fepvals->sc_r_power = 6.0;
550 if (file_version >= 15)
552 gmx_fio_do_real(fio, fepvals->sc_sigma);
556 fepvals->sc_sigma = 0.3;
560 if (file_version >= 71)
562 fepvals->sc_sigma_min = fepvals->sc_sigma;
566 fepvals->sc_sigma_min = 0;
569 if (file_version >= 79)
571 gmx_fio_do_int(fio, fepvals->bScCoul);
575 fepvals->bScCoul = TRUE;
577 if (file_version >= 64)
579 gmx_fio_do_int(fio, fepvals->nstdhdl);
583 fepvals->nstdhdl = 1;
586 if (file_version >= 73)
588 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
589 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
593 fepvals->separate_dhdl_file = esepdhdlfileYES;
594 fepvals->dhdl_derivatives = edhdlderivativesYES;
596 if (file_version >= 71)
598 gmx_fio_do_int(fio, fepvals->dh_hist_size);
599 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
603 fepvals->dh_hist_size = 0;
604 fepvals->dh_hist_spacing = 0.1;
606 if (file_version >= 79)
608 gmx_fio_do_int(fio, fepvals->bPrintEnergy);
612 fepvals->bPrintEnergy = FALSE;
615 /* handle lambda_neighbors */
616 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
618 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
619 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
620 (fepvals->init_lambda < 0) )
622 fepvals->lambda_start_n = (fepvals->init_fep_state -
623 fepvals->lambda_neighbors);
624 fepvals->lambda_stop_n = (fepvals->init_fep_state +
625 fepvals->lambda_neighbors + 1);
626 if (fepvals->lambda_start_n < 0)
628 fepvals->lambda_start_n = 0;;
630 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
632 fepvals->lambda_stop_n = fepvals->n_lambda;
637 fepvals->lambda_start_n = 0;
638 fepvals->lambda_stop_n = fepvals->n_lambda;
643 fepvals->lambda_start_n = 0;
644 fepvals->lambda_stop_n = fepvals->n_lambda;
648 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
652 if (file_version >= 95)
654 gmx_fio_do_int(fio, pull->ngroup);
656 gmx_fio_do_int(fio, pull->ncoord);
657 if (file_version < 95)
659 pull->ngroup = pull->ncoord + 1;
661 gmx_fio_do_int(fio, pull->eGeom);
662 gmx_fio_do_ivec(fio, pull->dim);
663 gmx_fio_do_real(fio, pull->cyl_r1);
664 gmx_fio_do_real(fio, pull->cyl_r0);
665 gmx_fio_do_real(fio, pull->constr_tol);
666 if (file_version >= 95)
668 gmx_fio_do_int(fio, pull->bPrintRef);
670 gmx_fio_do_int(fio, pull->nstxout);
671 gmx_fio_do_int(fio, pull->nstfout);
674 snew(pull->group, pull->ngroup);
675 snew(pull->coord, pull->ncoord);
677 if (file_version < 95)
679 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
680 if (pull->eGeom == epullgDIRPBC)
682 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
684 if (pull->eGeom > epullgDIRPBC)
689 for (g = 0; g < pull->ngroup; g++)
691 /* We read and ignore a pull coordinate for group 0 */
692 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[max(g-1, 0)],
693 bRead, file_version);
696 pull->coord[g-1].group[0] = 0;
697 pull->coord[g-1].group[1] = g;
701 pull->bPrintRef = (pull->group[0].nat > 0);
705 for (g = 0; g < pull->ngroup; g++)
707 do_pull_group(fio, &pull->group[g], bRead);
709 for (g = 0; g < pull->ncoord; g++)
711 do_pull_coord(fio, &pull->coord[g]);
717 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
721 gmx_fio_do_int(fio, rotg->eType);
722 gmx_fio_do_int(fio, rotg->bMassW);
723 gmx_fio_do_int(fio, rotg->nat);
726 snew(rotg->ind, rotg->nat);
728 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
731 snew(rotg->x_ref, rotg->nat);
733 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
734 gmx_fio_do_rvec(fio, rotg->vec);
735 gmx_fio_do_rvec(fio, rotg->pivot);
736 gmx_fio_do_real(fio, rotg->rate);
737 gmx_fio_do_real(fio, rotg->k);
738 gmx_fio_do_real(fio, rotg->slab_dist);
739 gmx_fio_do_real(fio, rotg->min_gaussian);
740 gmx_fio_do_real(fio, rotg->eps);
741 gmx_fio_do_int(fio, rotg->eFittype);
742 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
743 gmx_fio_do_real(fio, rotg->PotAngle_step);
746 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
750 gmx_fio_do_int(fio, rot->ngrp);
751 gmx_fio_do_int(fio, rot->nstrout);
752 gmx_fio_do_int(fio, rot->nstsout);
755 snew(rot->grp, rot->ngrp);
757 for (g = 0; g < rot->ngrp; g++)
759 do_rotgrp(fio, &rot->grp[g], bRead);
764 static void do_swapcoords(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead)
769 gmx_fio_do_int(fio, swap->nat);
770 gmx_fio_do_int(fio, swap->nat_sol);
771 for (j = 0; j < 2; j++)
773 gmx_fio_do_int(fio, swap->nat_split[j]);
774 gmx_fio_do_int(fio, swap->massw_split[j]);
776 gmx_fio_do_int(fio, swap->nstswap);
777 gmx_fio_do_int(fio, swap->nAverage);
778 gmx_fio_do_real(fio, swap->threshold);
779 gmx_fio_do_real(fio, swap->cyl0r);
780 gmx_fio_do_real(fio, swap->cyl0u);
781 gmx_fio_do_real(fio, swap->cyl0l);
782 gmx_fio_do_real(fio, swap->cyl1r);
783 gmx_fio_do_real(fio, swap->cyl1u);
784 gmx_fio_do_real(fio, swap->cyl1l);
788 snew(swap->ind, swap->nat);
789 snew(swap->ind_sol, swap->nat_sol);
790 for (j = 0; j < 2; j++)
792 snew(swap->ind_split[j], swap->nat_split[j]);
796 gmx_fio_ndo_int(fio, swap->ind, swap->nat);
797 gmx_fio_ndo_int(fio, swap->ind_sol, swap->nat_sol);
798 for (j = 0; j < 2; j++)
800 gmx_fio_ndo_int(fio, swap->ind_split[j], swap->nat_split[j]);
803 for (j = 0; j < eCompNR; j++)
805 gmx_fio_do_int(fio, swap->nanions[j]);
806 gmx_fio_do_int(fio, swap->ncations[j]);
812 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
813 int file_version, real *fudgeQQ)
815 int i, j, k, *tmp, idum = 0;
819 real zerotemptime, finish_t, init_temp, finish_temp;
821 if (file_version != tpx_version)
823 /* Give a warning about features that are not accessible */
824 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
825 file_version, tpx_version);
833 if (file_version == 0)
838 /* Basic inputrec stuff */
839 gmx_fio_do_int(fio, ir->eI);
840 if (file_version >= 62)
842 gmx_fio_do_int64(fio, ir->nsteps);
846 gmx_fio_do_int(fio, idum);
849 if (file_version > 25)
851 if (file_version >= 62)
853 gmx_fio_do_int64(fio, ir->init_step);
857 gmx_fio_do_int(fio, idum);
858 ir->init_step = idum;
866 if (file_version >= 58)
868 gmx_fio_do_int(fio, ir->simulation_part);
872 ir->simulation_part = 1;
875 if (file_version >= 67)
877 gmx_fio_do_int(fio, ir->nstcalcenergy);
881 ir->nstcalcenergy = 1;
883 if (file_version < 53)
885 /* The pbc info has been moved out of do_inputrec,
886 * since we always want it, also without reading the inputrec.
888 gmx_fio_do_int(fio, ir->ePBC);
889 if ((file_version <= 15) && (ir->ePBC == 2))
893 if (file_version >= 45)
895 gmx_fio_do_int(fio, ir->bPeriodicMols);
902 ir->bPeriodicMols = TRUE;
906 ir->bPeriodicMols = FALSE;
910 if (file_version >= 81)
912 gmx_fio_do_int(fio, ir->cutoff_scheme);
913 if (file_version < 94)
915 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
920 ir->cutoff_scheme = ecutsGROUP;
922 gmx_fio_do_int(fio, ir->ns_type);
923 gmx_fio_do_int(fio, ir->nstlist);
924 gmx_fio_do_int(fio, ir->ndelta);
925 if (file_version < 41)
927 gmx_fio_do_int(fio, idum);
928 gmx_fio_do_int(fio, idum);
930 if (file_version >= 45)
932 gmx_fio_do_real(fio, ir->rtpi);
938 gmx_fio_do_int(fio, ir->nstcomm);
939 if (file_version > 34)
941 gmx_fio_do_int(fio, ir->comm_mode);
943 else if (ir->nstcomm < 0)
945 ir->comm_mode = ecmANGULAR;
949 ir->comm_mode = ecmLINEAR;
951 ir->nstcomm = abs(ir->nstcomm);
953 if (file_version > 25)
955 gmx_fio_do_int(fio, ir->nstcheckpoint);
959 ir->nstcheckpoint = 0;
962 gmx_fio_do_int(fio, ir->nstcgsteep);
964 if (file_version >= 30)
966 gmx_fio_do_int(fio, ir->nbfgscorr);
973 gmx_fio_do_int(fio, ir->nstlog);
974 gmx_fio_do_int(fio, ir->nstxout);
975 gmx_fio_do_int(fio, ir->nstvout);
976 gmx_fio_do_int(fio, ir->nstfout);
977 gmx_fio_do_int(fio, ir->nstenergy);
978 gmx_fio_do_int(fio, ir->nstxout_compressed);
979 if (file_version >= 59)
981 gmx_fio_do_double(fio, ir->init_t);
982 gmx_fio_do_double(fio, ir->delta_t);
986 gmx_fio_do_real(fio, rdum);
988 gmx_fio_do_real(fio, rdum);
991 gmx_fio_do_real(fio, ir->x_compression_precision);
992 if (file_version < 19)
994 gmx_fio_do_int(fio, idum);
995 gmx_fio_do_int(fio, idum);
997 if (file_version < 18)
999 gmx_fio_do_int(fio, idum);
1001 if (file_version >= 81)
1003 gmx_fio_do_real(fio, ir->verletbuf_tol);
1007 ir->verletbuf_tol = 0;
1009 gmx_fio_do_real(fio, ir->rlist);
1010 if (file_version >= 67)
1012 gmx_fio_do_real(fio, ir->rlistlong);
1014 if (file_version >= 82 && file_version != 90)
1016 gmx_fio_do_int(fio, ir->nstcalclr);
1020 /* Calculate at NS steps */
1021 ir->nstcalclr = ir->nstlist;
1023 gmx_fio_do_int(fio, ir->coulombtype);
1024 if (file_version < 32 && ir->coulombtype == eelRF)
1026 ir->coulombtype = eelRF_NEC;
1028 if (file_version >= 81)
1030 gmx_fio_do_int(fio, ir->coulomb_modifier);
1034 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1036 gmx_fio_do_real(fio, ir->rcoulomb_switch);
1037 gmx_fio_do_real(fio, ir->rcoulomb);
1038 gmx_fio_do_int(fio, ir->vdwtype);
1039 if (file_version >= 81)
1041 gmx_fio_do_int(fio, ir->vdw_modifier);
1045 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1047 gmx_fio_do_real(fio, ir->rvdw_switch);
1048 gmx_fio_do_real(fio, ir->rvdw);
1049 if (file_version < 67)
1051 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1053 gmx_fio_do_int(fio, ir->eDispCorr);
1054 gmx_fio_do_real(fio, ir->epsilon_r);
1055 if (file_version >= 37)
1057 gmx_fio_do_real(fio, ir->epsilon_rf);
1061 if (EEL_RF(ir->coulombtype))
1063 ir->epsilon_rf = ir->epsilon_r;
1064 ir->epsilon_r = 1.0;
1068 ir->epsilon_rf = 1.0;
1071 if (file_version >= 29)
1073 gmx_fio_do_real(fio, ir->tabext);
1080 if (file_version > 25)
1082 gmx_fio_do_int(fio, ir->gb_algorithm);
1083 gmx_fio_do_int(fio, ir->nstgbradii);
1084 gmx_fio_do_real(fio, ir->rgbradii);
1085 gmx_fio_do_real(fio, ir->gb_saltconc);
1086 gmx_fio_do_int(fio, ir->implicit_solvent);
1090 ir->gb_algorithm = egbSTILL;
1093 ir->gb_saltconc = 0;
1094 ir->implicit_solvent = eisNO;
1096 if (file_version >= 55)
1098 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1099 gmx_fio_do_real(fio, ir->gb_obc_alpha);
1100 gmx_fio_do_real(fio, ir->gb_obc_beta);
1101 gmx_fio_do_real(fio, ir->gb_obc_gamma);
1102 if (file_version >= 60)
1104 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1105 gmx_fio_do_int(fio, ir->sa_algorithm);
1109 ir->gb_dielectric_offset = 0.009;
1110 ir->sa_algorithm = esaAPPROX;
1112 gmx_fio_do_real(fio, ir->sa_surface_tension);
1114 /* Override sa_surface_tension if it is not changed in the mpd-file */
1115 if (ir->sa_surface_tension < 0)
1117 if (ir->gb_algorithm == egbSTILL)
1119 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1121 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1123 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1130 /* Better use sensible values than insane (0.0) ones... */
1131 ir->gb_epsilon_solvent = 80;
1132 ir->gb_obc_alpha = 1.0;
1133 ir->gb_obc_beta = 0.8;
1134 ir->gb_obc_gamma = 4.85;
1135 ir->sa_surface_tension = 2.092;
1139 if (file_version >= 81)
1141 gmx_fio_do_real(fio, ir->fourier_spacing);
1145 ir->fourier_spacing = 0.0;
1147 gmx_fio_do_int(fio, ir->nkx);
1148 gmx_fio_do_int(fio, ir->nky);
1149 gmx_fio_do_int(fio, ir->nkz);
1150 gmx_fio_do_int(fio, ir->pme_order);
1151 gmx_fio_do_real(fio, ir->ewald_rtol);
1153 if (file_version >= 93)
1155 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1159 ir->ewald_rtol_lj = ir->ewald_rtol;
1162 if (file_version >= 24)
1164 gmx_fio_do_int(fio, ir->ewald_geometry);
1168 ir->ewald_geometry = eewg3D;
1171 if (file_version <= 17)
1173 ir->epsilon_surface = 0;
1174 if (file_version == 17)
1176 gmx_fio_do_int(fio, idum);
1181 gmx_fio_do_real(fio, ir->epsilon_surface);
1184 gmx_fio_do_gmx_bool(fio, ir->bOptFFT);
1186 if (file_version >= 93)
1188 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1190 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1191 gmx_fio_do_int(fio, ir->etc);
1192 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1193 * but the values 0 and 1 still mean no and
1194 * berendsen temperature coupling, respectively.
1196 if (file_version >= 79)
1198 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1200 if (file_version >= 71)
1202 gmx_fio_do_int(fio, ir->nsttcouple);
1206 ir->nsttcouple = ir->nstcalcenergy;
1208 if (file_version <= 15)
1210 gmx_fio_do_int(fio, idum);
1212 if (file_version <= 17)
1214 gmx_fio_do_int(fio, ir->epct);
1215 if (file_version <= 15)
1219 ir->epct = epctSURFACETENSION;
1221 gmx_fio_do_int(fio, idum);
1224 /* we have removed the NO alternative at the beginning */
1228 ir->epct = epctISOTROPIC;
1232 ir->epc = epcBERENDSEN;
1237 gmx_fio_do_int(fio, ir->epc);
1238 gmx_fio_do_int(fio, ir->epct);
1240 if (file_version >= 71)
1242 gmx_fio_do_int(fio, ir->nstpcouple);
1246 ir->nstpcouple = ir->nstcalcenergy;
1248 gmx_fio_do_real(fio, ir->tau_p);
1249 if (file_version <= 15)
1251 gmx_fio_do_rvec(fio, vdum);
1252 clear_mat(ir->ref_p);
1253 for (i = 0; i < DIM; i++)
1255 ir->ref_p[i][i] = vdum[i];
1260 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1261 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1262 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1264 if (file_version <= 15)
1266 gmx_fio_do_rvec(fio, vdum);
1267 clear_mat(ir->compress);
1268 for (i = 0; i < DIM; i++)
1270 ir->compress[i][i] = vdum[i];
1275 gmx_fio_do_rvec(fio, ir->compress[XX]);
1276 gmx_fio_do_rvec(fio, ir->compress[YY]);
1277 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1279 if (file_version >= 47)
1281 gmx_fio_do_int(fio, ir->refcoord_scaling);
1282 gmx_fio_do_rvec(fio, ir->posres_com);
1283 gmx_fio_do_rvec(fio, ir->posres_comB);
1287 ir->refcoord_scaling = erscNO;
1288 clear_rvec(ir->posres_com);
1289 clear_rvec(ir->posres_comB);
1291 if ((file_version > 25) && (file_version < 79))
1293 gmx_fio_do_int(fio, ir->andersen_seed);
1297 ir->andersen_seed = 0;
1299 if (file_version < 26)
1301 gmx_fio_do_gmx_bool(fio, bSimAnn);
1302 gmx_fio_do_real(fio, zerotemptime);
1305 if (file_version < 37)
1307 gmx_fio_do_real(fio, rdum);
1310 gmx_fio_do_real(fio, ir->shake_tol);
1311 if (file_version < 54)
1313 gmx_fio_do_real(fio, *fudgeQQ);
1316 gmx_fio_do_int(fio, ir->efep);
1317 if (file_version <= 14 && ir->efep != efepNO)
1321 do_fepvals(fio, ir->fepvals, bRead, file_version);
1323 if (file_version >= 79)
1325 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1328 ir->bSimTemp = TRUE;
1333 ir->bSimTemp = FALSE;
1337 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1340 if (file_version >= 79)
1342 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1345 ir->bExpanded = TRUE;
1349 ir->bExpanded = FALSE;
1354 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1356 if (file_version >= 57)
1358 gmx_fio_do_int(fio, ir->eDisre);
1360 gmx_fio_do_int(fio, ir->eDisreWeighting);
1361 if (file_version < 22)
1363 if (ir->eDisreWeighting == 0)
1365 ir->eDisreWeighting = edrwEqual;
1369 ir->eDisreWeighting = edrwConservative;
1372 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1373 gmx_fio_do_real(fio, ir->dr_fc);
1374 gmx_fio_do_real(fio, ir->dr_tau);
1375 gmx_fio_do_int(fio, ir->nstdisreout);
1376 if (file_version >= 22)
1378 gmx_fio_do_real(fio, ir->orires_fc);
1379 gmx_fio_do_real(fio, ir->orires_tau);
1380 gmx_fio_do_int(fio, ir->nstorireout);
1386 ir->nstorireout = 0;
1388 if (file_version >= 26 && file_version < 79)
1390 gmx_fio_do_real(fio, ir->dihre_fc);
1391 if (file_version < 56)
1393 gmx_fio_do_real(fio, rdum);
1394 gmx_fio_do_int(fio, idum);
1402 gmx_fio_do_real(fio, ir->em_stepsize);
1403 gmx_fio_do_real(fio, ir->em_tol);
1404 if (file_version >= 22)
1406 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1410 ir->bShakeSOR = TRUE;
1412 if (file_version >= 11)
1414 gmx_fio_do_int(fio, ir->niter);
1419 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1422 if (file_version >= 21)
1424 gmx_fio_do_real(fio, ir->fc_stepsize);
1428 ir->fc_stepsize = 0;
1430 gmx_fio_do_int(fio, ir->eConstrAlg);
1431 gmx_fio_do_int(fio, ir->nProjOrder);
1432 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1433 if (file_version <= 14)
1435 gmx_fio_do_int(fio, idum);
1437 if (file_version >= 26)
1439 gmx_fio_do_int(fio, ir->nLincsIter);
1444 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1447 if (file_version < 33)
1449 gmx_fio_do_real(fio, bd_temp);
1451 gmx_fio_do_real(fio, ir->bd_fric);
1452 if (file_version >= tpxv_Use64BitRandomSeed)
1454 gmx_fio_do_int64(fio, ir->ld_seed);
1458 gmx_fio_do_int(fio, idum);
1461 if (file_version >= 33)
1463 for (i = 0; i < DIM; i++)
1465 gmx_fio_do_rvec(fio, ir->deform[i]);
1470 for (i = 0; i < DIM; i++)
1472 clear_rvec(ir->deform[i]);
1475 if (file_version >= 14)
1477 gmx_fio_do_real(fio, ir->cos_accel);
1483 gmx_fio_do_int(fio, ir->userint1);
1484 gmx_fio_do_int(fio, ir->userint2);
1485 gmx_fio_do_int(fio, ir->userint3);
1486 gmx_fio_do_int(fio, ir->userint4);
1487 gmx_fio_do_real(fio, ir->userreal1);
1488 gmx_fio_do_real(fio, ir->userreal2);
1489 gmx_fio_do_real(fio, ir->userreal3);
1490 gmx_fio_do_real(fio, ir->userreal4);
1493 if (file_version >= 77)
1495 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1500 snew(ir->adress, 1);
1502 gmx_fio_do_int(fio, ir->adress->type);
1503 gmx_fio_do_real(fio, ir->adress->const_wf);
1504 gmx_fio_do_real(fio, ir->adress->ex_width);
1505 gmx_fio_do_real(fio, ir->adress->hy_width);
1506 gmx_fio_do_int(fio, ir->adress->icor);
1507 gmx_fio_do_int(fio, ir->adress->site);
1508 gmx_fio_do_rvec(fio, ir->adress->refs);
1509 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1510 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1511 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1512 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1516 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1518 if (ir->adress->n_tf_grps > 0)
1520 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1524 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1526 if (ir->adress->n_energy_grps > 0)
1528 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1534 ir->bAdress = FALSE;
1538 if (file_version >= 48)
1540 gmx_fio_do_int(fio, ir->ePull);
1541 if (ir->ePull != epullNO)
1547 do_pull(fio, ir->pull, bRead, file_version);
1552 ir->ePull = epullNO;
1555 /* Enforced rotation */
1556 if (file_version >= 74)
1558 gmx_fio_do_int(fio, ir->bRot);
1559 if (ir->bRot == TRUE)
1565 do_rot(fio, ir->rot, bRead);
1573 /* Interactive molecular dynamics */
1574 if (file_version >= tpxv_InteractiveMolecularDynamics)
1576 gmx_fio_do_int(fio, ir->bIMD);
1577 if (TRUE == ir->bIMD)
1583 do_imd(fio, ir->imd, bRead);
1588 /* We don't support IMD sessions for old .tpr files */
1593 gmx_fio_do_int(fio, ir->opts.ngtc);
1594 if (file_version >= 69)
1596 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1600 ir->opts.nhchainlength = 1;
1602 gmx_fio_do_int(fio, ir->opts.ngacc);
1603 gmx_fio_do_int(fio, ir->opts.ngfrz);
1604 gmx_fio_do_int(fio, ir->opts.ngener);
1608 snew(ir->opts.nrdf, ir->opts.ngtc);
1609 snew(ir->opts.ref_t, ir->opts.ngtc);
1610 snew(ir->opts.annealing, ir->opts.ngtc);
1611 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1612 snew(ir->opts.anneal_time, ir->opts.ngtc);
1613 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1614 snew(ir->opts.tau_t, ir->opts.ngtc);
1615 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1616 snew(ir->opts.acc, ir->opts.ngacc);
1617 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1619 if (ir->opts.ngtc > 0)
1621 if (bRead && file_version < 13)
1623 snew(tmp, ir->opts.ngtc);
1624 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1625 for (i = 0; i < ir->opts.ngtc; i++)
1627 ir->opts.nrdf[i] = tmp[i];
1633 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1635 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1636 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1637 if (file_version < 33 && ir->eI == eiBD)
1639 for (i = 0; i < ir->opts.ngtc; i++)
1641 ir->opts.tau_t[i] = bd_temp;
1645 if (ir->opts.ngfrz > 0)
1647 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1649 if (ir->opts.ngacc > 0)
1651 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1653 if (file_version >= 12)
1655 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1656 ir->opts.ngener*ir->opts.ngener);
1659 if (bRead && file_version < 26)
1661 for (i = 0; i < ir->opts.ngtc; i++)
1665 ir->opts.annealing[i] = eannSINGLE;
1666 ir->opts.anneal_npoints[i] = 2;
1667 snew(ir->opts.anneal_time[i], 2);
1668 snew(ir->opts.anneal_temp[i], 2);
1669 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1670 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1671 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1672 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1673 ir->opts.anneal_time[i][0] = ir->init_t;
1674 ir->opts.anneal_time[i][1] = finish_t;
1675 ir->opts.anneal_temp[i][0] = init_temp;
1676 ir->opts.anneal_temp[i][1] = finish_temp;
1680 ir->opts.annealing[i] = eannNO;
1681 ir->opts.anneal_npoints[i] = 0;
1687 /* file version 26 or later */
1688 /* First read the lists with annealing and npoints for each group */
1689 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1690 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1691 for (j = 0; j < (ir->opts.ngtc); j++)
1693 k = ir->opts.anneal_npoints[j];
1696 snew(ir->opts.anneal_time[j], k);
1697 snew(ir->opts.anneal_temp[j], k);
1699 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1700 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1704 if (file_version >= 45)
1706 gmx_fio_do_int(fio, ir->nwall);
1707 gmx_fio_do_int(fio, ir->wall_type);
1708 if (file_version >= 50)
1710 gmx_fio_do_real(fio, ir->wall_r_linpot);
1714 ir->wall_r_linpot = -1;
1716 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1717 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1718 gmx_fio_do_real(fio, ir->wall_density[0]);
1719 gmx_fio_do_real(fio, ir->wall_density[1]);
1720 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1726 ir->wall_atomtype[0] = -1;
1727 ir->wall_atomtype[1] = -1;
1728 ir->wall_density[0] = 0;
1729 ir->wall_density[1] = 0;
1730 ir->wall_ewald_zfac = 3;
1732 /* Cosine stuff for electric fields */
1733 for (j = 0; (j < DIM); j++)
1735 gmx_fio_do_int(fio, ir->ex[j].n);
1736 gmx_fio_do_int(fio, ir->et[j].n);
1739 snew(ir->ex[j].a, ir->ex[j].n);
1740 snew(ir->ex[j].phi, ir->ex[j].n);
1741 snew(ir->et[j].a, ir->et[j].n);
1742 snew(ir->et[j].phi, ir->et[j].n);
1744 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1745 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1746 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1747 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1751 if (file_version >= tpxv_ComputationalElectrophysiology)
1753 gmx_fio_do_int(fio, ir->eSwapCoords);
1754 if (ir->eSwapCoords != eswapNO)
1760 do_swapcoords(fio, ir->swap, bRead);
1765 if (file_version >= 39)
1767 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1768 gmx_fio_do_int(fio, ir->QMMMscheme);
1769 gmx_fio_do_real(fio, ir->scalefactor);
1770 gmx_fio_do_int(fio, ir->opts.ngQM);
1773 snew(ir->opts.QMmethod, ir->opts.ngQM);
1774 snew(ir->opts.QMbasis, ir->opts.ngQM);
1775 snew(ir->opts.QMcharge, ir->opts.ngQM);
1776 snew(ir->opts.QMmult, ir->opts.ngQM);
1777 snew(ir->opts.bSH, ir->opts.ngQM);
1778 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1779 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1780 snew(ir->opts.SAon, ir->opts.ngQM);
1781 snew(ir->opts.SAoff, ir->opts.ngQM);
1782 snew(ir->opts.SAsteps, ir->opts.ngQM);
1783 snew(ir->opts.bOPT, ir->opts.ngQM);
1784 snew(ir->opts.bTS, ir->opts.ngQM);
1786 if (ir->opts.ngQM > 0)
1788 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1789 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1790 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1791 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1792 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1793 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1794 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1795 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1796 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1797 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1798 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1799 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1801 /* end of QMMM stuff */
1806 static void do_harm(t_fileio *fio, t_iparams *iparams)
1808 gmx_fio_do_real(fio, iparams->harmonic.rA);
1809 gmx_fio_do_real(fio, iparams->harmonic.krA);
1810 gmx_fio_do_real(fio, iparams->harmonic.rB);
1811 gmx_fio_do_real(fio, iparams->harmonic.krB);
1814 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1815 gmx_bool bRead, int file_version)
1822 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1832 do_harm(fio, iparams);
1833 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1835 /* Correct incorrect storage of parameters */
1836 iparams->pdihs.phiB = iparams->pdihs.phiA;
1837 iparams->pdihs.cpB = iparams->pdihs.cpA;
1841 gmx_fio_do_real(fio, iparams->harmonic.rA);
1842 gmx_fio_do_real(fio, iparams->harmonic.krA);
1844 case F_LINEAR_ANGLES:
1845 gmx_fio_do_real(fio, iparams->linangle.klinA);
1846 gmx_fio_do_real(fio, iparams->linangle.aA);
1847 gmx_fio_do_real(fio, iparams->linangle.klinB);
1848 gmx_fio_do_real(fio, iparams->linangle.aB);
1851 gmx_fio_do_real(fio, iparams->fene.bm);
1852 gmx_fio_do_real(fio, iparams->fene.kb);
1856 gmx_fio_do_real(fio, iparams->restraint.lowA);
1857 gmx_fio_do_real(fio, iparams->restraint.up1A);
1858 gmx_fio_do_real(fio, iparams->restraint.up2A);
1859 gmx_fio_do_real(fio, iparams->restraint.kA);
1860 gmx_fio_do_real(fio, iparams->restraint.lowB);
1861 gmx_fio_do_real(fio, iparams->restraint.up1B);
1862 gmx_fio_do_real(fio, iparams->restraint.up2B);
1863 gmx_fio_do_real(fio, iparams->restraint.kB);
1869 gmx_fio_do_real(fio, iparams->tab.kA);
1870 gmx_fio_do_int(fio, iparams->tab.table);
1871 gmx_fio_do_real(fio, iparams->tab.kB);
1873 case F_CROSS_BOND_BONDS:
1874 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1875 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1876 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1878 case F_CROSS_BOND_ANGLES:
1879 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1880 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1881 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1882 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1884 case F_UREY_BRADLEY:
1885 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1886 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1887 gmx_fio_do_real(fio, iparams->u_b.r13A);
1888 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1889 if (file_version >= 79)
1891 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1892 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1893 gmx_fio_do_real(fio, iparams->u_b.r13B);
1894 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1898 iparams->u_b.thetaB = iparams->u_b.thetaA;
1899 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1900 iparams->u_b.r13B = iparams->u_b.r13A;
1901 iparams->u_b.kUBB = iparams->u_b.kUBA;
1904 case F_QUARTIC_ANGLES:
1905 gmx_fio_do_real(fio, iparams->qangle.theta);
1906 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1909 gmx_fio_do_real(fio, iparams->bham.a);
1910 gmx_fio_do_real(fio, iparams->bham.b);
1911 gmx_fio_do_real(fio, iparams->bham.c);
1914 gmx_fio_do_real(fio, iparams->morse.b0A);
1915 gmx_fio_do_real(fio, iparams->morse.cbA);
1916 gmx_fio_do_real(fio, iparams->morse.betaA);
1917 if (file_version >= 79)
1919 gmx_fio_do_real(fio, iparams->morse.b0B);
1920 gmx_fio_do_real(fio, iparams->morse.cbB);
1921 gmx_fio_do_real(fio, iparams->morse.betaB);
1925 iparams->morse.b0B = iparams->morse.b0A;
1926 iparams->morse.cbB = iparams->morse.cbA;
1927 iparams->morse.betaB = iparams->morse.betaA;
1931 gmx_fio_do_real(fio, iparams->cubic.b0);
1932 gmx_fio_do_real(fio, iparams->cubic.kb);
1933 gmx_fio_do_real(fio, iparams->cubic.kcub);
1937 case F_POLARIZATION:
1938 gmx_fio_do_real(fio, iparams->polarize.alpha);
1941 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1942 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1943 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1946 if (file_version < 31)
1948 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1950 gmx_fio_do_real(fio, iparams->wpol.al_x);
1951 gmx_fio_do_real(fio, iparams->wpol.al_y);
1952 gmx_fio_do_real(fio, iparams->wpol.al_z);
1953 gmx_fio_do_real(fio, iparams->wpol.rOH);
1954 gmx_fio_do_real(fio, iparams->wpol.rHH);
1955 gmx_fio_do_real(fio, iparams->wpol.rOD);
1958 gmx_fio_do_real(fio, iparams->thole.a);
1959 gmx_fio_do_real(fio, iparams->thole.alpha1);
1960 gmx_fio_do_real(fio, iparams->thole.alpha2);
1961 gmx_fio_do_real(fio, iparams->thole.rfac);
1964 gmx_fio_do_real(fio, iparams->lj.c6);
1965 gmx_fio_do_real(fio, iparams->lj.c12);
1968 gmx_fio_do_real(fio, iparams->lj14.c6A);
1969 gmx_fio_do_real(fio, iparams->lj14.c12A);
1970 gmx_fio_do_real(fio, iparams->lj14.c6B);
1971 gmx_fio_do_real(fio, iparams->lj14.c12B);
1974 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1975 gmx_fio_do_real(fio, iparams->ljc14.qi);
1976 gmx_fio_do_real(fio, iparams->ljc14.qj);
1977 gmx_fio_do_real(fio, iparams->ljc14.c6);
1978 gmx_fio_do_real(fio, iparams->ljc14.c12);
1980 case F_LJC_PAIRS_NB:
1981 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1982 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1983 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1984 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1990 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1991 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1992 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1994 /* Read the incorrectly stored multiplicity */
1995 gmx_fio_do_real(fio, iparams->harmonic.rB);
1996 gmx_fio_do_real(fio, iparams->harmonic.krB);
1997 iparams->pdihs.phiB = iparams->pdihs.phiA;
1998 iparams->pdihs.cpB = iparams->pdihs.cpA;
2002 gmx_fio_do_real(fio, iparams->pdihs.phiB);
2003 gmx_fio_do_real(fio, iparams->pdihs.cpB);
2004 gmx_fio_do_int(fio, iparams->pdihs.mult);
2008 gmx_fio_do_real(fio, iparams->pdihs.phiA);
2009 gmx_fio_do_real(fio, iparams->pdihs.cpA);
2012 gmx_fio_do_int(fio, iparams->disres.label);
2013 gmx_fio_do_int(fio, iparams->disres.type);
2014 gmx_fio_do_real(fio, iparams->disres.low);
2015 gmx_fio_do_real(fio, iparams->disres.up1);
2016 gmx_fio_do_real(fio, iparams->disres.up2);
2017 gmx_fio_do_real(fio, iparams->disres.kfac);
2020 gmx_fio_do_int(fio, iparams->orires.ex);
2021 gmx_fio_do_int(fio, iparams->orires.label);
2022 gmx_fio_do_int(fio, iparams->orires.power);
2023 gmx_fio_do_real(fio, iparams->orires.c);
2024 gmx_fio_do_real(fio, iparams->orires.obs);
2025 gmx_fio_do_real(fio, iparams->orires.kfac);
2028 if (file_version < 82)
2030 gmx_fio_do_int(fio, idum);
2031 gmx_fio_do_int(fio, idum);
2033 gmx_fio_do_real(fio, iparams->dihres.phiA);
2034 gmx_fio_do_real(fio, iparams->dihres.dphiA);
2035 gmx_fio_do_real(fio, iparams->dihres.kfacA);
2036 if (file_version >= 82)
2038 gmx_fio_do_real(fio, iparams->dihres.phiB);
2039 gmx_fio_do_real(fio, iparams->dihres.dphiB);
2040 gmx_fio_do_real(fio, iparams->dihres.kfacB);
2044 iparams->dihres.phiB = iparams->dihres.phiA;
2045 iparams->dihres.dphiB = iparams->dihres.dphiA;
2046 iparams->dihres.kfacB = iparams->dihres.kfacA;
2050 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
2051 gmx_fio_do_rvec(fio, iparams->posres.fcA);
2052 if (bRead && file_version < 27)
2054 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
2055 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
2059 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
2060 gmx_fio_do_rvec(fio, iparams->posres.fcB);
2064 gmx_fio_do_int(fio, iparams->fbposres.geom);
2065 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
2066 gmx_fio_do_real(fio, iparams->fbposres.r);
2067 gmx_fio_do_real(fio, iparams->fbposres.k);
2070 gmx_fio_ndo_real(fio, iparams->cbtdihs.cbtcA, NR_CBTDIHS);
2073 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2074 if (file_version >= 25)
2076 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2080 /* Fourier dihedrals are internally represented
2081 * as Ryckaert-Bellemans since those are faster to compute.
2083 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2084 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2088 gmx_fio_do_real(fio, iparams->constr.dA);
2089 gmx_fio_do_real(fio, iparams->constr.dB);
2092 gmx_fio_do_real(fio, iparams->settle.doh);
2093 gmx_fio_do_real(fio, iparams->settle.dhh);
2096 gmx_fio_do_real(fio, iparams->vsite.a);
2101 gmx_fio_do_real(fio, iparams->vsite.a);
2102 gmx_fio_do_real(fio, iparams->vsite.b);
2107 gmx_fio_do_real(fio, iparams->vsite.a);
2108 gmx_fio_do_real(fio, iparams->vsite.b);
2109 gmx_fio_do_real(fio, iparams->vsite.c);
2112 gmx_fio_do_int(fio, iparams->vsiten.n);
2113 gmx_fio_do_real(fio, iparams->vsiten.a);
2118 /* We got rid of some parameters in version 68 */
2119 if (bRead && file_version < 68)
2121 gmx_fio_do_real(fio, rdum);
2122 gmx_fio_do_real(fio, rdum);
2123 gmx_fio_do_real(fio, rdum);
2124 gmx_fio_do_real(fio, rdum);
2126 gmx_fio_do_real(fio, iparams->gb.sar);
2127 gmx_fio_do_real(fio, iparams->gb.st);
2128 gmx_fio_do_real(fio, iparams->gb.pi);
2129 gmx_fio_do_real(fio, iparams->gb.gbr);
2130 gmx_fio_do_real(fio, iparams->gb.bmlt);
2133 gmx_fio_do_int(fio, iparams->cmap.cmapA);
2134 gmx_fio_do_int(fio, iparams->cmap.cmapB);
2137 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
2138 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
2142 gmx_fio_unset_comment(fio);
2146 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
2153 gmx_fio_set_comment(fio, interaction_function[ftype].name);
2155 if (file_version < 44)
2157 for (i = 0; i < MAXNODES; i++)
2159 gmx_fio_do_int(fio, idum);
2162 gmx_fio_do_int(fio, ilist->nr);
2165 snew(ilist->iatoms, ilist->nr);
2167 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2170 gmx_fio_unset_comment(fio);
2174 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2175 gmx_bool bRead, int file_version)
2180 gmx_fio_do_int(fio, ffparams->atnr);
2181 if (file_version < 57)
2183 gmx_fio_do_int(fio, idum);
2185 gmx_fio_do_int(fio, ffparams->ntypes);
2188 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2189 ffparams->atnr, ffparams->ntypes);
2193 snew(ffparams->functype, ffparams->ntypes);
2194 snew(ffparams->iparams, ffparams->ntypes);
2196 /* Read/write all the function types */
2197 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2200 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2203 if (file_version >= 66)
2205 gmx_fio_do_double(fio, ffparams->reppow);
2209 ffparams->reppow = 12.0;
2212 if (file_version >= 57)
2214 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2217 /* Check whether all these function types are supported by the code.
2218 * In practice the code is backwards compatible, which means that the
2219 * numbering may have to be altered from old numbering to new numbering
2221 for (i = 0; (i < ffparams->ntypes); i++)
2225 /* Loop over file versions */
2226 for (k = 0; (k < NFTUPD); k++)
2228 /* Compare the read file_version to the update table */
2229 if ((file_version < ftupd[k].fvnr) &&
2230 (ffparams->functype[i] >= ftupd[k].ftype))
2232 ffparams->functype[i] += 1;
2235 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2236 i, ffparams->functype[i],
2237 interaction_function[ftupd[k].ftype].longname);
2244 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2248 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2253 static void add_settle_atoms(t_ilist *ilist)
2257 /* Settle used to only store the first atom: add the other two */
2258 srenew(ilist->iatoms, 2*ilist->nr);
2259 for (i = ilist->nr/2-1; i >= 0; i--)
2261 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2262 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2263 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2264 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2266 ilist->nr = 2*ilist->nr;
2269 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2272 int i, j, renum[F_NRE];
2276 for (j = 0; (j < F_NRE); j++)
2281 for (k = 0; k < NFTUPD; k++)
2283 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2292 ilist[j].iatoms = NULL;
2296 do_ilist(fio, &ilist[j], bRead, file_version, j);
2297 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2299 add_settle_atoms(&ilist[j]);
2303 if (bRead && gmx_debug_at)
2304 pr_ilist(debug,0,interaction_function[j].longname,
2305 functype,&ilist[j],TRUE);
2310 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2311 gmx_bool bRead, int file_version)
2313 do_ffparams(fio, ffparams, bRead, file_version);
2315 if (file_version >= 54)
2317 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2320 do_ilists(fio, molt->ilist, bRead, file_version);
2323 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2325 int i, idum, dum_nra, *dum_a;
2327 if (file_version < 44)
2329 for (i = 0; i < MAXNODES; i++)
2331 gmx_fio_do_int(fio, idum);
2334 gmx_fio_do_int(fio, block->nr);
2335 if (file_version < 51)
2337 gmx_fio_do_int(fio, dum_nra);
2341 if ((block->nalloc_index > 0) && (NULL != block->index))
2343 sfree(block->index);
2345 block->nalloc_index = block->nr+1;
2346 snew(block->index, block->nalloc_index);
2348 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2350 if (file_version < 51 && dum_nra > 0)
2352 snew(dum_a, dum_nra);
2353 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2358 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2363 if (file_version < 44)
2365 for (i = 0; i < MAXNODES; i++)
2367 gmx_fio_do_int(fio, idum);
2370 gmx_fio_do_int(fio, block->nr);
2371 gmx_fio_do_int(fio, block->nra);
2374 block->nalloc_index = block->nr+1;
2375 snew(block->index, block->nalloc_index);
2376 block->nalloc_a = block->nra;
2377 snew(block->a, block->nalloc_a);
2379 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2380 gmx_fio_ndo_int(fio, block->a, block->nra);
2383 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2384 int file_version, gmx_groups_t *groups, int atnr)
2388 gmx_fio_do_real(fio, atom->m);
2389 gmx_fio_do_real(fio, atom->q);
2390 gmx_fio_do_real(fio, atom->mB);
2391 gmx_fio_do_real(fio, atom->qB);
2392 gmx_fio_do_ushort(fio, atom->type);
2393 gmx_fio_do_ushort(fio, atom->typeB);
2394 gmx_fio_do_int(fio, atom->ptype);
2395 gmx_fio_do_int(fio, atom->resind);
2396 if (file_version >= 52)
2398 gmx_fio_do_int(fio, atom->atomnumber);
2402 atom->atomnumber = NOTSET;
2404 if (file_version < 23)
2408 else if (file_version < 39)
2417 if (file_version < 57)
2419 unsigned char uchar[egcNR];
2420 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2421 for (i = myngrp; (i < ngrp); i++)
2425 /* Copy the old data format to the groups struct */
2426 for (i = 0; i < ngrp; i++)
2428 groups->grpnr[i][atnr] = uchar[i];
2433 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2438 if (file_version < 23)
2442 else if (file_version < 39)
2451 for (j = 0; (j < ngrp); j++)
2455 gmx_fio_do_int(fio, grps[j].nr);
2458 snew(grps[j].nm_ind, grps[j].nr);
2460 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2465 snew(grps[j].nm_ind, grps[j].nr);
2470 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2476 gmx_fio_do_int(fio, ls);
2477 *nm = get_symtab_handle(symtab, ls);
2481 ls = lookup_symtab(symtab, *nm);
2482 gmx_fio_do_int(fio, ls);
2486 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2491 for (j = 0; (j < nstr); j++)
2493 do_symstr(fio, &(nm[j]), bRead, symtab);
2497 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2498 t_symtab *symtab, int file_version)
2502 for (j = 0; (j < n); j++)
2504 do_symstr(fio, &(ri[j].name), bRead, symtab);
2505 if (file_version >= 63)
2507 gmx_fio_do_int(fio, ri[j].nr);
2508 gmx_fio_do_uchar(fio, ri[j].ic);
2518 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2520 gmx_groups_t *groups)
2524 gmx_fio_do_int(fio, atoms->nr);
2525 gmx_fio_do_int(fio, atoms->nres);
2526 if (file_version < 57)
2528 gmx_fio_do_int(fio, groups->ngrpname);
2529 for (i = 0; i < egcNR; i++)
2531 groups->ngrpnr[i] = atoms->nr;
2532 snew(groups->grpnr[i], groups->ngrpnr[i]);
2537 snew(atoms->atom, atoms->nr);
2538 snew(atoms->atomname, atoms->nr);
2539 snew(atoms->atomtype, atoms->nr);
2540 snew(atoms->atomtypeB, atoms->nr);
2541 snew(atoms->resinfo, atoms->nres);
2542 if (file_version < 57)
2544 snew(groups->grpname, groups->ngrpname);
2546 atoms->pdbinfo = NULL;
2548 for (i = 0; (i < atoms->nr); i++)
2550 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2552 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2553 if (bRead && (file_version <= 20))
2555 for (i = 0; i < atoms->nr; i++)
2557 atoms->atomtype[i] = put_symtab(symtab, "?");
2558 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2563 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2564 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2566 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2568 if (file_version < 57)
2570 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2572 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2576 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2577 gmx_bool bRead, t_symtab *symtab,
2582 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2583 gmx_fio_do_int(fio, groups->ngrpname);
2586 snew(groups->grpname, groups->ngrpname);
2588 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2589 for (g = 0; g < egcNR; g++)
2591 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2592 if (groups->ngrpnr[g] == 0)
2596 groups->grpnr[g] = NULL;
2603 snew(groups->grpnr[g], groups->ngrpnr[g]);
2605 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2610 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2615 if (file_version > 25)
2617 gmx_fio_do_int(fio, atomtypes->nr);
2621 snew(atomtypes->radius, j);
2622 snew(atomtypes->vol, j);
2623 snew(atomtypes->surftens, j);
2624 snew(atomtypes->atomnumber, j);
2625 snew(atomtypes->gb_radius, j);
2626 snew(atomtypes->S_hct, j);
2628 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2629 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2630 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2631 if (file_version >= 40)
2633 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2635 if (file_version >= 60)
2637 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2638 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2643 /* File versions prior to 26 cannot do GBSA,
2644 * so they dont use this structure
2647 atomtypes->radius = NULL;
2648 atomtypes->vol = NULL;
2649 atomtypes->surftens = NULL;
2650 atomtypes->atomnumber = NULL;
2651 atomtypes->gb_radius = NULL;
2652 atomtypes->S_hct = NULL;
2656 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2662 gmx_fio_do_int(fio, symtab->nr);
2666 snew(symtab->symbuf, 1);
2667 symbuf = symtab->symbuf;
2668 symbuf->bufsize = nr;
2669 snew(symbuf->buf, nr);
2670 for (i = 0; (i < nr); i++)
2672 gmx_fio_do_string(fio, buf);
2673 symbuf->buf[i] = strdup(buf);
2678 symbuf = symtab->symbuf;
2679 while (symbuf != NULL)
2681 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2683 gmx_fio_do_string(fio, symbuf->buf[i]);
2686 symbuf = symbuf->next;
2690 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2695 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2697 int i, j, ngrid, gs, nelem;
2699 gmx_fio_do_int(fio, cmap_grid->ngrid);
2700 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2702 ngrid = cmap_grid->ngrid;
2703 gs = cmap_grid->grid_spacing;
2708 snew(cmap_grid->cmapdata, ngrid);
2710 for (i = 0; i < cmap_grid->ngrid; i++)
2712 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2716 for (i = 0; i < cmap_grid->ngrid; i++)
2718 for (j = 0; j < nelem; j++)
2720 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2721 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2722 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2723 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2729 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2731 int m, a, a0, a1, r;
2735 /* We always assign a new chain number, but save the chain id characters
2736 * for larger molecules.
2738 #define CHAIN_MIN_ATOMS 15
2742 for (m = 0; m < mols->nr; m++)
2744 a0 = mols->index[m];
2745 a1 = mols->index[m+1];
2746 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2755 for (a = a0; a < a1; a++)
2757 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2758 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2763 /* Blank out the chain id if there was only one chain */
2766 for (r = 0; r < atoms->nres; r++)
2768 atoms->resinfo[r].chainid = ' ';
2773 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2774 t_symtab *symtab, int file_version,
2775 gmx_groups_t *groups)
2779 if (file_version >= 57)
2781 do_symstr(fio, &(molt->name), bRead, symtab);
2784 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2786 if (bRead && gmx_debug_at)
2788 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2791 if (file_version >= 57)
2793 do_ilists(fio, molt->ilist, bRead, file_version);
2795 do_block(fio, &molt->cgs, bRead, file_version);
2796 if (bRead && gmx_debug_at)
2798 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2802 /* This used to be in the atoms struct */
2803 do_blocka(fio, &molt->excls, bRead, file_version);
2806 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2810 gmx_fio_do_int(fio, molb->type);
2811 gmx_fio_do_int(fio, molb->nmol);
2812 gmx_fio_do_int(fio, molb->natoms_mol);
2813 /* Position restraint coordinates */
2814 gmx_fio_do_int(fio, molb->nposres_xA);
2815 if (molb->nposres_xA > 0)
2819 snew(molb->posres_xA, molb->nposres_xA);
2821 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2823 gmx_fio_do_int(fio, molb->nposres_xB);
2824 if (molb->nposres_xB > 0)
2828 snew(molb->posres_xB, molb->nposres_xB);
2830 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2835 static t_block mtop_mols(gmx_mtop_t *mtop)
2841 for (mb = 0; mb < mtop->nmolblock; mb++)
2843 mols.nr += mtop->molblock[mb].nmol;
2845 mols.nalloc_index = mols.nr + 1;
2846 snew(mols.index, mols.nalloc_index);
2851 for (mb = 0; mb < mtop->nmolblock; mb++)
2853 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2855 a += mtop->molblock[mb].natoms_mol;
2864 static void add_posres_molblock(gmx_mtop_t *mtop)
2869 gmx_molblock_t *molb;
2872 /* posres reference positions are stored in ip->posres (if present) and
2873 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2874 posres.pos0A are identical to fbposres.pos0. */
2875 il = &mtop->moltype[0].ilist[F_POSRES];
2876 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2877 if (il->nr == 0 && ilfb->nr == 0)
2883 for (i = 0; i < il->nr; i += 2)
2885 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2886 am = max(am, il->iatoms[i+1]);
2887 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2888 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2889 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2894 /* This loop is required if we have only flat-bottomed posres:
2896 - bFE == FALSE (no B-state for flat-bottomed posres) */
2899 for (i = 0; i < ilfb->nr; i += 2)
2901 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2902 am = max(am, ilfb->iatoms[i+1]);
2905 /* Make the posres coordinate block end at a molecule end */
2907 while (am >= mtop->mols.index[mol+1])
2911 molb = &mtop->molblock[0];
2912 molb->nposres_xA = mtop->mols.index[mol+1];
2913 snew(molb->posres_xA, molb->nposres_xA);
2916 molb->nposres_xB = molb->nposres_xA;
2917 snew(molb->posres_xB, molb->nposres_xB);
2921 molb->nposres_xB = 0;
2923 for (i = 0; i < il->nr; i += 2)
2925 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2926 a = il->iatoms[i+1];
2927 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2928 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2929 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2932 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2933 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2934 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2939 /* If only flat-bottomed posres are present, take reference pos from them.
2940 Here: bFE == FALSE */
2941 for (i = 0; i < ilfb->nr; i += 2)
2943 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2944 a = ilfb->iatoms[i+1];
2945 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2946 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2947 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2952 static void set_disres_npair(gmx_mtop_t *mtop)
2959 ip = mtop->ffparams.iparams;
2961 for (mt = 0; mt < mtop->nmoltype; mt++)
2963 il = &mtop->moltype[mt].ilist[F_DISRES];
2968 for (i = 0; i < il->nr; i += 3)
2971 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2973 ip[a[i]].disres.npair = npair;
2981 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2991 do_symtab(fio, &(mtop->symtab), bRead);
2994 pr_symtab(debug, 0, "symtab", &mtop->symtab);
2997 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2999 if (file_version >= 57)
3001 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
3003 gmx_fio_do_int(fio, mtop->nmoltype);
3011 snew(mtop->moltype, mtop->nmoltype);
3012 if (file_version < 57)
3014 mtop->moltype[0].name = mtop->name;
3017 for (mt = 0; mt < mtop->nmoltype; mt++)
3019 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
3023 if (file_version >= 57)
3025 gmx_fio_do_int(fio, mtop->nmolblock);
3029 mtop->nmolblock = 1;
3033 snew(mtop->molblock, mtop->nmolblock);
3035 if (file_version >= 57)
3037 for (mb = 0; mb < mtop->nmolblock; mb++)
3039 do_molblock(fio, &mtop->molblock[mb], bRead);
3041 gmx_fio_do_int(fio, mtop->natoms);
3045 mtop->molblock[0].type = 0;
3046 mtop->molblock[0].nmol = 1;
3047 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
3048 mtop->molblock[0].nposres_xA = 0;
3049 mtop->molblock[0].nposres_xB = 0;
3052 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
3055 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
3058 if (file_version < 57)
3060 /* Debug statements are inside do_idef */
3061 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
3062 mtop->natoms = mtop->moltype[0].atoms.nr;
3065 if (file_version >= 65)
3067 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
3071 mtop->ffparams.cmap_grid.ngrid = 0;
3072 mtop->ffparams.cmap_grid.grid_spacing = 0;
3073 mtop->ffparams.cmap_grid.cmapdata = NULL;
3076 if (file_version >= 57)
3078 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
3081 if (file_version < 57)
3083 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
3084 if (bRead && gmx_debug_at)
3086 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
3088 do_block(fio, &mtop->mols, bRead, file_version);
3089 /* Add the posres coordinates to the molblock */
3090 add_posres_molblock(mtop);
3094 if (file_version >= 57)
3096 done_block(&mtop->mols);
3097 mtop->mols = mtop_mols(mtop);
3101 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
3105 if (file_version < 51)
3107 /* Here used to be the shake blocks */
3108 do_blocka(fio, &dumb, bRead, file_version);
3121 close_symtab(&(mtop->symtab));
3125 /* If TopOnlyOK is TRUE then we can read even future versions
3126 * of tpx files, provided the file_generation hasn't changed.
3127 * If it is FALSE, we need the inputrecord too, and bail out
3128 * if the file is newer than the program.
3130 * The version and generation if the topology (see top of this file)
3131 * are returned in the two last arguments.
3133 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
3135 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
3136 gmx_bool TopOnlyOK, int *file_version,
3137 int *file_generation)
3140 char file_tag[STRLEN];
3147 gmx_fio_checktype(fio);
3148 gmx_fio_setdebug(fio, bDebugMode());
3150 /* NEW! XDR tpb file */
3151 precision = sizeof(real);
3154 gmx_fio_do_string(fio, buf);
3155 if (strncmp(buf, "VERSION", 7))
3157 gmx_fatal(FARGS, "Can not read file %s,\n"
3158 " this file is from a Gromacs version which is older than 2.0\n"
3159 " Make a new one with grompp or use a gro or pdb file, if possible",
3160 gmx_fio_getname(fio));
3162 gmx_fio_do_int(fio, precision);
3163 bDouble = (precision == sizeof(double));
3164 if ((precision != sizeof(float)) && !bDouble)
3166 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3167 "instead of %d or %d",
3168 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3170 gmx_fio_setprecision(fio, bDouble);
3171 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3172 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3176 gmx_fio_write_string(fio, GromacsVersion());
3177 bDouble = (precision == sizeof(double));
3178 gmx_fio_setprecision(fio, bDouble);
3179 gmx_fio_do_int(fio, precision);
3181 sprintf(file_tag, "%s", tpx_tag);
3182 fgen = tpx_generation;
3185 /* Check versions! */
3186 gmx_fio_do_int(fio, fver);
3188 /* This is for backward compatibility with development versions 77-79
3189 * where the tag was, mistakenly, placed before the generation,
3190 * which would cause a segv instead of a proper error message
3191 * when reading the topology only from tpx with <77 code.
3193 if (fver >= 77 && fver <= 79)
3195 gmx_fio_do_string(fio, file_tag);
3200 gmx_fio_do_int(fio, fgen);
3209 gmx_fio_do_string(fio, file_tag);
3215 /* Versions before 77 don't have the tag, set it to release */
3216 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3219 if (strcmp(file_tag, tpx_tag) != 0)
3221 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3224 /* We only support reading tpx files with the same tag as the code
3225 * or tpx files with the release tag and with lower version number.
3227 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
3229 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3230 gmx_fio_getname(fio), fver, file_tag,
3231 tpx_version, tpx_tag);
3236 if (file_version != NULL)
3238 *file_version = fver;
3240 if (file_generation != NULL)
3242 *file_generation = fgen;
3246 if ((fver <= tpx_incompatible_version) ||
3247 ((fver > tpx_version) && !TopOnlyOK) ||
3248 (fgen > tpx_generation) ||
3249 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3251 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3252 gmx_fio_getname(fio), fver, tpx_version);
3255 do_section(fio, eitemHEADER, bRead);
3256 gmx_fio_do_int(fio, tpx->natoms);
3259 gmx_fio_do_int(fio, tpx->ngtc);
3267 gmx_fio_do_int(fio, idum);
3268 gmx_fio_do_real(fio, rdum);
3270 /*a better decision will eventually (5.0 or later) need to be made
3271 on how to treat the alchemical state of the system, which can now
3272 vary through a simulation, and cannot be completely described
3273 though a single lambda variable, or even a single state
3274 index. Eventually, should probably be a vector. MRS*/
3277 gmx_fio_do_int(fio, tpx->fep_state);
3279 gmx_fio_do_real(fio, tpx->lambda);
3280 gmx_fio_do_int(fio, tpx->bIr);
3281 gmx_fio_do_int(fio, tpx->bTop);
3282 gmx_fio_do_int(fio, tpx->bX);
3283 gmx_fio_do_int(fio, tpx->bV);
3284 gmx_fio_do_int(fio, tpx->bF);
3285 gmx_fio_do_int(fio, tpx->bBox);
3287 if ((fgen > tpx_generation))
3289 /* This can only happen if TopOnlyOK=TRUE */
3294 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3295 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3296 gmx_bool bXVallocated)
3302 int file_version, file_generation;
3306 gmx_bool bPeriodicMols;
3310 tpx.natoms = state->natoms;
3311 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3312 tpx.fep_state = state->fep_state;
3313 tpx.lambda = state->lambda[efptFEP];
3314 tpx.bIr = (ir != NULL);
3315 tpx.bTop = (mtop != NULL);
3316 tpx.bX = (state->x != NULL);
3317 tpx.bV = (state->v != NULL);
3318 tpx.bF = (f != NULL);
3322 TopOnlyOK = (ir == NULL);
3324 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3329 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3330 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3335 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3336 state->natoms = tpx.natoms;
3337 state->nalloc = tpx.natoms;
3343 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3347 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3349 do_test(fio, tpx.bBox, state->box);
3350 do_section(fio, eitemBOX, bRead);
3353 gmx_fio_ndo_rvec(fio, state->box, DIM);
3354 if (file_version >= 51)
3356 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3360 /* We initialize box_rel after reading the inputrec */
3361 clear_mat(state->box_rel);
3363 if (file_version >= 28)
3365 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3366 if (file_version < 56)
3369 gmx_fio_ndo_rvec(fio, mdum, DIM);
3374 if (state->ngtc > 0 && file_version >= 28)
3377 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3378 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3379 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3380 snew(dumv, state->ngtc);
3381 if (file_version < 69)
3383 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3385 /* These used to be the Berendsen tcoupl_lambda's */
3386 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3390 /* Prior to tpx version 26, the inputrec was here.
3391 * I moved it to enable partial forward-compatibility
3392 * for analysis/viewer programs.
3394 if (file_version < 26)
3396 do_test(fio, tpx.bIr, ir);
3397 do_section(fio, eitemIR, bRead);
3402 do_inputrec(fio, ir, bRead, file_version,
3403 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3406 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3411 do_inputrec(fio, &dum_ir, bRead, file_version,
3412 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3415 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3417 done_inputrec(&dum_ir);
3423 do_test(fio, tpx.bTop, mtop);
3424 do_section(fio, eitemTOP, bRead);
3429 do_mtop(fio, mtop, bRead, file_version);
3433 do_mtop(fio, &dum_top, bRead, file_version);
3434 done_mtop(&dum_top, TRUE);
3437 do_test(fio, tpx.bX, state->x);
3438 do_section(fio, eitemX, bRead);
3443 state->flags |= (1<<estX);
3445 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3448 do_test(fio, tpx.bV, state->v);
3449 do_section(fio, eitemV, bRead);
3454 state->flags |= (1<<estV);
3456 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3459 do_test(fio, tpx.bF, f);
3460 do_section(fio, eitemF, bRead);
3463 gmx_fio_ndo_rvec(fio, f, state->natoms);
3466 /* Starting with tpx version 26, we have the inputrec
3467 * at the end of the file, so we can ignore it
3468 * if the file is never than the software (but still the
3469 * same generation - see comments at the top of this file.
3474 bPeriodicMols = FALSE;
3475 if (file_version >= 26)
3477 do_test(fio, tpx.bIr, ir);
3478 do_section(fio, eitemIR, bRead);
3481 if (file_version >= 53)
3483 /* Removed the pbc info from do_inputrec, since we always want it */
3487 bPeriodicMols = ir->bPeriodicMols;
3489 gmx_fio_do_int(fio, ePBC);
3490 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3492 if (file_generation <= tpx_generation && ir)
3494 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3497 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3499 if (file_version < 51)
3501 set_box_rel(ir, state);
3503 if (file_version < 53)
3506 bPeriodicMols = ir->bPeriodicMols;
3509 if (bRead && ir && file_version >= 53)
3511 /* We need to do this after do_inputrec, since that initializes ir */
3513 ir->bPeriodicMols = bPeriodicMols;
3522 if (state->ngtc == 0)
3524 /* Reading old version without tcoupl state data: set it */
3525 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3527 if (tpx.bTop && mtop)
3529 if (file_version < 57)
3531 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3533 ir->eDisre = edrSimple;
3537 ir->eDisre = edrNone;
3540 set_disres_npair(mtop);
3544 if (tpx.bTop && mtop)
3546 gmx_mtop_finalize(mtop);
3549 if (file_version >= 57)
3553 env = getenv("GMX_NOCHARGEGROUPS");
3556 sscanf(env, "%d", &ienv);
3557 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3562 "Will make single atomic charge groups in non-solvent%s\n",
3563 ienv > 1 ? " and solvent" : "");
3564 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3566 fprintf(stderr, "\n");
3574 /************************************************************
3576 * The following routines are the exported ones
3578 ************************************************************/
3580 t_fileio *open_tpx(const char *fn, const char *mode)
3582 return gmx_fio_open(fn, mode);
3585 void close_tpx(t_fileio *fio)
3590 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3591 int *file_version, int *file_generation)
3595 fio = open_tpx(fn, "r");
3596 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3600 void write_tpx_state(const char *fn,
3601 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3605 fio = open_tpx(fn, "w");
3606 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3610 void read_tpx_state(const char *fn,
3611 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3615 fio = open_tpx(fn, "r");
3616 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3620 int read_tpx(const char *fn,
3621 t_inputrec *ir, matrix box, int *natoms,
3622 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3630 fio = open_tpx(fn, "r");
3631 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3633 *natoms = state.natoms;
3636 copy_mat(state.box, box);
3645 int read_tpx_top(const char *fn,
3646 t_inputrec *ir, matrix box, int *natoms,
3647 rvec *x, rvec *v, rvec *f, t_topology *top)
3653 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3655 *top = gmx_mtop_t_to_t_topology(&mtop);
3660 gmx_bool fn2bTPX(const char *file)
3662 switch (fn2ftp(file))
3673 static void done_gmx_groups_t(gmx_groups_t *g)
3677 for (i = 0; (i < egcNR); i++)
3679 if (NULL != g->grps[i].nm_ind)
3681 sfree(g->grps[i].nm_ind);
3682 g->grps[i].nm_ind = NULL;
3684 if (NULL != g->grpnr[i])
3690 /* The contents of this array is in symtab, don't free it here */
3694 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3695 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3698 int natoms, i, version, generation;
3699 gmx_bool bTop, bXNULL = FALSE;
3701 t_topology *topconv;
3704 bTop = fn2bTPX(infile);
3708 read_tpxheader(infile, &header, TRUE, &version, &generation);
3711 snew(*x, header.natoms);
3715 snew(*v, header.natoms);
3718 *ePBC = read_tpx(infile, NULL, box, &natoms,
3719 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3720 *top = gmx_mtop_t_to_t_topology(mtop);
3721 /* In this case we need to throw away the group data too */
3722 done_gmx_groups_t(&mtop->groups);
3724 strcpy(title, *top->name);
3725 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3729 get_stx_coordnum(infile, &natoms);
3730 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3741 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3749 aps = gmx_atomprop_init();
3750 for (i = 0; (i < natoms); i++)
3752 if (!gmx_atomprop_query(aps, epropMass,
3753 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3754 *top->atoms.atomname[i],
3755 &(top->atoms.atom[i].m)))
3759 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3760 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3761 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3762 *top->atoms.atomname[i]);
3766 gmx_atomprop_destroy(aps);
3768 top->idef.ntypes = -1;