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.
39 /* This file is completely threadsafe - keep it that way! */
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/legacyheaders/names.h"
46 #include "gromacs/utility/futil.h"
50 #include "gromacs/legacyheaders/txtdump.h"
52 #include "gromacs/legacyheaders/copyrite.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/topology/atomprop.h"
56 #include "gromacs/topology/block.h"
57 #include "gromacs/topology/mtop_util.h"
58 #include "gromacs/topology/symtab.h"
59 #include "gromacs/topology/topology.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/smalloc.h"
64 #define TPX_TAG_RELEASE "release"
66 /*! \brief Tag string for the file format written to run input files
67 * written by this version of the code.
69 * Change this if you want to change the run input format in a feature
70 * branch. This ensures that there will not be different run input
71 * formats around which cannot be distinguished, while not causing
72 * problems rebasing the feature branch onto upstream changes. When
73 * merging with mainstream GROMACS, set this tag string back to
74 * TPX_TAG_RELEASE, and instead add an element to tpxv and set
75 * tpx_version to that.
77 static const char *tpx_tag = TPX_TAG_RELEASE;
79 /*! \brief Enum of values that describe the contents of a tpr file
80 * whose format matches a version number
82 * The enum helps the code be more self-documenting and ensure merges
83 * do not silently resolve when two patches make the same bump. When
84 * adding new functionality, add a new element to the end of this
85 * enumeration, change the definition of tpx_version, and write code
86 * below that does the right thing according to the value of
89 tpxv_ComputationalElectrophysiology = 96, /**< support for ion/water position swaps (computational electrophysiology) */
90 tpxv_Use64BitRandomSeed, /**< change ld_seed from int to gmx_int64_t */
91 tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, /**< potentials for supporting coarse-grained force fields */
92 tpxv_InteractiveMolecularDynamics, /**< interactive molecular dynamics (IMD) */
93 tpxv_RemoveObsoleteParameters1 /**< remove optimize_fft, dihre_fc, nstcheckpoint */
96 /*! \brief Version number of the file format written to run input
97 * files by this version of the code.
99 * The tpx_version number should be increased whenever the file format
100 * in the main development branch changes, generally to the highest
101 * value present in tpxv. Backward compatibility for reading old run
102 * input files is maintained by checking this version number against
103 * that of the file and then using the correct code path.
105 * When developing a feature branch that needs to change the run input
106 * file format, change tpx_tag instead. */
107 static const int tpx_version = tpxv_RemoveObsoleteParameters1;
110 /* This number should only be increased when you edit the TOPOLOGY section
111 * or the HEADER of the tpx format.
112 * This way we can maintain forward compatibility too for all analysis tools
113 * and/or external programs that only need to know the atom/residue names,
114 * charges, and bond connectivity.
116 * It first appeared in tpx version 26, when I also moved the inputrecord
117 * to the end of the tpx file, so we can just skip it if we only
120 * In particular, it must be increased when adding new elements to
121 * ftupd, so that old code can read new .tpr files.
123 static const int tpx_generation = 26;
125 /* This number should be the most recent backwards incompatible version
126 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
128 static const int tpx_incompatible_version = 9;
132 /* Struct used to maintain tpx compatibility when function types are added */
134 int fvnr; /* file version number in which the function type first appeared */
135 int ftype; /* function type */
139 * The entries should be ordered in:
140 * 1. ascending file version number
141 * 2. ascending function type number
143 /*static const t_ftupd ftupd[] = {
144 { 20, F_CUBICBONDS },
148 { 22, F_DISRESVIOL },
154 { 26, F_DIHRESVIOL },
155 { 30, F_CROSS_BOND_BONDS },
156 { 30, F_CROSS_BOND_ANGLES },
157 { 30, F_UREY_BRADLEY },
158 { 30, F_POLARIZATION },
162 * The entries should be ordered in:
163 * 1. ascending function type number
164 * 2. ascending file version number
166 /* question; what is the purpose of the commented code above? */
167 static const t_ftupd ftupd[] = {
168 { 20, F_CUBICBONDS },
173 { 43, F_TABBONDSNC },
174 { 70, F_RESTRBONDS },
175 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRANGLES },
176 { 76, F_LINEAR_ANGLES },
177 { 30, F_CROSS_BOND_BONDS },
178 { 30, F_CROSS_BOND_ANGLES },
179 { 30, F_UREY_BRADLEY },
180 { 34, F_QUARTIC_ANGLES },
182 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_RESTRDIHS },
183 { tpxv_RestrictedBendingAndCombinedAngleTorsionPotentials, F_CBTDIHS },
192 { 72, F_NPSOLVATION },
194 { 41, F_LJC_PAIRS_NB },
197 { 32, F_COUL_RECIP },
200 { 30, F_POLARIZATION },
203 { 22, F_DISRESVIOL },
207 { 26, F_DIHRESVIOL },
212 { 46, F_ECONSERVED },
213 { 69, F_VTEMP_NOLONGERUSED},
215 { 54, F_DVDL_CONSTR },
216 { 76, F_ANHARM_POL },
219 { 79, F_DVDL_BONDED, },
220 { 79, F_DVDL_RESTRAINT },
221 { 79, F_DVDL_TEMPERATURE },
223 #define NFTUPD asize(ftupd)
225 /* Needed for backward compatibility */
228 /**************************************************************
230 * Now the higer level routines that do io of the structures and arrays
232 **************************************************************/
233 static void do_pullgrp_tpx_pre95(t_fileio *fio,
242 gmx_fio_do_int(fio, pgrp->nat);
245 snew(pgrp->ind, pgrp->nat);
247 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
248 gmx_fio_do_int(fio, pgrp->nweight);
251 snew(pgrp->weight, pgrp->nweight);
253 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
254 gmx_fio_do_int(fio, pgrp->pbcatom);
255 gmx_fio_do_rvec(fio, pcrd->vec);
256 clear_rvec(pcrd->origin);
257 gmx_fio_do_rvec(fio, tmp);
259 gmx_fio_do_real(fio, pcrd->rate);
260 gmx_fio_do_real(fio, pcrd->k);
261 if (file_version >= 56)
263 gmx_fio_do_real(fio, pcrd->kB);
271 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
275 gmx_fio_do_int(fio, pgrp->nat);
278 snew(pgrp->ind, pgrp->nat);
280 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
281 gmx_fio_do_int(fio, pgrp->nweight);
284 snew(pgrp->weight, pgrp->nweight);
286 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
287 gmx_fio_do_int(fio, pgrp->pbcatom);
290 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd)
294 gmx_fio_do_int(fio, pcrd->group[0]);
295 gmx_fio_do_int(fio, pcrd->group[1]);
296 gmx_fio_do_rvec(fio, pcrd->origin);
297 gmx_fio_do_rvec(fio, pcrd->vec);
298 gmx_fio_do_real(fio, pcrd->init);
299 gmx_fio_do_real(fio, pcrd->rate);
300 gmx_fio_do_real(fio, pcrd->k);
301 gmx_fio_do_real(fio, pcrd->kB);
304 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
306 /* i is used in the ndo_double macro*/
310 int n_lambda = fepvals->n_lambda;
312 /* reset the lambda calculation window */
313 fepvals->lambda_start_n = 0;
314 fepvals->lambda_stop_n = n_lambda;
315 if (file_version >= 79)
321 snew(expand->init_lambda_weights, n_lambda);
323 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
324 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
327 gmx_fio_do_int(fio, expand->nstexpanded);
328 gmx_fio_do_int(fio, expand->elmcmove);
329 gmx_fio_do_int(fio, expand->elamstats);
330 gmx_fio_do_int(fio, expand->lmc_repeats);
331 gmx_fio_do_int(fio, expand->gibbsdeltalam);
332 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
333 gmx_fio_do_int(fio, expand->lmc_seed);
334 gmx_fio_do_real(fio, expand->mc_temp);
335 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
336 gmx_fio_do_int(fio, expand->nstTij);
337 gmx_fio_do_int(fio, expand->minvarmin);
338 gmx_fio_do_int(fio, expand->c_range);
339 gmx_fio_do_real(fio, expand->wl_scale);
340 gmx_fio_do_real(fio, expand->wl_ratio);
341 gmx_fio_do_real(fio, expand->init_wl_delta);
342 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
343 gmx_fio_do_int(fio, expand->elmceq);
344 gmx_fio_do_int(fio, expand->equil_steps);
345 gmx_fio_do_int(fio, expand->equil_samples);
346 gmx_fio_do_int(fio, expand->equil_n_at_lam);
347 gmx_fio_do_real(fio, expand->equil_wl_delta);
348 gmx_fio_do_real(fio, expand->equil_ratio);
352 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
355 if (file_version >= 79)
357 gmx_fio_do_int(fio, simtemp->eSimTempScale);
358 gmx_fio_do_real(fio, simtemp->simtemp_high);
359 gmx_fio_do_real(fio, simtemp->simtemp_low);
364 snew(simtemp->temperatures, n_lambda);
366 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
371 static void do_imd(t_fileio *fio, t_IMD *imd, gmx_bool bRead)
373 gmx_fio_do_int(fio, imd->nat);
376 snew(imd->ind, imd->nat);
378 gmx_fio_ndo_int(fio, imd->ind, imd->nat);
381 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
383 /* i is defined in the ndo_double macro; use g to iterate. */
388 /* free energy values */
390 if (file_version >= 79)
392 gmx_fio_do_int(fio, fepvals->init_fep_state);
393 gmx_fio_do_double(fio, fepvals->init_lambda);
394 gmx_fio_do_double(fio, fepvals->delta_lambda);
396 else if (file_version >= 59)
398 gmx_fio_do_double(fio, fepvals->init_lambda);
399 gmx_fio_do_double(fio, fepvals->delta_lambda);
403 gmx_fio_do_real(fio, rdum);
404 fepvals->init_lambda = rdum;
405 gmx_fio_do_real(fio, rdum);
406 fepvals->delta_lambda = rdum;
408 if (file_version >= 79)
410 gmx_fio_do_int(fio, fepvals->n_lambda);
413 snew(fepvals->all_lambda, efptNR);
415 for (g = 0; g < efptNR; g++)
417 if (fepvals->n_lambda > 0)
421 snew(fepvals->all_lambda[g], fepvals->n_lambda);
423 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
424 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
426 else if (fepvals->init_lambda >= 0)
428 fepvals->separate_dvdl[efptFEP] = TRUE;
432 else if (file_version >= 64)
434 gmx_fio_do_int(fio, fepvals->n_lambda);
439 snew(fepvals->all_lambda, efptNR);
440 /* still allocate the all_lambda array's contents. */
441 for (g = 0; g < efptNR; g++)
443 if (fepvals->n_lambda > 0)
445 snew(fepvals->all_lambda[g], fepvals->n_lambda);
449 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
451 if (fepvals->init_lambda >= 0)
455 fepvals->separate_dvdl[efptFEP] = TRUE;
459 /* copy the contents of the efptFEP lambda component to all
460 the other components */
461 for (g = 0; g < efptNR; g++)
463 for (h = 0; h < fepvals->n_lambda; h++)
467 fepvals->all_lambda[g][h] =
468 fepvals->all_lambda[efptFEP][h];
477 fepvals->n_lambda = 0;
478 fepvals->all_lambda = NULL;
479 if (fepvals->init_lambda >= 0)
481 fepvals->separate_dvdl[efptFEP] = TRUE;
484 if (file_version >= 13)
486 gmx_fio_do_real(fio, fepvals->sc_alpha);
490 fepvals->sc_alpha = 0;
492 if (file_version >= 38)
494 gmx_fio_do_int(fio, fepvals->sc_power);
498 fepvals->sc_power = 2;
500 if (file_version >= 79)
502 gmx_fio_do_real(fio, fepvals->sc_r_power);
506 fepvals->sc_r_power = 6.0;
508 if (file_version >= 15)
510 gmx_fio_do_real(fio, fepvals->sc_sigma);
514 fepvals->sc_sigma = 0.3;
518 if (file_version >= 71)
520 fepvals->sc_sigma_min = fepvals->sc_sigma;
524 fepvals->sc_sigma_min = 0;
527 if (file_version >= 79)
529 gmx_fio_do_int(fio, fepvals->bScCoul);
533 fepvals->bScCoul = TRUE;
535 if (file_version >= 64)
537 gmx_fio_do_int(fio, fepvals->nstdhdl);
541 fepvals->nstdhdl = 1;
544 if (file_version >= 73)
546 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
547 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
551 fepvals->separate_dhdl_file = esepdhdlfileYES;
552 fepvals->dhdl_derivatives = edhdlderivativesYES;
554 if (file_version >= 71)
556 gmx_fio_do_int(fio, fepvals->dh_hist_size);
557 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
561 fepvals->dh_hist_size = 0;
562 fepvals->dh_hist_spacing = 0.1;
564 if (file_version >= 79)
566 gmx_fio_do_int(fio, fepvals->edHdLPrintEnergy);
570 fepvals->edHdLPrintEnergy = edHdLPrintEnergyNO;
573 /* handle lambda_neighbors */
574 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
576 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
577 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
578 (fepvals->init_lambda < 0) )
580 fepvals->lambda_start_n = (fepvals->init_fep_state -
581 fepvals->lambda_neighbors);
582 fepvals->lambda_stop_n = (fepvals->init_fep_state +
583 fepvals->lambda_neighbors + 1);
584 if (fepvals->lambda_start_n < 0)
586 fepvals->lambda_start_n = 0;;
588 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
590 fepvals->lambda_stop_n = fepvals->n_lambda;
595 fepvals->lambda_start_n = 0;
596 fepvals->lambda_stop_n = fepvals->n_lambda;
601 fepvals->lambda_start_n = 0;
602 fepvals->lambda_stop_n = fepvals->n_lambda;
606 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
610 if (file_version >= 95)
612 gmx_fio_do_int(fio, pull->ngroup);
614 gmx_fio_do_int(fio, pull->ncoord);
615 if (file_version < 95)
617 pull->ngroup = pull->ncoord + 1;
619 gmx_fio_do_int(fio, pull->eGeom);
620 gmx_fio_do_ivec(fio, pull->dim);
621 gmx_fio_do_real(fio, pull->cyl_r1);
622 gmx_fio_do_real(fio, pull->cyl_r0);
623 gmx_fio_do_real(fio, pull->constr_tol);
624 if (file_version >= 95)
626 gmx_fio_do_int(fio, pull->bPrintRef);
628 gmx_fio_do_int(fio, pull->nstxout);
629 gmx_fio_do_int(fio, pull->nstfout);
632 snew(pull->group, pull->ngroup);
633 snew(pull->coord, pull->ncoord);
635 if (file_version < 95)
637 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
638 if (pull->eGeom == epullgDIRPBC)
640 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
642 if (pull->eGeom > epullgDIRPBC)
647 for (g = 0; g < pull->ngroup; g++)
649 /* We read and ignore a pull coordinate for group 0 */
650 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[max(g-1, 0)],
651 bRead, file_version);
654 pull->coord[g-1].group[0] = 0;
655 pull->coord[g-1].group[1] = g;
659 pull->bPrintRef = (pull->group[0].nat > 0);
663 for (g = 0; g < pull->ngroup; g++)
665 do_pull_group(fio, &pull->group[g], bRead);
667 for (g = 0; g < pull->ncoord; g++)
669 do_pull_coord(fio, &pull->coord[g]);
675 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
679 gmx_fio_do_int(fio, rotg->eType);
680 gmx_fio_do_int(fio, rotg->bMassW);
681 gmx_fio_do_int(fio, rotg->nat);
684 snew(rotg->ind, rotg->nat);
686 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
689 snew(rotg->x_ref, rotg->nat);
691 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
692 gmx_fio_do_rvec(fio, rotg->vec);
693 gmx_fio_do_rvec(fio, rotg->pivot);
694 gmx_fio_do_real(fio, rotg->rate);
695 gmx_fio_do_real(fio, rotg->k);
696 gmx_fio_do_real(fio, rotg->slab_dist);
697 gmx_fio_do_real(fio, rotg->min_gaussian);
698 gmx_fio_do_real(fio, rotg->eps);
699 gmx_fio_do_int(fio, rotg->eFittype);
700 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
701 gmx_fio_do_real(fio, rotg->PotAngle_step);
704 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
708 gmx_fio_do_int(fio, rot->ngrp);
709 gmx_fio_do_int(fio, rot->nstrout);
710 gmx_fio_do_int(fio, rot->nstsout);
713 snew(rot->grp, rot->ngrp);
715 for (g = 0; g < rot->ngrp; g++)
717 do_rotgrp(fio, &rot->grp[g], bRead);
722 static void do_swapcoords(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead)
727 gmx_fio_do_int(fio, swap->nat);
728 gmx_fio_do_int(fio, swap->nat_sol);
729 for (j = 0; j < 2; j++)
731 gmx_fio_do_int(fio, swap->nat_split[j]);
732 gmx_fio_do_int(fio, swap->massw_split[j]);
734 gmx_fio_do_int(fio, swap->nstswap);
735 gmx_fio_do_int(fio, swap->nAverage);
736 gmx_fio_do_real(fio, swap->threshold);
737 gmx_fio_do_real(fio, swap->cyl0r);
738 gmx_fio_do_real(fio, swap->cyl0u);
739 gmx_fio_do_real(fio, swap->cyl0l);
740 gmx_fio_do_real(fio, swap->cyl1r);
741 gmx_fio_do_real(fio, swap->cyl1u);
742 gmx_fio_do_real(fio, swap->cyl1l);
746 snew(swap->ind, swap->nat);
747 snew(swap->ind_sol, swap->nat_sol);
748 for (j = 0; j < 2; j++)
750 snew(swap->ind_split[j], swap->nat_split[j]);
754 gmx_fio_ndo_int(fio, swap->ind, swap->nat);
755 gmx_fio_ndo_int(fio, swap->ind_sol, swap->nat_sol);
756 for (j = 0; j < 2; j++)
758 gmx_fio_ndo_int(fio, swap->ind_split[j], swap->nat_split[j]);
761 for (j = 0; j < eCompNR; j++)
763 gmx_fio_do_int(fio, swap->nanions[j]);
764 gmx_fio_do_int(fio, swap->ncations[j]);
770 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
771 int file_version, real *fudgeQQ)
773 int i, j, k, *tmp, idum = 0;
776 gmx_bool bSimAnn, bdum = 0;
777 real zerotemptime, finish_t, init_temp, finish_temp;
779 if (file_version != tpx_version)
781 /* Give a warning about features that are not accessible */
782 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
783 file_version, tpx_version);
791 if (file_version == 0)
796 /* Basic inputrec stuff */
797 gmx_fio_do_int(fio, ir->eI);
798 if (file_version >= 62)
800 gmx_fio_do_int64(fio, ir->nsteps);
804 gmx_fio_do_int(fio, idum);
808 if (file_version > 25)
810 if (file_version >= 62)
812 gmx_fio_do_int64(fio, ir->init_step);
816 gmx_fio_do_int(fio, idum);
817 ir->init_step = idum;
825 if (file_version >= 58)
827 gmx_fio_do_int(fio, ir->simulation_part);
831 ir->simulation_part = 1;
834 if (file_version >= 67)
836 gmx_fio_do_int(fio, ir->nstcalcenergy);
840 ir->nstcalcenergy = 1;
842 if (file_version < 53)
844 /* The pbc info has been moved out of do_inputrec,
845 * since we always want it, also without reading the inputrec.
847 gmx_fio_do_int(fio, ir->ePBC);
848 if ((file_version <= 15) && (ir->ePBC == 2))
852 if (file_version >= 45)
854 gmx_fio_do_int(fio, ir->bPeriodicMols);
861 ir->bPeriodicMols = TRUE;
865 ir->bPeriodicMols = FALSE;
869 if (file_version >= 81)
871 gmx_fio_do_int(fio, ir->cutoff_scheme);
872 if (file_version < 94)
874 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
879 ir->cutoff_scheme = ecutsGROUP;
881 gmx_fio_do_int(fio, ir->ns_type);
882 gmx_fio_do_int(fio, ir->nstlist);
883 gmx_fio_do_int(fio, idum); /* used to be ndelta; not used anymore */
884 if (file_version < 41)
886 gmx_fio_do_int(fio, idum);
887 gmx_fio_do_int(fio, idum);
889 if (file_version >= 45)
891 gmx_fio_do_real(fio, ir->rtpi);
897 gmx_fio_do_int(fio, ir->nstcomm);
898 if (file_version > 34)
900 gmx_fio_do_int(fio, ir->comm_mode);
902 else if (ir->nstcomm < 0)
904 ir->comm_mode = ecmANGULAR;
908 ir->comm_mode = ecmLINEAR;
910 ir->nstcomm = abs(ir->nstcomm);
912 /* ignore nstcheckpoint */
913 if (file_version > 25 && file_version < tpxv_RemoveObsoleteParameters1)
915 gmx_fio_do_int(fio, idum);
918 gmx_fio_do_int(fio, ir->nstcgsteep);
920 if (file_version >= 30)
922 gmx_fio_do_int(fio, ir->nbfgscorr);
929 gmx_fio_do_int(fio, ir->nstlog);
930 gmx_fio_do_int(fio, ir->nstxout);
931 gmx_fio_do_int(fio, ir->nstvout);
932 gmx_fio_do_int(fio, ir->nstfout);
933 gmx_fio_do_int(fio, ir->nstenergy);
934 gmx_fio_do_int(fio, ir->nstxout_compressed);
935 if (file_version >= 59)
937 gmx_fio_do_double(fio, ir->init_t);
938 gmx_fio_do_double(fio, ir->delta_t);
942 gmx_fio_do_real(fio, rdum);
944 gmx_fio_do_real(fio, rdum);
947 gmx_fio_do_real(fio, ir->x_compression_precision);
948 if (file_version < 19)
950 gmx_fio_do_int(fio, idum);
951 gmx_fio_do_int(fio, idum);
953 if (file_version < 18)
955 gmx_fio_do_int(fio, idum);
957 if (file_version >= 81)
959 gmx_fio_do_real(fio, ir->verletbuf_tol);
963 ir->verletbuf_tol = 0;
965 gmx_fio_do_real(fio, ir->rlist);
966 if (file_version >= 67)
968 gmx_fio_do_real(fio, ir->rlistlong);
970 if (file_version >= 82 && file_version != 90)
972 gmx_fio_do_int(fio, ir->nstcalclr);
976 /* Calculate at NS steps */
977 ir->nstcalclr = ir->nstlist;
979 gmx_fio_do_int(fio, ir->coulombtype);
980 if (file_version < 32 && ir->coulombtype == eelRF)
982 ir->coulombtype = eelRF_NEC;
984 if (file_version >= 81)
986 gmx_fio_do_int(fio, ir->coulomb_modifier);
990 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
992 gmx_fio_do_real(fio, ir->rcoulomb_switch);
993 gmx_fio_do_real(fio, ir->rcoulomb);
994 gmx_fio_do_int(fio, ir->vdwtype);
995 if (file_version >= 81)
997 gmx_fio_do_int(fio, ir->vdw_modifier);
1001 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1003 gmx_fio_do_real(fio, ir->rvdw_switch);
1004 gmx_fio_do_real(fio, ir->rvdw);
1005 if (file_version < 67)
1007 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1009 gmx_fio_do_int(fio, ir->eDispCorr);
1010 gmx_fio_do_real(fio, ir->epsilon_r);
1011 if (file_version >= 37)
1013 gmx_fio_do_real(fio, ir->epsilon_rf);
1017 if (EEL_RF(ir->coulombtype))
1019 ir->epsilon_rf = ir->epsilon_r;
1020 ir->epsilon_r = 1.0;
1024 ir->epsilon_rf = 1.0;
1027 if (file_version >= 29)
1029 gmx_fio_do_real(fio, ir->tabext);
1036 if (file_version > 25)
1038 gmx_fio_do_int(fio, ir->gb_algorithm);
1039 gmx_fio_do_int(fio, ir->nstgbradii);
1040 gmx_fio_do_real(fio, ir->rgbradii);
1041 gmx_fio_do_real(fio, ir->gb_saltconc);
1042 gmx_fio_do_int(fio, ir->implicit_solvent);
1046 ir->gb_algorithm = egbSTILL;
1049 ir->gb_saltconc = 0;
1050 ir->implicit_solvent = eisNO;
1052 if (file_version >= 55)
1054 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1055 gmx_fio_do_real(fio, ir->gb_obc_alpha);
1056 gmx_fio_do_real(fio, ir->gb_obc_beta);
1057 gmx_fio_do_real(fio, ir->gb_obc_gamma);
1058 if (file_version >= 60)
1060 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1061 gmx_fio_do_int(fio, ir->sa_algorithm);
1065 ir->gb_dielectric_offset = 0.009;
1066 ir->sa_algorithm = esaAPPROX;
1068 gmx_fio_do_real(fio, ir->sa_surface_tension);
1070 /* Override sa_surface_tension if it is not changed in the mpd-file */
1071 if (ir->sa_surface_tension < 0)
1073 if (ir->gb_algorithm == egbSTILL)
1075 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1077 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1079 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1086 /* Better use sensible values than insane (0.0) ones... */
1087 ir->gb_epsilon_solvent = 80;
1088 ir->gb_obc_alpha = 1.0;
1089 ir->gb_obc_beta = 0.8;
1090 ir->gb_obc_gamma = 4.85;
1091 ir->sa_surface_tension = 2.092;
1095 if (file_version >= 81)
1097 gmx_fio_do_real(fio, ir->fourier_spacing);
1101 ir->fourier_spacing = 0.0;
1103 gmx_fio_do_int(fio, ir->nkx);
1104 gmx_fio_do_int(fio, ir->nky);
1105 gmx_fio_do_int(fio, ir->nkz);
1106 gmx_fio_do_int(fio, ir->pme_order);
1107 gmx_fio_do_real(fio, ir->ewald_rtol);
1109 if (file_version >= 93)
1111 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1115 ir->ewald_rtol_lj = ir->ewald_rtol;
1118 if (file_version >= 24)
1120 gmx_fio_do_int(fio, ir->ewald_geometry);
1124 ir->ewald_geometry = eewg3D;
1127 if (file_version <= 17)
1129 ir->epsilon_surface = 0;
1130 if (file_version == 17)
1132 gmx_fio_do_int(fio, idum);
1137 gmx_fio_do_real(fio, ir->epsilon_surface);
1140 /* ignore bOptFFT */
1141 if (file_version < tpxv_RemoveObsoleteParameters1)
1143 gmx_fio_do_gmx_bool(fio, bdum);
1146 if (file_version >= 93)
1148 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1150 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1151 gmx_fio_do_int(fio, ir->etc);
1152 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1153 * but the values 0 and 1 still mean no and
1154 * berendsen temperature coupling, respectively.
1156 if (file_version >= 79)
1158 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1160 if (file_version >= 71)
1162 gmx_fio_do_int(fio, ir->nsttcouple);
1166 ir->nsttcouple = ir->nstcalcenergy;
1168 if (file_version <= 15)
1170 gmx_fio_do_int(fio, idum);
1172 if (file_version <= 17)
1174 gmx_fio_do_int(fio, ir->epct);
1175 if (file_version <= 15)
1179 ir->epct = epctSURFACETENSION;
1181 gmx_fio_do_int(fio, idum);
1184 /* we have removed the NO alternative at the beginning */
1188 ir->epct = epctISOTROPIC;
1192 ir->epc = epcBERENDSEN;
1197 gmx_fio_do_int(fio, ir->epc);
1198 gmx_fio_do_int(fio, ir->epct);
1200 if (file_version >= 71)
1202 gmx_fio_do_int(fio, ir->nstpcouple);
1206 ir->nstpcouple = ir->nstcalcenergy;
1208 gmx_fio_do_real(fio, ir->tau_p);
1209 if (file_version <= 15)
1211 gmx_fio_do_rvec(fio, vdum);
1212 clear_mat(ir->ref_p);
1213 for (i = 0; i < DIM; i++)
1215 ir->ref_p[i][i] = vdum[i];
1220 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1221 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1222 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1224 if (file_version <= 15)
1226 gmx_fio_do_rvec(fio, vdum);
1227 clear_mat(ir->compress);
1228 for (i = 0; i < DIM; i++)
1230 ir->compress[i][i] = vdum[i];
1235 gmx_fio_do_rvec(fio, ir->compress[XX]);
1236 gmx_fio_do_rvec(fio, ir->compress[YY]);
1237 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1239 if (file_version >= 47)
1241 gmx_fio_do_int(fio, ir->refcoord_scaling);
1242 gmx_fio_do_rvec(fio, ir->posres_com);
1243 gmx_fio_do_rvec(fio, ir->posres_comB);
1247 ir->refcoord_scaling = erscNO;
1248 clear_rvec(ir->posres_com);
1249 clear_rvec(ir->posres_comB);
1251 if ((file_version > 25) && (file_version < 79))
1253 gmx_fio_do_int(fio, ir->andersen_seed);
1257 ir->andersen_seed = 0;
1259 if (file_version < 26)
1261 gmx_fio_do_gmx_bool(fio, bSimAnn);
1262 gmx_fio_do_real(fio, zerotemptime);
1265 if (file_version < 37)
1267 gmx_fio_do_real(fio, rdum);
1270 gmx_fio_do_real(fio, ir->shake_tol);
1271 if (file_version < 54)
1273 gmx_fio_do_real(fio, *fudgeQQ);
1276 gmx_fio_do_int(fio, ir->efep);
1277 if (file_version <= 14 && ir->efep != efepNO)
1281 do_fepvals(fio, ir->fepvals, bRead, file_version);
1283 if (file_version >= 79)
1285 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1288 ir->bSimTemp = TRUE;
1293 ir->bSimTemp = FALSE;
1297 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1300 if (file_version >= 79)
1302 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1305 ir->bExpanded = TRUE;
1309 ir->bExpanded = FALSE;
1314 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1316 if (file_version >= 57)
1318 gmx_fio_do_int(fio, ir->eDisre);
1320 gmx_fio_do_int(fio, ir->eDisreWeighting);
1321 if (file_version < 22)
1323 if (ir->eDisreWeighting == 0)
1325 ir->eDisreWeighting = edrwEqual;
1329 ir->eDisreWeighting = edrwConservative;
1332 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1333 gmx_fio_do_real(fio, ir->dr_fc);
1334 gmx_fio_do_real(fio, ir->dr_tau);
1335 gmx_fio_do_int(fio, ir->nstdisreout);
1336 if (file_version >= 22)
1338 gmx_fio_do_real(fio, ir->orires_fc);
1339 gmx_fio_do_real(fio, ir->orires_tau);
1340 gmx_fio_do_int(fio, ir->nstorireout);
1346 ir->nstorireout = 0;
1349 /* ignore dihre_fc */
1350 if (file_version >= 26 && file_version < 79)
1352 gmx_fio_do_real(fio, rdum);
1353 if (file_version < 56)
1355 gmx_fio_do_real(fio, rdum);
1356 gmx_fio_do_int(fio, idum);
1360 gmx_fio_do_real(fio, ir->em_stepsize);
1361 gmx_fio_do_real(fio, ir->em_tol);
1362 if (file_version >= 22)
1364 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1368 ir->bShakeSOR = TRUE;
1370 if (file_version >= 11)
1372 gmx_fio_do_int(fio, ir->niter);
1377 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1380 if (file_version >= 21)
1382 gmx_fio_do_real(fio, ir->fc_stepsize);
1386 ir->fc_stepsize = 0;
1388 gmx_fio_do_int(fio, ir->eConstrAlg);
1389 gmx_fio_do_int(fio, ir->nProjOrder);
1390 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1391 if (file_version <= 14)
1393 gmx_fio_do_int(fio, idum);
1395 if (file_version >= 26)
1397 gmx_fio_do_int(fio, ir->nLincsIter);
1402 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1405 if (file_version < 33)
1407 gmx_fio_do_real(fio, bd_temp);
1409 gmx_fio_do_real(fio, ir->bd_fric);
1410 if (file_version >= tpxv_Use64BitRandomSeed)
1412 gmx_fio_do_int64(fio, ir->ld_seed);
1416 gmx_fio_do_int(fio, idum);
1419 if (file_version >= 33)
1421 for (i = 0; i < DIM; i++)
1423 gmx_fio_do_rvec(fio, ir->deform[i]);
1428 for (i = 0; i < DIM; i++)
1430 clear_rvec(ir->deform[i]);
1433 if (file_version >= 14)
1435 gmx_fio_do_real(fio, ir->cos_accel);
1441 gmx_fio_do_int(fio, ir->userint1);
1442 gmx_fio_do_int(fio, ir->userint2);
1443 gmx_fio_do_int(fio, ir->userint3);
1444 gmx_fio_do_int(fio, ir->userint4);
1445 gmx_fio_do_real(fio, ir->userreal1);
1446 gmx_fio_do_real(fio, ir->userreal2);
1447 gmx_fio_do_real(fio, ir->userreal3);
1448 gmx_fio_do_real(fio, ir->userreal4);
1451 if (file_version >= 77)
1453 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1458 snew(ir->adress, 1);
1460 gmx_fio_do_int(fio, ir->adress->type);
1461 gmx_fio_do_real(fio, ir->adress->const_wf);
1462 gmx_fio_do_real(fio, ir->adress->ex_width);
1463 gmx_fio_do_real(fio, ir->adress->hy_width);
1464 gmx_fio_do_int(fio, ir->adress->icor);
1465 gmx_fio_do_int(fio, ir->adress->site);
1466 gmx_fio_do_rvec(fio, ir->adress->refs);
1467 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1468 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1469 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1470 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1474 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1476 if (ir->adress->n_tf_grps > 0)
1478 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1482 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1484 if (ir->adress->n_energy_grps > 0)
1486 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1492 ir->bAdress = FALSE;
1496 if (file_version >= 48)
1498 gmx_fio_do_int(fio, ir->ePull);
1499 if (ir->ePull != epullNO)
1505 do_pull(fio, ir->pull, bRead, file_version);
1510 ir->ePull = epullNO;
1513 /* Enforced rotation */
1514 if (file_version >= 74)
1516 gmx_fio_do_int(fio, ir->bRot);
1517 if (ir->bRot == TRUE)
1523 do_rot(fio, ir->rot, bRead);
1531 /* Interactive molecular dynamics */
1532 if (file_version >= tpxv_InteractiveMolecularDynamics)
1534 gmx_fio_do_int(fio, ir->bIMD);
1535 if (TRUE == ir->bIMD)
1541 do_imd(fio, ir->imd, bRead);
1546 /* We don't support IMD sessions for old .tpr files */
1551 gmx_fio_do_int(fio, ir->opts.ngtc);
1552 if (file_version >= 69)
1554 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1558 ir->opts.nhchainlength = 1;
1560 gmx_fio_do_int(fio, ir->opts.ngacc);
1561 gmx_fio_do_int(fio, ir->opts.ngfrz);
1562 gmx_fio_do_int(fio, ir->opts.ngener);
1566 snew(ir->opts.nrdf, ir->opts.ngtc);
1567 snew(ir->opts.ref_t, ir->opts.ngtc);
1568 snew(ir->opts.annealing, ir->opts.ngtc);
1569 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1570 snew(ir->opts.anneal_time, ir->opts.ngtc);
1571 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1572 snew(ir->opts.tau_t, ir->opts.ngtc);
1573 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1574 snew(ir->opts.acc, ir->opts.ngacc);
1575 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1577 if (ir->opts.ngtc > 0)
1579 if (bRead && file_version < 13)
1581 snew(tmp, ir->opts.ngtc);
1582 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1583 for (i = 0; i < ir->opts.ngtc; i++)
1585 ir->opts.nrdf[i] = tmp[i];
1591 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1593 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1594 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1595 if (file_version < 33 && ir->eI == eiBD)
1597 for (i = 0; i < ir->opts.ngtc; i++)
1599 ir->opts.tau_t[i] = bd_temp;
1603 if (ir->opts.ngfrz > 0)
1605 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1607 if (ir->opts.ngacc > 0)
1609 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1611 if (file_version >= 12)
1613 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1614 ir->opts.ngener*ir->opts.ngener);
1617 if (bRead && file_version < 26)
1619 for (i = 0; i < ir->opts.ngtc; i++)
1623 ir->opts.annealing[i] = eannSINGLE;
1624 ir->opts.anneal_npoints[i] = 2;
1625 snew(ir->opts.anneal_time[i], 2);
1626 snew(ir->opts.anneal_temp[i], 2);
1627 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1628 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1629 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1630 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1631 ir->opts.anneal_time[i][0] = ir->init_t;
1632 ir->opts.anneal_time[i][1] = finish_t;
1633 ir->opts.anneal_temp[i][0] = init_temp;
1634 ir->opts.anneal_temp[i][1] = finish_temp;
1638 ir->opts.annealing[i] = eannNO;
1639 ir->opts.anneal_npoints[i] = 0;
1645 /* file version 26 or later */
1646 /* First read the lists with annealing and npoints for each group */
1647 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1648 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1649 for (j = 0; j < (ir->opts.ngtc); j++)
1651 k = ir->opts.anneal_npoints[j];
1654 snew(ir->opts.anneal_time[j], k);
1655 snew(ir->opts.anneal_temp[j], k);
1657 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1658 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1662 if (file_version >= 45)
1664 gmx_fio_do_int(fio, ir->nwall);
1665 gmx_fio_do_int(fio, ir->wall_type);
1666 if (file_version >= 50)
1668 gmx_fio_do_real(fio, ir->wall_r_linpot);
1672 ir->wall_r_linpot = -1;
1674 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1675 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1676 gmx_fio_do_real(fio, ir->wall_density[0]);
1677 gmx_fio_do_real(fio, ir->wall_density[1]);
1678 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1684 ir->wall_atomtype[0] = -1;
1685 ir->wall_atomtype[1] = -1;
1686 ir->wall_density[0] = 0;
1687 ir->wall_density[1] = 0;
1688 ir->wall_ewald_zfac = 3;
1690 /* Cosine stuff for electric fields */
1691 for (j = 0; (j < DIM); j++)
1693 gmx_fio_do_int(fio, ir->ex[j].n);
1694 gmx_fio_do_int(fio, ir->et[j].n);
1697 snew(ir->ex[j].a, ir->ex[j].n);
1698 snew(ir->ex[j].phi, ir->ex[j].n);
1699 snew(ir->et[j].a, ir->et[j].n);
1700 snew(ir->et[j].phi, ir->et[j].n);
1702 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1703 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1704 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1705 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1709 if (file_version >= tpxv_ComputationalElectrophysiology)
1711 gmx_fio_do_int(fio, ir->eSwapCoords);
1712 if (ir->eSwapCoords != eswapNO)
1718 do_swapcoords(fio, ir->swap, bRead);
1723 if (file_version >= 39)
1725 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1726 gmx_fio_do_int(fio, ir->QMMMscheme);
1727 gmx_fio_do_real(fio, ir->scalefactor);
1728 gmx_fio_do_int(fio, ir->opts.ngQM);
1731 snew(ir->opts.QMmethod, ir->opts.ngQM);
1732 snew(ir->opts.QMbasis, ir->opts.ngQM);
1733 snew(ir->opts.QMcharge, ir->opts.ngQM);
1734 snew(ir->opts.QMmult, ir->opts.ngQM);
1735 snew(ir->opts.bSH, ir->opts.ngQM);
1736 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1737 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1738 snew(ir->opts.SAon, ir->opts.ngQM);
1739 snew(ir->opts.SAoff, ir->opts.ngQM);
1740 snew(ir->opts.SAsteps, ir->opts.ngQM);
1741 snew(ir->opts.bOPT, ir->opts.ngQM);
1742 snew(ir->opts.bTS, ir->opts.ngQM);
1744 if (ir->opts.ngQM > 0)
1746 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1747 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1748 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1749 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1750 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1751 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1752 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1753 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1754 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1755 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1756 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1757 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1759 /* end of QMMM stuff */
1764 static void do_harm(t_fileio *fio, t_iparams *iparams)
1766 gmx_fio_do_real(fio, iparams->harmonic.rA);
1767 gmx_fio_do_real(fio, iparams->harmonic.krA);
1768 gmx_fio_do_real(fio, iparams->harmonic.rB);
1769 gmx_fio_do_real(fio, iparams->harmonic.krB);
1772 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1773 gmx_bool bRead, int file_version)
1780 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1790 do_harm(fio, iparams);
1791 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1793 /* Correct incorrect storage of parameters */
1794 iparams->pdihs.phiB = iparams->pdihs.phiA;
1795 iparams->pdihs.cpB = iparams->pdihs.cpA;
1799 gmx_fio_do_real(fio, iparams->harmonic.rA);
1800 gmx_fio_do_real(fio, iparams->harmonic.krA);
1802 case F_LINEAR_ANGLES:
1803 gmx_fio_do_real(fio, iparams->linangle.klinA);
1804 gmx_fio_do_real(fio, iparams->linangle.aA);
1805 gmx_fio_do_real(fio, iparams->linangle.klinB);
1806 gmx_fio_do_real(fio, iparams->linangle.aB);
1809 gmx_fio_do_real(fio, iparams->fene.bm);
1810 gmx_fio_do_real(fio, iparams->fene.kb);
1814 gmx_fio_do_real(fio, iparams->restraint.lowA);
1815 gmx_fio_do_real(fio, iparams->restraint.up1A);
1816 gmx_fio_do_real(fio, iparams->restraint.up2A);
1817 gmx_fio_do_real(fio, iparams->restraint.kA);
1818 gmx_fio_do_real(fio, iparams->restraint.lowB);
1819 gmx_fio_do_real(fio, iparams->restraint.up1B);
1820 gmx_fio_do_real(fio, iparams->restraint.up2B);
1821 gmx_fio_do_real(fio, iparams->restraint.kB);
1827 gmx_fio_do_real(fio, iparams->tab.kA);
1828 gmx_fio_do_int(fio, iparams->tab.table);
1829 gmx_fio_do_real(fio, iparams->tab.kB);
1831 case F_CROSS_BOND_BONDS:
1832 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1833 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1834 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1836 case F_CROSS_BOND_ANGLES:
1837 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1838 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1839 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1840 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1842 case F_UREY_BRADLEY:
1843 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1844 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1845 gmx_fio_do_real(fio, iparams->u_b.r13A);
1846 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1847 if (file_version >= 79)
1849 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1850 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1851 gmx_fio_do_real(fio, iparams->u_b.r13B);
1852 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1856 iparams->u_b.thetaB = iparams->u_b.thetaA;
1857 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1858 iparams->u_b.r13B = iparams->u_b.r13A;
1859 iparams->u_b.kUBB = iparams->u_b.kUBA;
1862 case F_QUARTIC_ANGLES:
1863 gmx_fio_do_real(fio, iparams->qangle.theta);
1864 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1867 gmx_fio_do_real(fio, iparams->bham.a);
1868 gmx_fio_do_real(fio, iparams->bham.b);
1869 gmx_fio_do_real(fio, iparams->bham.c);
1872 gmx_fio_do_real(fio, iparams->morse.b0A);
1873 gmx_fio_do_real(fio, iparams->morse.cbA);
1874 gmx_fio_do_real(fio, iparams->morse.betaA);
1875 if (file_version >= 79)
1877 gmx_fio_do_real(fio, iparams->morse.b0B);
1878 gmx_fio_do_real(fio, iparams->morse.cbB);
1879 gmx_fio_do_real(fio, iparams->morse.betaB);
1883 iparams->morse.b0B = iparams->morse.b0A;
1884 iparams->morse.cbB = iparams->morse.cbA;
1885 iparams->morse.betaB = iparams->morse.betaA;
1889 gmx_fio_do_real(fio, iparams->cubic.b0);
1890 gmx_fio_do_real(fio, iparams->cubic.kb);
1891 gmx_fio_do_real(fio, iparams->cubic.kcub);
1895 case F_POLARIZATION:
1896 gmx_fio_do_real(fio, iparams->polarize.alpha);
1899 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1900 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1901 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1904 if (file_version < 31)
1906 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1908 gmx_fio_do_real(fio, iparams->wpol.al_x);
1909 gmx_fio_do_real(fio, iparams->wpol.al_y);
1910 gmx_fio_do_real(fio, iparams->wpol.al_z);
1911 gmx_fio_do_real(fio, iparams->wpol.rOH);
1912 gmx_fio_do_real(fio, iparams->wpol.rHH);
1913 gmx_fio_do_real(fio, iparams->wpol.rOD);
1916 gmx_fio_do_real(fio, iparams->thole.a);
1917 gmx_fio_do_real(fio, iparams->thole.alpha1);
1918 gmx_fio_do_real(fio, iparams->thole.alpha2);
1919 gmx_fio_do_real(fio, iparams->thole.rfac);
1922 gmx_fio_do_real(fio, iparams->lj.c6);
1923 gmx_fio_do_real(fio, iparams->lj.c12);
1926 gmx_fio_do_real(fio, iparams->lj14.c6A);
1927 gmx_fio_do_real(fio, iparams->lj14.c12A);
1928 gmx_fio_do_real(fio, iparams->lj14.c6B);
1929 gmx_fio_do_real(fio, iparams->lj14.c12B);
1932 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1933 gmx_fio_do_real(fio, iparams->ljc14.qi);
1934 gmx_fio_do_real(fio, iparams->ljc14.qj);
1935 gmx_fio_do_real(fio, iparams->ljc14.c6);
1936 gmx_fio_do_real(fio, iparams->ljc14.c12);
1938 case F_LJC_PAIRS_NB:
1939 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1940 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1941 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1942 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1948 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1949 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1950 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1952 /* Read the incorrectly stored multiplicity */
1953 gmx_fio_do_real(fio, iparams->harmonic.rB);
1954 gmx_fio_do_real(fio, iparams->harmonic.krB);
1955 iparams->pdihs.phiB = iparams->pdihs.phiA;
1956 iparams->pdihs.cpB = iparams->pdihs.cpA;
1960 gmx_fio_do_real(fio, iparams->pdihs.phiB);
1961 gmx_fio_do_real(fio, iparams->pdihs.cpB);
1962 gmx_fio_do_int(fio, iparams->pdihs.mult);
1966 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1967 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1970 gmx_fio_do_int(fio, iparams->disres.label);
1971 gmx_fio_do_int(fio, iparams->disres.type);
1972 gmx_fio_do_real(fio, iparams->disres.low);
1973 gmx_fio_do_real(fio, iparams->disres.up1);
1974 gmx_fio_do_real(fio, iparams->disres.up2);
1975 gmx_fio_do_real(fio, iparams->disres.kfac);
1978 gmx_fio_do_int(fio, iparams->orires.ex);
1979 gmx_fio_do_int(fio, iparams->orires.label);
1980 gmx_fio_do_int(fio, iparams->orires.power);
1981 gmx_fio_do_real(fio, iparams->orires.c);
1982 gmx_fio_do_real(fio, iparams->orires.obs);
1983 gmx_fio_do_real(fio, iparams->orires.kfac);
1986 if (file_version < 82)
1988 gmx_fio_do_int(fio, idum);
1989 gmx_fio_do_int(fio, idum);
1991 gmx_fio_do_real(fio, iparams->dihres.phiA);
1992 gmx_fio_do_real(fio, iparams->dihres.dphiA);
1993 gmx_fio_do_real(fio, iparams->dihres.kfacA);
1994 if (file_version >= 82)
1996 gmx_fio_do_real(fio, iparams->dihres.phiB);
1997 gmx_fio_do_real(fio, iparams->dihres.dphiB);
1998 gmx_fio_do_real(fio, iparams->dihres.kfacB);
2002 iparams->dihres.phiB = iparams->dihres.phiA;
2003 iparams->dihres.dphiB = iparams->dihres.dphiA;
2004 iparams->dihres.kfacB = iparams->dihres.kfacA;
2008 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
2009 gmx_fio_do_rvec(fio, iparams->posres.fcA);
2010 if (bRead && file_version < 27)
2012 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
2013 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
2017 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
2018 gmx_fio_do_rvec(fio, iparams->posres.fcB);
2022 gmx_fio_do_int(fio, iparams->fbposres.geom);
2023 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
2024 gmx_fio_do_real(fio, iparams->fbposres.r);
2025 gmx_fio_do_real(fio, iparams->fbposres.k);
2028 gmx_fio_ndo_real(fio, iparams->cbtdihs.cbtcA, NR_CBTDIHS);
2031 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2032 if (file_version >= 25)
2034 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2038 /* Fourier dihedrals are internally represented
2039 * as Ryckaert-Bellemans since those are faster to compute.
2041 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2042 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2046 gmx_fio_do_real(fio, iparams->constr.dA);
2047 gmx_fio_do_real(fio, iparams->constr.dB);
2050 gmx_fio_do_real(fio, iparams->settle.doh);
2051 gmx_fio_do_real(fio, iparams->settle.dhh);
2054 gmx_fio_do_real(fio, iparams->vsite.a);
2059 gmx_fio_do_real(fio, iparams->vsite.a);
2060 gmx_fio_do_real(fio, iparams->vsite.b);
2065 gmx_fio_do_real(fio, iparams->vsite.a);
2066 gmx_fio_do_real(fio, iparams->vsite.b);
2067 gmx_fio_do_real(fio, iparams->vsite.c);
2070 gmx_fio_do_int(fio, iparams->vsiten.n);
2071 gmx_fio_do_real(fio, iparams->vsiten.a);
2076 /* We got rid of some parameters in version 68 */
2077 if (bRead && file_version < 68)
2079 gmx_fio_do_real(fio, rdum);
2080 gmx_fio_do_real(fio, rdum);
2081 gmx_fio_do_real(fio, rdum);
2082 gmx_fio_do_real(fio, rdum);
2084 gmx_fio_do_real(fio, iparams->gb.sar);
2085 gmx_fio_do_real(fio, iparams->gb.st);
2086 gmx_fio_do_real(fio, iparams->gb.pi);
2087 gmx_fio_do_real(fio, iparams->gb.gbr);
2088 gmx_fio_do_real(fio, iparams->gb.bmlt);
2091 gmx_fio_do_int(fio, iparams->cmap.cmapA);
2092 gmx_fio_do_int(fio, iparams->cmap.cmapB);
2095 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
2096 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
2100 gmx_fio_unset_comment(fio);
2104 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
2111 gmx_fio_set_comment(fio, interaction_function[ftype].name);
2113 if (file_version < 44)
2115 for (i = 0; i < MAXNODES; i++)
2117 gmx_fio_do_int(fio, idum);
2120 gmx_fio_do_int(fio, ilist->nr);
2123 snew(ilist->iatoms, ilist->nr);
2125 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2128 gmx_fio_unset_comment(fio);
2132 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2133 gmx_bool bRead, int file_version)
2138 gmx_fio_do_int(fio, ffparams->atnr);
2139 if (file_version < 57)
2141 gmx_fio_do_int(fio, idum);
2143 gmx_fio_do_int(fio, ffparams->ntypes);
2146 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2147 ffparams->atnr, ffparams->ntypes);
2151 snew(ffparams->functype, ffparams->ntypes);
2152 snew(ffparams->iparams, ffparams->ntypes);
2154 /* Read/write all the function types */
2155 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2158 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2161 if (file_version >= 66)
2163 gmx_fio_do_double(fio, ffparams->reppow);
2167 ffparams->reppow = 12.0;
2170 if (file_version >= 57)
2172 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2175 /* Check whether all these function types are supported by the code.
2176 * In practice the code is backwards compatible, which means that the
2177 * numbering may have to be altered from old numbering to new numbering
2179 for (i = 0; (i < ffparams->ntypes); i++)
2183 /* Loop over file versions */
2184 for (k = 0; (k < NFTUPD); k++)
2186 /* Compare the read file_version to the update table */
2187 if ((file_version < ftupd[k].fvnr) &&
2188 (ffparams->functype[i] >= ftupd[k].ftype))
2190 ffparams->functype[i] += 1;
2193 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2194 i, ffparams->functype[i],
2195 interaction_function[ftupd[k].ftype].longname);
2202 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2206 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2211 static void add_settle_atoms(t_ilist *ilist)
2215 /* Settle used to only store the first atom: add the other two */
2216 srenew(ilist->iatoms, 2*ilist->nr);
2217 for (i = ilist->nr/2-1; i >= 0; i--)
2219 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2220 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2221 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2222 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2224 ilist->nr = 2*ilist->nr;
2227 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2230 int i, j, renum[F_NRE];
2234 for (j = 0; (j < F_NRE); j++)
2239 for (k = 0; k < NFTUPD; k++)
2241 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2250 ilist[j].iatoms = NULL;
2254 do_ilist(fio, &ilist[j], bRead, file_version, j);
2255 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2257 add_settle_atoms(&ilist[j]);
2261 if (bRead && gmx_debug_at)
2262 pr_ilist(debug,0,interaction_function[j].longname,
2263 functype,&ilist[j],TRUE);
2268 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2269 gmx_bool bRead, int file_version)
2271 do_ffparams(fio, ffparams, bRead, file_version);
2273 if (file_version >= 54)
2275 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2278 do_ilists(fio, molt->ilist, bRead, file_version);
2281 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2283 int i, idum, dum_nra, *dum_a;
2285 if (file_version < 44)
2287 for (i = 0; i < MAXNODES; i++)
2289 gmx_fio_do_int(fio, idum);
2292 gmx_fio_do_int(fio, block->nr);
2293 if (file_version < 51)
2295 gmx_fio_do_int(fio, dum_nra);
2299 if ((block->nalloc_index > 0) && (NULL != block->index))
2301 sfree(block->index);
2303 block->nalloc_index = block->nr+1;
2304 snew(block->index, block->nalloc_index);
2306 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2308 if (file_version < 51 && dum_nra > 0)
2310 snew(dum_a, dum_nra);
2311 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2316 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2321 if (file_version < 44)
2323 for (i = 0; i < MAXNODES; i++)
2325 gmx_fio_do_int(fio, idum);
2328 gmx_fio_do_int(fio, block->nr);
2329 gmx_fio_do_int(fio, block->nra);
2332 block->nalloc_index = block->nr+1;
2333 snew(block->index, block->nalloc_index);
2334 block->nalloc_a = block->nra;
2335 snew(block->a, block->nalloc_a);
2337 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2338 gmx_fio_ndo_int(fio, block->a, block->nra);
2341 /* This is a primitive routine to make it possible to translate atomic numbers
2342 * to element names when reading TPR files, without making the Gromacs library
2343 * directory a dependency on mdrun (which is the case if we need elements.dat).
2346 atomicnumber_to_element(int atomicnumber)
2350 /* This does not have to be complete, so we only include elements likely
2351 * to occur in PDB files.
2353 switch (atomicnumber)
2355 case 1: p = "H"; break;
2356 case 5: p = "B"; break;
2357 case 6: p = "C"; break;
2358 case 7: p = "N"; break;
2359 case 8: p = "O"; break;
2360 case 9: p = "F"; break;
2361 case 11: p = "Na"; break;
2362 case 12: p = "Mg"; break;
2363 case 15: p = "P"; break;
2364 case 16: p = "S"; break;
2365 case 17: p = "Cl"; break;
2366 case 18: p = "Ar"; break;
2367 case 19: p = "K"; break;
2368 case 20: p = "Ca"; break;
2369 case 25: p = "Mn"; break;
2370 case 26: p = "Fe"; break;
2371 case 28: p = "Ni"; break;
2372 case 29: p = "Cu"; break;
2373 case 30: p = "Zn"; break;
2374 case 35: p = "Br"; break;
2375 case 47: p = "Ag"; break;
2376 default: p = ""; break;
2382 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2383 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);
2401 /* Set element string from atomic number if present.
2402 * This routine returns an empty string if the name is not found.
2404 strncpy(atom->elem, atomicnumber_to_element(atom->atomnumber), 4);
2405 /* avoid warnings about potentially unterminated string */
2406 atom->elem[3] = '\0';
2411 atom->atomnumber = NOTSET;
2413 if (file_version < 23)
2417 else if (file_version < 39)
2426 if (file_version < 57)
2428 unsigned char uchar[egcNR];
2429 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2430 for (i = myngrp; (i < ngrp); i++)
2434 /* Copy the old data format to the groups struct */
2435 for (i = 0; i < ngrp; i++)
2437 groups->grpnr[i][atnr] = uchar[i];
2442 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2447 if (file_version < 23)
2451 else if (file_version < 39)
2460 for (j = 0; (j < ngrp); j++)
2464 gmx_fio_do_int(fio, grps[j].nr);
2467 snew(grps[j].nm_ind, grps[j].nr);
2469 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2474 snew(grps[j].nm_ind, grps[j].nr);
2479 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2485 gmx_fio_do_int(fio, ls);
2486 *nm = get_symtab_handle(symtab, ls);
2490 ls = lookup_symtab(symtab, *nm);
2491 gmx_fio_do_int(fio, ls);
2495 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2500 for (j = 0; (j < nstr); j++)
2502 do_symstr(fio, &(nm[j]), bRead, symtab);
2506 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2507 t_symtab *symtab, int file_version)
2511 for (j = 0; (j < n); j++)
2513 do_symstr(fio, &(ri[j].name), bRead, symtab);
2514 if (file_version >= 63)
2516 gmx_fio_do_int(fio, ri[j].nr);
2517 gmx_fio_do_uchar(fio, ri[j].ic);
2527 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2529 gmx_groups_t *groups)
2533 gmx_fio_do_int(fio, atoms->nr);
2534 gmx_fio_do_int(fio, atoms->nres);
2535 if (file_version < 57)
2537 gmx_fio_do_int(fio, groups->ngrpname);
2538 for (i = 0; i < egcNR; i++)
2540 groups->ngrpnr[i] = atoms->nr;
2541 snew(groups->grpnr[i], groups->ngrpnr[i]);
2546 snew(atoms->atom, atoms->nr);
2547 snew(atoms->atomname, atoms->nr);
2548 snew(atoms->atomtype, atoms->nr);
2549 snew(atoms->atomtypeB, atoms->nr);
2550 snew(atoms->resinfo, atoms->nres);
2551 if (file_version < 57)
2553 snew(groups->grpname, groups->ngrpname);
2555 atoms->pdbinfo = NULL;
2557 for (i = 0; (i < atoms->nr); i++)
2559 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2561 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2562 if (bRead && (file_version <= 20))
2564 for (i = 0; i < atoms->nr; i++)
2566 atoms->atomtype[i] = put_symtab(symtab, "?");
2567 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2572 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2573 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2575 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2577 if (file_version < 57)
2579 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2581 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2585 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2586 gmx_bool bRead, t_symtab *symtab,
2591 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2592 gmx_fio_do_int(fio, groups->ngrpname);
2595 snew(groups->grpname, groups->ngrpname);
2597 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2598 for (g = 0; g < egcNR; g++)
2600 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2601 if (groups->ngrpnr[g] == 0)
2605 groups->grpnr[g] = NULL;
2612 snew(groups->grpnr[g], groups->ngrpnr[g]);
2614 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2619 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2624 if (file_version > 25)
2626 gmx_fio_do_int(fio, atomtypes->nr);
2630 snew(atomtypes->radius, j);
2631 snew(atomtypes->vol, j);
2632 snew(atomtypes->surftens, j);
2633 snew(atomtypes->atomnumber, j);
2634 snew(atomtypes->gb_radius, j);
2635 snew(atomtypes->S_hct, j);
2637 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2638 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2639 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2640 if (file_version >= 40)
2642 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2644 if (file_version >= 60)
2646 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2647 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2652 /* File versions prior to 26 cannot do GBSA,
2653 * so they dont use this structure
2656 atomtypes->radius = NULL;
2657 atomtypes->vol = NULL;
2658 atomtypes->surftens = NULL;
2659 atomtypes->atomnumber = NULL;
2660 atomtypes->gb_radius = NULL;
2661 atomtypes->S_hct = NULL;
2665 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2671 gmx_fio_do_int(fio, symtab->nr);
2675 snew(symtab->symbuf, 1);
2676 symbuf = symtab->symbuf;
2677 symbuf->bufsize = nr;
2678 snew(symbuf->buf, nr);
2679 for (i = 0; (i < nr); i++)
2681 gmx_fio_do_string(fio, buf);
2682 symbuf->buf[i] = gmx_strdup(buf);
2687 symbuf = symtab->symbuf;
2688 while (symbuf != NULL)
2690 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2692 gmx_fio_do_string(fio, symbuf->buf[i]);
2695 symbuf = symbuf->next;
2699 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2704 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2706 int i, j, ngrid, gs, nelem;
2708 gmx_fio_do_int(fio, cmap_grid->ngrid);
2709 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2711 ngrid = cmap_grid->ngrid;
2712 gs = cmap_grid->grid_spacing;
2717 snew(cmap_grid->cmapdata, ngrid);
2719 for (i = 0; i < cmap_grid->ngrid; i++)
2721 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2725 for (i = 0; i < cmap_grid->ngrid; i++)
2727 for (j = 0; j < nelem; j++)
2729 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2730 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2731 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2732 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2738 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2740 int m, a, a0, a1, r;
2744 /* We always assign a new chain number, but save the chain id characters
2745 * for larger molecules.
2747 #define CHAIN_MIN_ATOMS 15
2751 for (m = 0; m < mols->nr; m++)
2753 a0 = mols->index[m];
2754 a1 = mols->index[m+1];
2755 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2764 for (a = a0; a < a1; a++)
2766 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2767 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2772 /* Blank out the chain id if there was only one chain */
2775 for (r = 0; r < atoms->nres; r++)
2777 atoms->resinfo[r].chainid = ' ';
2782 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2783 t_symtab *symtab, int file_version,
2784 gmx_groups_t *groups)
2788 if (file_version >= 57)
2790 do_symstr(fio, &(molt->name), bRead, symtab);
2793 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2795 if (bRead && gmx_debug_at)
2797 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2800 if (file_version >= 57)
2802 do_ilists(fio, molt->ilist, bRead, file_version);
2804 do_block(fio, &molt->cgs, bRead, file_version);
2805 if (bRead && gmx_debug_at)
2807 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2811 /* This used to be in the atoms struct */
2812 do_blocka(fio, &molt->excls, bRead, file_version);
2815 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2819 gmx_fio_do_int(fio, molb->type);
2820 gmx_fio_do_int(fio, molb->nmol);
2821 gmx_fio_do_int(fio, molb->natoms_mol);
2822 /* Position restraint coordinates */
2823 gmx_fio_do_int(fio, molb->nposres_xA);
2824 if (molb->nposres_xA > 0)
2828 snew(molb->posres_xA, molb->nposres_xA);
2830 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2832 gmx_fio_do_int(fio, molb->nposres_xB);
2833 if (molb->nposres_xB > 0)
2837 snew(molb->posres_xB, molb->nposres_xB);
2839 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2844 static t_block mtop_mols(gmx_mtop_t *mtop)
2850 for (mb = 0; mb < mtop->nmolblock; mb++)
2852 mols.nr += mtop->molblock[mb].nmol;
2854 mols.nalloc_index = mols.nr + 1;
2855 snew(mols.index, mols.nalloc_index);
2860 for (mb = 0; mb < mtop->nmolblock; mb++)
2862 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2864 a += mtop->molblock[mb].natoms_mol;
2873 static void add_posres_molblock(gmx_mtop_t *mtop)
2878 gmx_molblock_t *molb;
2881 /* posres reference positions are stored in ip->posres (if present) and
2882 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2883 posres.pos0A are identical to fbposres.pos0. */
2884 il = &mtop->moltype[0].ilist[F_POSRES];
2885 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2886 if (il->nr == 0 && ilfb->nr == 0)
2892 for (i = 0; i < il->nr; i += 2)
2894 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2895 am = max(am, il->iatoms[i+1]);
2896 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2897 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2898 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2903 /* This loop is required if we have only flat-bottomed posres:
2905 - bFE == FALSE (no B-state for flat-bottomed posres) */
2908 for (i = 0; i < ilfb->nr; i += 2)
2910 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2911 am = max(am, ilfb->iatoms[i+1]);
2914 /* Make the posres coordinate block end at a molecule end */
2916 while (am >= mtop->mols.index[mol+1])
2920 molb = &mtop->molblock[0];
2921 molb->nposres_xA = mtop->mols.index[mol+1];
2922 snew(molb->posres_xA, molb->nposres_xA);
2925 molb->nposres_xB = molb->nposres_xA;
2926 snew(molb->posres_xB, molb->nposres_xB);
2930 molb->nposres_xB = 0;
2932 for (i = 0; i < il->nr; i += 2)
2934 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2935 a = il->iatoms[i+1];
2936 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2937 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2938 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2941 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2942 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2943 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2948 /* If only flat-bottomed posres are present, take reference pos from them.
2949 Here: bFE == FALSE */
2950 for (i = 0; i < ilfb->nr; i += 2)
2952 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2953 a = ilfb->iatoms[i+1];
2954 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2955 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2956 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2961 static void set_disres_npair(gmx_mtop_t *mtop)
2968 ip = mtop->ffparams.iparams;
2970 for (mt = 0; mt < mtop->nmoltype; mt++)
2972 il = &mtop->moltype[mt].ilist[F_DISRES];
2977 for (i = 0; i < il->nr; i += 3)
2980 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2982 ip[a[i]].disres.npair = npair;
2990 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
3000 do_symtab(fio, &(mtop->symtab), bRead);
3003 pr_symtab(debug, 0, "symtab", &mtop->symtab);
3006 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
3008 if (file_version >= 57)
3010 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
3012 gmx_fio_do_int(fio, mtop->nmoltype);
3020 snew(mtop->moltype, mtop->nmoltype);
3021 if (file_version < 57)
3023 mtop->moltype[0].name = mtop->name;
3026 for (mt = 0; mt < mtop->nmoltype; mt++)
3028 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
3032 if (file_version >= 57)
3034 gmx_fio_do_int(fio, mtop->nmolblock);
3038 mtop->nmolblock = 1;
3042 snew(mtop->molblock, mtop->nmolblock);
3044 if (file_version >= 57)
3046 for (mb = 0; mb < mtop->nmolblock; mb++)
3048 do_molblock(fio, &mtop->molblock[mb], bRead);
3050 gmx_fio_do_int(fio, mtop->natoms);
3054 mtop->molblock[0].type = 0;
3055 mtop->molblock[0].nmol = 1;
3056 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
3057 mtop->molblock[0].nposres_xA = 0;
3058 mtop->molblock[0].nposres_xB = 0;
3061 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
3064 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
3067 if (file_version < 57)
3069 /* Debug statements are inside do_idef */
3070 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
3071 mtop->natoms = mtop->moltype[0].atoms.nr;
3074 if (file_version >= 65)
3076 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
3080 mtop->ffparams.cmap_grid.ngrid = 0;
3081 mtop->ffparams.cmap_grid.grid_spacing = 0;
3082 mtop->ffparams.cmap_grid.cmapdata = NULL;
3085 if (file_version >= 57)
3087 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
3090 if (file_version < 57)
3092 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
3093 if (bRead && gmx_debug_at)
3095 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
3097 do_block(fio, &mtop->mols, bRead, file_version);
3098 /* Add the posres coordinates to the molblock */
3099 add_posres_molblock(mtop);
3103 if (file_version >= 57)
3105 done_block(&mtop->mols);
3106 mtop->mols = mtop_mols(mtop);
3110 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
3114 if (file_version < 51)
3116 /* Here used to be the shake blocks */
3117 do_blocka(fio, &dumb, bRead, file_version);
3130 close_symtab(&(mtop->symtab));
3134 /* If TopOnlyOK is TRUE then we can read even future versions
3135 * of tpx files, provided the file_generation hasn't changed.
3136 * If it is FALSE, we need the inputrecord too, and bail out
3137 * if the file is newer than the program.
3139 * The version and generation if the topology (see top of this file)
3140 * are returned in the two last arguments.
3142 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
3144 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
3145 gmx_bool TopOnlyOK, int *file_version,
3146 int *file_generation)
3149 char file_tag[STRLEN];
3156 gmx_fio_checktype(fio);
3157 gmx_fio_setdebug(fio, bDebugMode());
3159 /* XDR binary topology file */
3160 precision = sizeof(real);
3163 gmx_fio_do_string(fio, buf);
3164 if (strncmp(buf, "VERSION", 7))
3166 gmx_fatal(FARGS, "Can not read file %s,\n"
3167 " this file is from a Gromacs version which is older than 2.0\n"
3168 " Make a new one with grompp or use a gro or pdb file, if possible",
3169 gmx_fio_getname(fio));
3171 gmx_fio_do_int(fio, precision);
3172 bDouble = (precision == sizeof(double));
3173 if ((precision != sizeof(float)) && !bDouble)
3175 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3176 "instead of %d or %d",
3177 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3179 gmx_fio_setprecision(fio, bDouble);
3180 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3181 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3185 gmx_fio_write_string(fio, GromacsVersion());
3186 bDouble = (precision == sizeof(double));
3187 gmx_fio_setprecision(fio, bDouble);
3188 gmx_fio_do_int(fio, precision);
3190 sprintf(file_tag, "%s", tpx_tag);
3191 fgen = tpx_generation;
3194 /* Check versions! */
3195 gmx_fio_do_int(fio, fver);
3197 /* This is for backward compatibility with development versions 77-79
3198 * where the tag was, mistakenly, placed before the generation,
3199 * which would cause a segv instead of a proper error message
3200 * when reading the topology only from tpx with <77 code.
3202 if (fver >= 77 && fver <= 79)
3204 gmx_fio_do_string(fio, file_tag);
3209 gmx_fio_do_int(fio, fgen);
3218 gmx_fio_do_string(fio, file_tag);
3224 /* Versions before 77 don't have the tag, set it to release */
3225 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3228 if (strcmp(file_tag, tpx_tag) != 0)
3230 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3233 /* We only support reading tpx files with the same tag as the code
3234 * or tpx files with the release tag and with lower version number.
3236 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
3238 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3239 gmx_fio_getname(fio), fver, file_tag,
3240 tpx_version, tpx_tag);
3245 if (file_version != NULL)
3247 *file_version = fver;
3249 if (file_generation != NULL)
3251 *file_generation = fgen;
3255 if ((fver <= tpx_incompatible_version) ||
3256 ((fver > tpx_version) && !TopOnlyOK) ||
3257 (fgen > tpx_generation) ||
3258 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3260 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3261 gmx_fio_getname(fio), fver, tpx_version);
3264 gmx_fio_do_int(fio, tpx->natoms);
3267 gmx_fio_do_int(fio, tpx->ngtc);
3275 gmx_fio_do_int(fio, idum);
3276 gmx_fio_do_real(fio, rdum);
3278 /*a better decision will eventually (5.0 or later) need to be made
3279 on how to treat the alchemical state of the system, which can now
3280 vary through a simulation, and cannot be completely described
3281 though a single lambda variable, or even a single state
3282 index. Eventually, should probably be a vector. MRS*/
3285 gmx_fio_do_int(fio, tpx->fep_state);
3287 gmx_fio_do_real(fio, tpx->lambda);
3288 gmx_fio_do_int(fio, tpx->bIr);
3289 gmx_fio_do_int(fio, tpx->bTop);
3290 gmx_fio_do_int(fio, tpx->bX);
3291 gmx_fio_do_int(fio, tpx->bV);
3292 gmx_fio_do_int(fio, tpx->bF);
3293 gmx_fio_do_int(fio, tpx->bBox);
3295 if ((fgen > tpx_generation))
3297 /* This can only happen if TopOnlyOK=TRUE */
3302 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3303 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3304 gmx_bool bXVallocated)
3310 int file_version, file_generation;
3314 gmx_bool bPeriodicMols;
3318 tpx.natoms = state->natoms;
3319 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3320 tpx.fep_state = state->fep_state;
3321 tpx.lambda = state->lambda[efptFEP];
3322 tpx.bIr = (ir != NULL);
3323 tpx.bTop = (mtop != NULL);
3324 tpx.bX = (state->x != NULL);
3325 tpx.bV = (state->v != NULL);
3326 tpx.bF = (f != NULL);
3330 TopOnlyOK = (ir == NULL);
3332 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3337 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3338 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3343 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3344 state->natoms = tpx.natoms;
3345 state->nalloc = tpx.natoms;
3351 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3355 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3357 do_test(fio, tpx.bBox, state->box);
3360 gmx_fio_ndo_rvec(fio, state->box, DIM);
3361 if (file_version >= 51)
3363 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3367 /* We initialize box_rel after reading the inputrec */
3368 clear_mat(state->box_rel);
3370 if (file_version >= 28)
3372 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3373 if (file_version < 56)
3376 gmx_fio_ndo_rvec(fio, mdum, DIM);
3381 if (state->ngtc > 0 && file_version >= 28)
3384 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3385 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3386 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3387 snew(dumv, state->ngtc);
3388 if (file_version < 69)
3390 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3392 /* These used to be the Berendsen tcoupl_lambda's */
3393 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3397 /* Prior to tpx version 26, the inputrec was here.
3398 * I moved it to enable partial forward-compatibility
3399 * for analysis/viewer programs.
3401 if (file_version < 26)
3403 do_test(fio, tpx.bIr, ir);
3408 do_inputrec(fio, ir, bRead, file_version,
3409 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3412 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3417 do_inputrec(fio, &dum_ir, bRead, file_version,
3418 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3421 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3423 done_inputrec(&dum_ir);
3429 do_test(fio, tpx.bTop, mtop);
3434 do_mtop(fio, mtop, bRead, file_version);
3438 do_mtop(fio, &dum_top, bRead, file_version);
3439 done_mtop(&dum_top, TRUE);
3442 do_test(fio, tpx.bX, state->x);
3447 state->flags |= (1<<estX);
3449 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3452 do_test(fio, tpx.bV, state->v);
3457 state->flags |= (1<<estV);
3459 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3462 do_test(fio, tpx.bF, f);
3465 gmx_fio_ndo_rvec(fio, f, state->natoms);
3468 /* Starting with tpx version 26, we have the inputrec
3469 * at the end of the file, so we can ignore it
3470 * if the file is never than the software (but still the
3471 * same generation - see comments at the top of this file.
3476 bPeriodicMols = FALSE;
3477 if (file_version >= 26)
3479 do_test(fio, tpx.bIr, ir);
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 return (efTPR == fn2ftp(file));
3666 static void done_gmx_groups_t(gmx_groups_t *g)
3670 for (i = 0; (i < egcNR); i++)
3672 if (NULL != g->grps[i].nm_ind)
3674 sfree(g->grps[i].nm_ind);
3675 g->grps[i].nm_ind = NULL;
3677 if (NULL != g->grpnr[i])
3683 /* The contents of this array is in symtab, don't free it here */
3687 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3688 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3691 int natoms, i, version, generation;
3692 gmx_bool bTop, bXNULL = FALSE;
3694 t_topology *topconv;
3697 bTop = fn2bTPX(infile);
3701 read_tpxheader(infile, &header, TRUE, &version, &generation);
3704 snew(*x, header.natoms);
3708 snew(*v, header.natoms);
3711 *ePBC = read_tpx(infile, NULL, box, &natoms,
3712 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3713 *top = gmx_mtop_t_to_t_topology(mtop);
3714 /* In this case we need to throw away the group data too */
3715 done_gmx_groups_t(&mtop->groups);
3717 strcpy(title, *top->name);
3718 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3722 get_stx_coordnum(infile, &natoms);
3723 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3734 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3742 aps = gmx_atomprop_init();
3743 for (i = 0; (i < natoms); i++)
3745 if (!gmx_atomprop_query(aps, epropMass,
3746 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3747 *top->atoms.atomname[i],
3748 &(top->atoms.atom[i].m)))
3752 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3753 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3754 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3755 *top->atoms.atomname[i]);
3759 gmx_atomprop_destroy(aps);
3761 top->idef.ntypes = -1;