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! */
47 #include "gmx_fatal.h"
60 #include "mtop_util.h"
62 #define TPX_TAG_RELEASE "release"
64 /*! \brief Tag string for the file format written to run input files
65 * written by this version of the code.
67 * Change this if you want to change the run input format in a feature
68 * branch. This ensures that there will not be different run input
69 * formats around which cannot be distinguished, while not causing
70 * problems rebasing the feature branch onto upstream changes. When
71 * merging with mainstream GROMACS, set this tag string back to
72 * TPX_TAG_RELEASE, and instead add an element to tpxv and set
73 * tpx_version to that.
75 static const char *tpx_tag = TPX_TAG_RELEASE;
77 /*! \brief Enum of values that describe the contents of a tpr file
78 * whose format matches a version number
80 * The enum helps the code be more self-documenting and ensure merges
81 * do not silently resolve when two patches make the same bump. When
82 * adding new functionality, add a new element to the end of this
83 * enumeration, change the definition of tpx_version, and write code
84 * below that does the right thing according to the value of
87 tpxv_ComputationalElectrophysiology = 96, /**< support for ion/water position swaps (computational electrophysiology) */
88 tpxv_Use64BitRandomSeed /**< change ld_seed from int to gmx_int64_t */
91 /*! \brief Version number of the file format written to run input
92 * files by this version of the code.
94 * The tpx_version number should be increased whenever the file format
95 * in the main development branch changes, generally to the highest
96 * value present in tpxv. Backward compatibility for reading old run
97 * input files is maintained by checking this version number against
98 * that of the file and then using the correct code path.
100 * When developing a feature branch that needs to change the run input
101 * file format, change tpx_tag instead. */
102 static const int tpx_version = tpxv_Use64BitRandomSeed;
105 /* This number should only be increased when you edit the TOPOLOGY section
106 * or the HEADER of the tpx format.
107 * This way we can maintain forward compatibility too for all analysis tools
108 * and/or external programs that only need to know the atom/residue names,
109 * charges, and bond connectivity.
111 * It first appeared in tpx version 26, when I also moved the inputrecord
112 * to the end of the tpx file, so we can just skip it if we only
115 static const int tpx_generation = 25;
117 /* This number should be the most recent backwards incompatible version
118 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
120 static const int tpx_incompatible_version = 9;
124 /* Struct used to maintain tpx compatibility when function types are added */
126 int fvnr; /* file version number in which the function type first appeared */
127 int ftype; /* function type */
131 * The entries should be ordered in:
132 * 1. ascending file version number
133 * 2. ascending function type number
135 /*static const t_ftupd ftupd[] = {
136 { 20, F_CUBICBONDS },
140 { 22, F_DISRESVIOL },
146 { 26, F_DIHRESVIOL },
147 { 30, F_CROSS_BOND_BONDS },
148 { 30, F_CROSS_BOND_ANGLES },
149 { 30, F_UREY_BRADLEY },
150 { 30, F_POLARIZATION },
154 * The entries should be ordered in:
155 * 1. ascending function type number
156 * 2. ascending file version number
158 /* question; what is the purpose of the commented code above? */
159 static const t_ftupd ftupd[] = {
160 { 20, F_CUBICBONDS },
165 { 43, F_TABBONDSNC },
166 { 70, F_RESTRBONDS },
167 { 76, F_LINEAR_ANGLES },
168 { 30, F_CROSS_BOND_BONDS },
169 { 30, F_CROSS_BOND_ANGLES },
170 { 30, F_UREY_BRADLEY },
171 { 34, F_QUARTIC_ANGLES },
181 { 72, F_NPSOLVATION },
183 { 41, F_LJC_PAIRS_NB },
186 { 32, F_COUL_RECIP },
189 { 30, F_POLARIZATION },
192 { 22, F_DISRESVIOL },
196 { 26, F_DIHRESVIOL },
201 { 46, F_ECONSERVED },
202 { 69, F_VTEMP_NOLONGERUSED},
204 { 54, F_DVDL_CONSTR },
205 { 76, F_ANHARM_POL },
208 { 79, F_DVDL_BONDED, },
209 { 79, F_DVDL_RESTRAINT },
210 { 79, F_DVDL_TEMPERATURE },
212 #define NFTUPD asize(ftupd)
214 /* Needed for backward compatibility */
217 static void _do_section(t_fileio *fio, int key, gmx_bool bRead, const char *src,
223 if (gmx_fio_getftp(fio) == efTPA)
227 gmx_fio_write_string(fio, itemstr[key]);
228 bDbg = gmx_fio_getdebug(fio);
229 gmx_fio_setdebug(fio, FALSE);
230 gmx_fio_write_string(fio, comment_str[key]);
231 gmx_fio_setdebug(fio, bDbg);
235 if (gmx_fio_getdebug(fio))
237 fprintf(stderr, "Looking for section %s (%s, %d)",
238 itemstr[key], src, line);
243 gmx_fio_do_string(fio, buf);
245 while ((gmx_strcasecmp(buf, itemstr[key]) != 0));
247 if (gmx_strcasecmp(buf, itemstr[key]) != 0)
249 gmx_fatal(FARGS, "\nCould not find section heading %s", itemstr[key]);
251 else if (gmx_fio_getdebug(fio))
253 fprintf(stderr, " and found it\n");
259 #define do_section(fio, key, bRead) _do_section(fio, key, bRead, __FILE__, __LINE__)
261 /**************************************************************
263 * Now the higer level routines that do io of the structures and arrays
265 **************************************************************/
266 static void do_pullgrp_tpx_pre95(t_fileio *fio,
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);
288 gmx_fio_do_rvec(fio, pcrd->vec);
289 clear_rvec(pcrd->origin);
290 gmx_fio_do_rvec(fio, tmp);
292 gmx_fio_do_real(fio, pcrd->rate);
293 gmx_fio_do_real(fio, pcrd->k);
294 if (file_version >= 56)
296 gmx_fio_do_real(fio, pcrd->kB);
304 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
308 gmx_fio_do_int(fio, pgrp->nat);
311 snew(pgrp->ind, pgrp->nat);
313 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
314 gmx_fio_do_int(fio, pgrp->nweight);
317 snew(pgrp->weight, pgrp->nweight);
319 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
320 gmx_fio_do_int(fio, pgrp->pbcatom);
323 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd)
327 gmx_fio_do_int(fio, pcrd->group[0]);
328 gmx_fio_do_int(fio, pcrd->group[1]);
329 gmx_fio_do_rvec(fio, pcrd->origin);
330 gmx_fio_do_rvec(fio, pcrd->vec);
331 gmx_fio_do_real(fio, pcrd->init);
332 gmx_fio_do_real(fio, pcrd->rate);
333 gmx_fio_do_real(fio, pcrd->k);
334 gmx_fio_do_real(fio, pcrd->kB);
337 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
339 /* i is used in the ndo_double macro*/
343 int n_lambda = fepvals->n_lambda;
345 /* reset the lambda calculation window */
346 fepvals->lambda_start_n = 0;
347 fepvals->lambda_stop_n = n_lambda;
348 if (file_version >= 79)
354 snew(expand->init_lambda_weights, n_lambda);
356 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
357 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
360 gmx_fio_do_int(fio, expand->nstexpanded);
361 gmx_fio_do_int(fio, expand->elmcmove);
362 gmx_fio_do_int(fio, expand->elamstats);
363 gmx_fio_do_int(fio, expand->lmc_repeats);
364 gmx_fio_do_int(fio, expand->gibbsdeltalam);
365 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
366 gmx_fio_do_int(fio, expand->lmc_seed);
367 gmx_fio_do_real(fio, expand->mc_temp);
368 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
369 gmx_fio_do_int(fio, expand->nstTij);
370 gmx_fio_do_int(fio, expand->minvarmin);
371 gmx_fio_do_int(fio, expand->c_range);
372 gmx_fio_do_real(fio, expand->wl_scale);
373 gmx_fio_do_real(fio, expand->wl_ratio);
374 gmx_fio_do_real(fio, expand->init_wl_delta);
375 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
376 gmx_fio_do_int(fio, expand->elmceq);
377 gmx_fio_do_int(fio, expand->equil_steps);
378 gmx_fio_do_int(fio, expand->equil_samples);
379 gmx_fio_do_int(fio, expand->equil_n_at_lam);
380 gmx_fio_do_real(fio, expand->equil_wl_delta);
381 gmx_fio_do_real(fio, expand->equil_ratio);
385 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
388 if (file_version >= 79)
390 gmx_fio_do_int(fio, simtemp->eSimTempScale);
391 gmx_fio_do_real(fio, simtemp->simtemp_high);
392 gmx_fio_do_real(fio, simtemp->simtemp_low);
397 snew(simtemp->temperatures, n_lambda);
399 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
404 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
406 /* i is defined in the ndo_double macro; use g to iterate. */
411 /* free energy values */
413 if (file_version >= 79)
415 gmx_fio_do_int(fio, fepvals->init_fep_state);
416 gmx_fio_do_double(fio, fepvals->init_lambda);
417 gmx_fio_do_double(fio, fepvals->delta_lambda);
419 else if (file_version >= 59)
421 gmx_fio_do_double(fio, fepvals->init_lambda);
422 gmx_fio_do_double(fio, fepvals->delta_lambda);
426 gmx_fio_do_real(fio, rdum);
427 fepvals->init_lambda = rdum;
428 gmx_fio_do_real(fio, rdum);
429 fepvals->delta_lambda = rdum;
431 if (file_version >= 79)
433 gmx_fio_do_int(fio, fepvals->n_lambda);
436 snew(fepvals->all_lambda, efptNR);
438 for (g = 0; g < efptNR; g++)
440 if (fepvals->n_lambda > 0)
444 snew(fepvals->all_lambda[g], fepvals->n_lambda);
446 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
447 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
449 else if (fepvals->init_lambda >= 0)
451 fepvals->separate_dvdl[efptFEP] = TRUE;
455 else if (file_version >= 64)
457 gmx_fio_do_int(fio, fepvals->n_lambda);
462 snew(fepvals->all_lambda, efptNR);
463 /* still allocate the all_lambda array's contents. */
464 for (g = 0; g < efptNR; g++)
466 if (fepvals->n_lambda > 0)
468 snew(fepvals->all_lambda[g], fepvals->n_lambda);
472 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
474 if (fepvals->init_lambda >= 0)
478 fepvals->separate_dvdl[efptFEP] = TRUE;
482 /* copy the contents of the efptFEP lambda component to all
483 the other components */
484 for (g = 0; g < efptNR; g++)
486 for (h = 0; h < fepvals->n_lambda; h++)
490 fepvals->all_lambda[g][h] =
491 fepvals->all_lambda[efptFEP][h];
500 fepvals->n_lambda = 0;
501 fepvals->all_lambda = NULL;
502 if (fepvals->init_lambda >= 0)
504 fepvals->separate_dvdl[efptFEP] = TRUE;
507 if (file_version >= 13)
509 gmx_fio_do_real(fio, fepvals->sc_alpha);
513 fepvals->sc_alpha = 0;
515 if (file_version >= 38)
517 gmx_fio_do_int(fio, fepvals->sc_power);
521 fepvals->sc_power = 2;
523 if (file_version >= 79)
525 gmx_fio_do_real(fio, fepvals->sc_r_power);
529 fepvals->sc_r_power = 6.0;
531 if (file_version >= 15)
533 gmx_fio_do_real(fio, fepvals->sc_sigma);
537 fepvals->sc_sigma = 0.3;
541 if (file_version >= 71)
543 fepvals->sc_sigma_min = fepvals->sc_sigma;
547 fepvals->sc_sigma_min = 0;
550 if (file_version >= 79)
552 gmx_fio_do_int(fio, fepvals->bScCoul);
556 fepvals->bScCoul = TRUE;
558 if (file_version >= 64)
560 gmx_fio_do_int(fio, fepvals->nstdhdl);
564 fepvals->nstdhdl = 1;
567 if (file_version >= 73)
569 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
570 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
574 fepvals->separate_dhdl_file = esepdhdlfileYES;
575 fepvals->dhdl_derivatives = edhdlderivativesYES;
577 if (file_version >= 71)
579 gmx_fio_do_int(fio, fepvals->dh_hist_size);
580 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
584 fepvals->dh_hist_size = 0;
585 fepvals->dh_hist_spacing = 0.1;
587 if (file_version >= 79)
589 gmx_fio_do_int(fio, fepvals->bPrintEnergy);
593 fepvals->bPrintEnergy = FALSE;
596 /* handle lambda_neighbors */
597 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
599 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
600 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
601 (fepvals->init_lambda < 0) )
603 fepvals->lambda_start_n = (fepvals->init_fep_state -
604 fepvals->lambda_neighbors);
605 fepvals->lambda_stop_n = (fepvals->init_fep_state +
606 fepvals->lambda_neighbors + 1);
607 if (fepvals->lambda_start_n < 0)
609 fepvals->lambda_start_n = 0;;
611 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
613 fepvals->lambda_stop_n = fepvals->n_lambda;
618 fepvals->lambda_start_n = 0;
619 fepvals->lambda_stop_n = fepvals->n_lambda;
624 fepvals->lambda_start_n = 0;
625 fepvals->lambda_stop_n = fepvals->n_lambda;
629 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
633 if (file_version >= 95)
635 gmx_fio_do_int(fio, pull->ngroup);
637 gmx_fio_do_int(fio, pull->ncoord);
638 if (file_version < 95)
640 pull->ngroup = pull->ncoord + 1;
642 gmx_fio_do_int(fio, pull->eGeom);
643 gmx_fio_do_ivec(fio, pull->dim);
644 gmx_fio_do_real(fio, pull->cyl_r1);
645 gmx_fio_do_real(fio, pull->cyl_r0);
646 gmx_fio_do_real(fio, pull->constr_tol);
647 if (file_version >= 95)
649 gmx_fio_do_int(fio, pull->bPrintRef);
651 gmx_fio_do_int(fio, pull->nstxout);
652 gmx_fio_do_int(fio, pull->nstfout);
655 snew(pull->group, pull->ngroup);
656 snew(pull->coord, pull->ncoord);
658 if (file_version < 95)
660 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
661 if (pull->eGeom == epullgDIRPBC)
663 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
665 if (pull->eGeom > epullgDIRPBC)
670 for (g = 0; g < pull->ngroup; g++)
672 /* We read and ignore a pull coordinate for group 0 */
673 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[max(g-1, 0)],
674 bRead, file_version);
677 pull->coord[g-1].group[0] = 0;
678 pull->coord[g-1].group[1] = g;
682 pull->bPrintRef = (pull->group[0].nat > 0);
686 for (g = 0; g < pull->ngroup; g++)
688 do_pull_group(fio, &pull->group[g], bRead);
690 for (g = 0; g < pull->ncoord; g++)
692 do_pull_coord(fio, &pull->coord[g]);
698 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
702 gmx_fio_do_int(fio, rotg->eType);
703 gmx_fio_do_int(fio, rotg->bMassW);
704 gmx_fio_do_int(fio, rotg->nat);
707 snew(rotg->ind, rotg->nat);
709 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
712 snew(rotg->x_ref, rotg->nat);
714 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
715 gmx_fio_do_rvec(fio, rotg->vec);
716 gmx_fio_do_rvec(fio, rotg->pivot);
717 gmx_fio_do_real(fio, rotg->rate);
718 gmx_fio_do_real(fio, rotg->k);
719 gmx_fio_do_real(fio, rotg->slab_dist);
720 gmx_fio_do_real(fio, rotg->min_gaussian);
721 gmx_fio_do_real(fio, rotg->eps);
722 gmx_fio_do_int(fio, rotg->eFittype);
723 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
724 gmx_fio_do_real(fio, rotg->PotAngle_step);
727 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
731 gmx_fio_do_int(fio, rot->ngrp);
732 gmx_fio_do_int(fio, rot->nstrout);
733 gmx_fio_do_int(fio, rot->nstsout);
736 snew(rot->grp, rot->ngrp);
738 for (g = 0; g < rot->ngrp; g++)
740 do_rotgrp(fio, &rot->grp[g], bRead);
745 static void do_swapcoords(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead)
750 gmx_fio_do_int(fio, swap->nat);
751 gmx_fio_do_int(fio, swap->nat_sol);
752 for (j = 0; j < 2; j++)
754 gmx_fio_do_int(fio, swap->nat_split[j]);
755 gmx_fio_do_int(fio, swap->massw_split[j]);
757 gmx_fio_do_int(fio, swap->nstswap);
758 gmx_fio_do_int(fio, swap->nAverage);
759 gmx_fio_do_real(fio, swap->threshold);
760 gmx_fio_do_real(fio, swap->cyl0r);
761 gmx_fio_do_real(fio, swap->cyl0u);
762 gmx_fio_do_real(fio, swap->cyl0l);
763 gmx_fio_do_real(fio, swap->cyl1r);
764 gmx_fio_do_real(fio, swap->cyl1u);
765 gmx_fio_do_real(fio, swap->cyl1l);
769 snew(swap->ind, swap->nat);
770 snew(swap->ind_sol, swap->nat_sol);
771 for (j = 0; j < 2; j++)
773 snew(swap->ind_split[j], swap->nat_split[j]);
777 gmx_fio_ndo_int(fio, swap->ind, swap->nat);
778 gmx_fio_ndo_int(fio, swap->ind_sol, swap->nat_sol);
779 for (j = 0; j < 2; j++)
781 gmx_fio_ndo_int(fio, swap->ind_split[j], swap->nat_split[j]);
784 for (j = 0; j < eCompNR; j++)
786 gmx_fio_do_int(fio, swap->nanions[j]);
787 gmx_fio_do_int(fio, swap->ncations[j]);
793 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
794 int file_version, real *fudgeQQ)
796 int i, j, k, *tmp, idum = 0;
800 real zerotemptime, finish_t, init_temp, finish_temp;
802 if (file_version != tpx_version)
804 /* Give a warning about features that are not accessible */
805 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
806 file_version, tpx_version);
814 if (file_version == 0)
819 /* Basic inputrec stuff */
820 gmx_fio_do_int(fio, ir->eI);
821 if (file_version >= 62)
823 gmx_fio_do_int64(fio, ir->nsteps);
827 gmx_fio_do_int(fio, idum);
830 if (file_version > 25)
832 if (file_version >= 62)
834 gmx_fio_do_int64(fio, ir->init_step);
838 gmx_fio_do_int(fio, idum);
839 ir->init_step = idum;
847 if (file_version >= 58)
849 gmx_fio_do_int(fio, ir->simulation_part);
853 ir->simulation_part = 1;
856 if (file_version >= 67)
858 gmx_fio_do_int(fio, ir->nstcalcenergy);
862 ir->nstcalcenergy = 1;
864 if (file_version < 53)
866 /* The pbc info has been moved out of do_inputrec,
867 * since we always want it, also without reading the inputrec.
869 gmx_fio_do_int(fio, ir->ePBC);
870 if ((file_version <= 15) && (ir->ePBC == 2))
874 if (file_version >= 45)
876 gmx_fio_do_int(fio, ir->bPeriodicMols);
883 ir->bPeriodicMols = TRUE;
887 ir->bPeriodicMols = FALSE;
891 if (file_version >= 81)
893 gmx_fio_do_int(fio, ir->cutoff_scheme);
894 if (file_version < 94)
896 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
901 ir->cutoff_scheme = ecutsGROUP;
903 gmx_fio_do_int(fio, ir->ns_type);
904 gmx_fio_do_int(fio, ir->nstlist);
905 gmx_fio_do_int(fio, ir->ndelta);
906 if (file_version < 41)
908 gmx_fio_do_int(fio, idum);
909 gmx_fio_do_int(fio, idum);
911 if (file_version >= 45)
913 gmx_fio_do_real(fio, ir->rtpi);
919 gmx_fio_do_int(fio, ir->nstcomm);
920 if (file_version > 34)
922 gmx_fio_do_int(fio, ir->comm_mode);
924 else if (ir->nstcomm < 0)
926 ir->comm_mode = ecmANGULAR;
930 ir->comm_mode = ecmLINEAR;
932 ir->nstcomm = abs(ir->nstcomm);
934 if (file_version > 25)
936 gmx_fio_do_int(fio, ir->nstcheckpoint);
940 ir->nstcheckpoint = 0;
943 gmx_fio_do_int(fio, ir->nstcgsteep);
945 if (file_version >= 30)
947 gmx_fio_do_int(fio, ir->nbfgscorr);
954 gmx_fio_do_int(fio, ir->nstlog);
955 gmx_fio_do_int(fio, ir->nstxout);
956 gmx_fio_do_int(fio, ir->nstvout);
957 gmx_fio_do_int(fio, ir->nstfout);
958 gmx_fio_do_int(fio, ir->nstenergy);
959 gmx_fio_do_int(fio, ir->nstxout_compressed);
960 if (file_version >= 59)
962 gmx_fio_do_double(fio, ir->init_t);
963 gmx_fio_do_double(fio, ir->delta_t);
967 gmx_fio_do_real(fio, rdum);
969 gmx_fio_do_real(fio, rdum);
972 gmx_fio_do_real(fio, ir->x_compression_precision);
973 if (file_version < 19)
975 gmx_fio_do_int(fio, idum);
976 gmx_fio_do_int(fio, idum);
978 if (file_version < 18)
980 gmx_fio_do_int(fio, idum);
982 if (file_version >= 81)
984 gmx_fio_do_real(fio, ir->verletbuf_tol);
988 ir->verletbuf_tol = 0;
990 gmx_fio_do_real(fio, ir->rlist);
991 if (file_version >= 67)
993 gmx_fio_do_real(fio, ir->rlistlong);
995 if (file_version >= 82 && file_version != 90)
997 gmx_fio_do_int(fio, ir->nstcalclr);
1001 /* Calculate at NS steps */
1002 ir->nstcalclr = ir->nstlist;
1004 gmx_fio_do_int(fio, ir->coulombtype);
1005 if (file_version < 32 && ir->coulombtype == eelRF)
1007 ir->coulombtype = eelRF_NEC;
1009 if (file_version >= 81)
1011 gmx_fio_do_int(fio, ir->coulomb_modifier);
1015 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1017 gmx_fio_do_real(fio, ir->rcoulomb_switch);
1018 gmx_fio_do_real(fio, ir->rcoulomb);
1019 gmx_fio_do_int(fio, ir->vdwtype);
1020 if (file_version >= 81)
1022 gmx_fio_do_int(fio, ir->vdw_modifier);
1026 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1028 gmx_fio_do_real(fio, ir->rvdw_switch);
1029 gmx_fio_do_real(fio, ir->rvdw);
1030 if (file_version < 67)
1032 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1034 gmx_fio_do_int(fio, ir->eDispCorr);
1035 gmx_fio_do_real(fio, ir->epsilon_r);
1036 if (file_version >= 37)
1038 gmx_fio_do_real(fio, ir->epsilon_rf);
1042 if (EEL_RF(ir->coulombtype))
1044 ir->epsilon_rf = ir->epsilon_r;
1045 ir->epsilon_r = 1.0;
1049 ir->epsilon_rf = 1.0;
1052 if (file_version >= 29)
1054 gmx_fio_do_real(fio, ir->tabext);
1061 if (file_version > 25)
1063 gmx_fio_do_int(fio, ir->gb_algorithm);
1064 gmx_fio_do_int(fio, ir->nstgbradii);
1065 gmx_fio_do_real(fio, ir->rgbradii);
1066 gmx_fio_do_real(fio, ir->gb_saltconc);
1067 gmx_fio_do_int(fio, ir->implicit_solvent);
1071 ir->gb_algorithm = egbSTILL;
1074 ir->gb_saltconc = 0;
1075 ir->implicit_solvent = eisNO;
1077 if (file_version >= 55)
1079 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1080 gmx_fio_do_real(fio, ir->gb_obc_alpha);
1081 gmx_fio_do_real(fio, ir->gb_obc_beta);
1082 gmx_fio_do_real(fio, ir->gb_obc_gamma);
1083 if (file_version >= 60)
1085 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1086 gmx_fio_do_int(fio, ir->sa_algorithm);
1090 ir->gb_dielectric_offset = 0.009;
1091 ir->sa_algorithm = esaAPPROX;
1093 gmx_fio_do_real(fio, ir->sa_surface_tension);
1095 /* Override sa_surface_tension if it is not changed in the mpd-file */
1096 if (ir->sa_surface_tension < 0)
1098 if (ir->gb_algorithm == egbSTILL)
1100 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1102 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1104 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1111 /* Better use sensible values than insane (0.0) ones... */
1112 ir->gb_epsilon_solvent = 80;
1113 ir->gb_obc_alpha = 1.0;
1114 ir->gb_obc_beta = 0.8;
1115 ir->gb_obc_gamma = 4.85;
1116 ir->sa_surface_tension = 2.092;
1120 if (file_version >= 81)
1122 gmx_fio_do_real(fio, ir->fourier_spacing);
1126 ir->fourier_spacing = 0.0;
1128 gmx_fio_do_int(fio, ir->nkx);
1129 gmx_fio_do_int(fio, ir->nky);
1130 gmx_fio_do_int(fio, ir->nkz);
1131 gmx_fio_do_int(fio, ir->pme_order);
1132 gmx_fio_do_real(fio, ir->ewald_rtol);
1134 if (file_version >= 93)
1136 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1140 ir->ewald_rtol_lj = ir->ewald_rtol;
1143 if (file_version >= 24)
1145 gmx_fio_do_int(fio, ir->ewald_geometry);
1149 ir->ewald_geometry = eewg3D;
1152 if (file_version <= 17)
1154 ir->epsilon_surface = 0;
1155 if (file_version == 17)
1157 gmx_fio_do_int(fio, idum);
1162 gmx_fio_do_real(fio, ir->epsilon_surface);
1165 gmx_fio_do_gmx_bool(fio, ir->bOptFFT);
1167 if (file_version >= 93)
1169 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1171 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1172 gmx_fio_do_int(fio, ir->etc);
1173 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1174 * but the values 0 and 1 still mean no and
1175 * berendsen temperature coupling, respectively.
1177 if (file_version >= 79)
1179 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1181 if (file_version >= 71)
1183 gmx_fio_do_int(fio, ir->nsttcouple);
1187 ir->nsttcouple = ir->nstcalcenergy;
1189 if (file_version <= 15)
1191 gmx_fio_do_int(fio, idum);
1193 if (file_version <= 17)
1195 gmx_fio_do_int(fio, ir->epct);
1196 if (file_version <= 15)
1200 ir->epct = epctSURFACETENSION;
1202 gmx_fio_do_int(fio, idum);
1205 /* we have removed the NO alternative at the beginning */
1209 ir->epct = epctISOTROPIC;
1213 ir->epc = epcBERENDSEN;
1218 gmx_fio_do_int(fio, ir->epc);
1219 gmx_fio_do_int(fio, ir->epct);
1221 if (file_version >= 71)
1223 gmx_fio_do_int(fio, ir->nstpcouple);
1227 ir->nstpcouple = ir->nstcalcenergy;
1229 gmx_fio_do_real(fio, ir->tau_p);
1230 if (file_version <= 15)
1232 gmx_fio_do_rvec(fio, vdum);
1233 clear_mat(ir->ref_p);
1234 for (i = 0; i < DIM; i++)
1236 ir->ref_p[i][i] = vdum[i];
1241 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1242 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1243 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1245 if (file_version <= 15)
1247 gmx_fio_do_rvec(fio, vdum);
1248 clear_mat(ir->compress);
1249 for (i = 0; i < DIM; i++)
1251 ir->compress[i][i] = vdum[i];
1256 gmx_fio_do_rvec(fio, ir->compress[XX]);
1257 gmx_fio_do_rvec(fio, ir->compress[YY]);
1258 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1260 if (file_version >= 47)
1262 gmx_fio_do_int(fio, ir->refcoord_scaling);
1263 gmx_fio_do_rvec(fio, ir->posres_com);
1264 gmx_fio_do_rvec(fio, ir->posres_comB);
1268 ir->refcoord_scaling = erscNO;
1269 clear_rvec(ir->posres_com);
1270 clear_rvec(ir->posres_comB);
1272 if ((file_version > 25) && (file_version < 79))
1274 gmx_fio_do_int(fio, ir->andersen_seed);
1278 ir->andersen_seed = 0;
1280 if (file_version < 26)
1282 gmx_fio_do_gmx_bool(fio, bSimAnn);
1283 gmx_fio_do_real(fio, zerotemptime);
1286 if (file_version < 37)
1288 gmx_fio_do_real(fio, rdum);
1291 gmx_fio_do_real(fio, ir->shake_tol);
1292 if (file_version < 54)
1294 gmx_fio_do_real(fio, *fudgeQQ);
1297 gmx_fio_do_int(fio, ir->efep);
1298 if (file_version <= 14 && ir->efep != efepNO)
1302 do_fepvals(fio, ir->fepvals, bRead, file_version);
1304 if (file_version >= 79)
1306 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1309 ir->bSimTemp = TRUE;
1314 ir->bSimTemp = FALSE;
1318 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1321 if (file_version >= 79)
1323 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1326 ir->bExpanded = TRUE;
1330 ir->bExpanded = FALSE;
1335 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1337 if (file_version >= 57)
1339 gmx_fio_do_int(fio, ir->eDisre);
1341 gmx_fio_do_int(fio, ir->eDisreWeighting);
1342 if (file_version < 22)
1344 if (ir->eDisreWeighting == 0)
1346 ir->eDisreWeighting = edrwEqual;
1350 ir->eDisreWeighting = edrwConservative;
1353 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1354 gmx_fio_do_real(fio, ir->dr_fc);
1355 gmx_fio_do_real(fio, ir->dr_tau);
1356 gmx_fio_do_int(fio, ir->nstdisreout);
1357 if (file_version >= 22)
1359 gmx_fio_do_real(fio, ir->orires_fc);
1360 gmx_fio_do_real(fio, ir->orires_tau);
1361 gmx_fio_do_int(fio, ir->nstorireout);
1367 ir->nstorireout = 0;
1369 if (file_version >= 26 && file_version < 79)
1371 gmx_fio_do_real(fio, ir->dihre_fc);
1372 if (file_version < 56)
1374 gmx_fio_do_real(fio, rdum);
1375 gmx_fio_do_int(fio, idum);
1383 gmx_fio_do_real(fio, ir->em_stepsize);
1384 gmx_fio_do_real(fio, ir->em_tol);
1385 if (file_version >= 22)
1387 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1391 ir->bShakeSOR = TRUE;
1393 if (file_version >= 11)
1395 gmx_fio_do_int(fio, ir->niter);
1400 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1403 if (file_version >= 21)
1405 gmx_fio_do_real(fio, ir->fc_stepsize);
1409 ir->fc_stepsize = 0;
1411 gmx_fio_do_int(fio, ir->eConstrAlg);
1412 gmx_fio_do_int(fio, ir->nProjOrder);
1413 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1414 if (file_version <= 14)
1416 gmx_fio_do_int(fio, idum);
1418 if (file_version >= 26)
1420 gmx_fio_do_int(fio, ir->nLincsIter);
1425 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1428 if (file_version < 33)
1430 gmx_fio_do_real(fio, bd_temp);
1432 gmx_fio_do_real(fio, ir->bd_fric);
1433 if (file_version >= tpxv_Use64BitRandomSeed)
1435 gmx_fio_do_int64(fio, ir->ld_seed);
1439 gmx_fio_do_int(fio, idum);
1442 if (file_version >= 33)
1444 for (i = 0; i < DIM; i++)
1446 gmx_fio_do_rvec(fio, ir->deform[i]);
1451 for (i = 0; i < DIM; i++)
1453 clear_rvec(ir->deform[i]);
1456 if (file_version >= 14)
1458 gmx_fio_do_real(fio, ir->cos_accel);
1464 gmx_fio_do_int(fio, ir->userint1);
1465 gmx_fio_do_int(fio, ir->userint2);
1466 gmx_fio_do_int(fio, ir->userint3);
1467 gmx_fio_do_int(fio, ir->userint4);
1468 gmx_fio_do_real(fio, ir->userreal1);
1469 gmx_fio_do_real(fio, ir->userreal2);
1470 gmx_fio_do_real(fio, ir->userreal3);
1471 gmx_fio_do_real(fio, ir->userreal4);
1474 if (file_version >= 77)
1476 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1481 snew(ir->adress, 1);
1483 gmx_fio_do_int(fio, ir->adress->type);
1484 gmx_fio_do_real(fio, ir->adress->const_wf);
1485 gmx_fio_do_real(fio, ir->adress->ex_width);
1486 gmx_fio_do_real(fio, ir->adress->hy_width);
1487 gmx_fio_do_int(fio, ir->adress->icor);
1488 gmx_fio_do_int(fio, ir->adress->site);
1489 gmx_fio_do_rvec(fio, ir->adress->refs);
1490 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1491 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1492 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1493 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1497 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1499 if (ir->adress->n_tf_grps > 0)
1501 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1505 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1507 if (ir->adress->n_energy_grps > 0)
1509 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1515 ir->bAdress = FALSE;
1519 if (file_version >= 48)
1521 gmx_fio_do_int(fio, ir->ePull);
1522 if (ir->ePull != epullNO)
1528 do_pull(fio, ir->pull, bRead, file_version);
1533 ir->ePull = epullNO;
1536 /* Enforced rotation */
1537 if (file_version >= 74)
1539 gmx_fio_do_int(fio, ir->bRot);
1540 if (ir->bRot == TRUE)
1546 do_rot(fio, ir->rot, bRead);
1555 gmx_fio_do_int(fio, ir->opts.ngtc);
1556 if (file_version >= 69)
1558 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1562 ir->opts.nhchainlength = 1;
1564 gmx_fio_do_int(fio, ir->opts.ngacc);
1565 gmx_fio_do_int(fio, ir->opts.ngfrz);
1566 gmx_fio_do_int(fio, ir->opts.ngener);
1570 snew(ir->opts.nrdf, ir->opts.ngtc);
1571 snew(ir->opts.ref_t, ir->opts.ngtc);
1572 snew(ir->opts.annealing, ir->opts.ngtc);
1573 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1574 snew(ir->opts.anneal_time, ir->opts.ngtc);
1575 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1576 snew(ir->opts.tau_t, ir->opts.ngtc);
1577 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1578 snew(ir->opts.acc, ir->opts.ngacc);
1579 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1581 if (ir->opts.ngtc > 0)
1583 if (bRead && file_version < 13)
1585 snew(tmp, ir->opts.ngtc);
1586 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1587 for (i = 0; i < ir->opts.ngtc; i++)
1589 ir->opts.nrdf[i] = tmp[i];
1595 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1597 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1598 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1599 if (file_version < 33 && ir->eI == eiBD)
1601 for (i = 0; i < ir->opts.ngtc; i++)
1603 ir->opts.tau_t[i] = bd_temp;
1607 if (ir->opts.ngfrz > 0)
1609 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1611 if (ir->opts.ngacc > 0)
1613 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1615 if (file_version >= 12)
1617 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1618 ir->opts.ngener*ir->opts.ngener);
1621 if (bRead && file_version < 26)
1623 for (i = 0; i < ir->opts.ngtc; i++)
1627 ir->opts.annealing[i] = eannSINGLE;
1628 ir->opts.anneal_npoints[i] = 2;
1629 snew(ir->opts.anneal_time[i], 2);
1630 snew(ir->opts.anneal_temp[i], 2);
1631 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1632 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1633 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1634 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1635 ir->opts.anneal_time[i][0] = ir->init_t;
1636 ir->opts.anneal_time[i][1] = finish_t;
1637 ir->opts.anneal_temp[i][0] = init_temp;
1638 ir->opts.anneal_temp[i][1] = finish_temp;
1642 ir->opts.annealing[i] = eannNO;
1643 ir->opts.anneal_npoints[i] = 0;
1649 /* file version 26 or later */
1650 /* First read the lists with annealing and npoints for each group */
1651 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1652 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1653 for (j = 0; j < (ir->opts.ngtc); j++)
1655 k = ir->opts.anneal_npoints[j];
1658 snew(ir->opts.anneal_time[j], k);
1659 snew(ir->opts.anneal_temp[j], k);
1661 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1662 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1666 if (file_version >= 45)
1668 gmx_fio_do_int(fio, ir->nwall);
1669 gmx_fio_do_int(fio, ir->wall_type);
1670 if (file_version >= 50)
1672 gmx_fio_do_real(fio, ir->wall_r_linpot);
1676 ir->wall_r_linpot = -1;
1678 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1679 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1680 gmx_fio_do_real(fio, ir->wall_density[0]);
1681 gmx_fio_do_real(fio, ir->wall_density[1]);
1682 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1688 ir->wall_atomtype[0] = -1;
1689 ir->wall_atomtype[1] = -1;
1690 ir->wall_density[0] = 0;
1691 ir->wall_density[1] = 0;
1692 ir->wall_ewald_zfac = 3;
1694 /* Cosine stuff for electric fields */
1695 for (j = 0; (j < DIM); j++)
1697 gmx_fio_do_int(fio, ir->ex[j].n);
1698 gmx_fio_do_int(fio, ir->et[j].n);
1701 snew(ir->ex[j].a, ir->ex[j].n);
1702 snew(ir->ex[j].phi, ir->ex[j].n);
1703 snew(ir->et[j].a, ir->et[j].n);
1704 snew(ir->et[j].phi, ir->et[j].n);
1706 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1707 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1708 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1709 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1713 if (file_version >= tpxv_ComputationalElectrophysiology)
1715 gmx_fio_do_int(fio, ir->eSwapCoords);
1716 if (ir->eSwapCoords != eswapNO)
1722 do_swapcoords(fio, ir->swap, bRead);
1727 if (file_version >= 39)
1729 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1730 gmx_fio_do_int(fio, ir->QMMMscheme);
1731 gmx_fio_do_real(fio, ir->scalefactor);
1732 gmx_fio_do_int(fio, ir->opts.ngQM);
1735 snew(ir->opts.QMmethod, ir->opts.ngQM);
1736 snew(ir->opts.QMbasis, ir->opts.ngQM);
1737 snew(ir->opts.QMcharge, ir->opts.ngQM);
1738 snew(ir->opts.QMmult, ir->opts.ngQM);
1739 snew(ir->opts.bSH, ir->opts.ngQM);
1740 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1741 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1742 snew(ir->opts.SAon, ir->opts.ngQM);
1743 snew(ir->opts.SAoff, ir->opts.ngQM);
1744 snew(ir->opts.SAsteps, ir->opts.ngQM);
1745 snew(ir->opts.bOPT, ir->opts.ngQM);
1746 snew(ir->opts.bTS, ir->opts.ngQM);
1748 if (ir->opts.ngQM > 0)
1750 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1751 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1752 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1753 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1754 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1755 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1756 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1757 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1758 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1759 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1760 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1761 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1763 /* end of QMMM stuff */
1768 static void do_harm(t_fileio *fio, t_iparams *iparams)
1770 gmx_fio_do_real(fio, iparams->harmonic.rA);
1771 gmx_fio_do_real(fio, iparams->harmonic.krA);
1772 gmx_fio_do_real(fio, iparams->harmonic.rB);
1773 gmx_fio_do_real(fio, iparams->harmonic.krB);
1776 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1777 gmx_bool bRead, int file_version)
1784 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1794 do_harm(fio, iparams);
1795 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1797 /* Correct incorrect storage of parameters */
1798 iparams->pdihs.phiB = iparams->pdihs.phiA;
1799 iparams->pdihs.cpB = iparams->pdihs.cpA;
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);
1813 gmx_fio_do_real(fio, iparams->restraint.lowA);
1814 gmx_fio_do_real(fio, iparams->restraint.up1A);
1815 gmx_fio_do_real(fio, iparams->restraint.up2A);
1816 gmx_fio_do_real(fio, iparams->restraint.kA);
1817 gmx_fio_do_real(fio, iparams->restraint.lowB);
1818 gmx_fio_do_real(fio, iparams->restraint.up1B);
1819 gmx_fio_do_real(fio, iparams->restraint.up2B);
1820 gmx_fio_do_real(fio, iparams->restraint.kB);
1826 gmx_fio_do_real(fio, iparams->tab.kA);
1827 gmx_fio_do_int(fio, iparams->tab.table);
1828 gmx_fio_do_real(fio, iparams->tab.kB);
1830 case F_CROSS_BOND_BONDS:
1831 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1832 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1833 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1835 case F_CROSS_BOND_ANGLES:
1836 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1837 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1838 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1839 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1841 case F_UREY_BRADLEY:
1842 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1843 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1844 gmx_fio_do_real(fio, iparams->u_b.r13A);
1845 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1846 if (file_version >= 79)
1848 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1849 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1850 gmx_fio_do_real(fio, iparams->u_b.r13B);
1851 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1855 iparams->u_b.thetaB = iparams->u_b.thetaA;
1856 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1857 iparams->u_b.r13B = iparams->u_b.r13A;
1858 iparams->u_b.kUBB = iparams->u_b.kUBA;
1861 case F_QUARTIC_ANGLES:
1862 gmx_fio_do_real(fio, iparams->qangle.theta);
1863 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1866 gmx_fio_do_real(fio, iparams->bham.a);
1867 gmx_fio_do_real(fio, iparams->bham.b);
1868 gmx_fio_do_real(fio, iparams->bham.c);
1871 gmx_fio_do_real(fio, iparams->morse.b0A);
1872 gmx_fio_do_real(fio, iparams->morse.cbA);
1873 gmx_fio_do_real(fio, iparams->morse.betaA);
1874 if (file_version >= 79)
1876 gmx_fio_do_real(fio, iparams->morse.b0B);
1877 gmx_fio_do_real(fio, iparams->morse.cbB);
1878 gmx_fio_do_real(fio, iparams->morse.betaB);
1882 iparams->morse.b0B = iparams->morse.b0A;
1883 iparams->morse.cbB = iparams->morse.cbA;
1884 iparams->morse.betaB = iparams->morse.betaA;
1888 gmx_fio_do_real(fio, iparams->cubic.b0);
1889 gmx_fio_do_real(fio, iparams->cubic.kb);
1890 gmx_fio_do_real(fio, iparams->cubic.kcub);
1894 case F_POLARIZATION:
1895 gmx_fio_do_real(fio, iparams->polarize.alpha);
1898 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1899 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1900 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1903 if (file_version < 31)
1905 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1907 gmx_fio_do_real(fio, iparams->wpol.al_x);
1908 gmx_fio_do_real(fio, iparams->wpol.al_y);
1909 gmx_fio_do_real(fio, iparams->wpol.al_z);
1910 gmx_fio_do_real(fio, iparams->wpol.rOH);
1911 gmx_fio_do_real(fio, iparams->wpol.rHH);
1912 gmx_fio_do_real(fio, iparams->wpol.rOD);
1915 gmx_fio_do_real(fio, iparams->thole.a);
1916 gmx_fio_do_real(fio, iparams->thole.alpha1);
1917 gmx_fio_do_real(fio, iparams->thole.alpha2);
1918 gmx_fio_do_real(fio, iparams->thole.rfac);
1921 gmx_fio_do_real(fio, iparams->lj.c6);
1922 gmx_fio_do_real(fio, iparams->lj.c12);
1925 gmx_fio_do_real(fio, iparams->lj14.c6A);
1926 gmx_fio_do_real(fio, iparams->lj14.c12A);
1927 gmx_fio_do_real(fio, iparams->lj14.c6B);
1928 gmx_fio_do_real(fio, iparams->lj14.c12B);
1931 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1932 gmx_fio_do_real(fio, iparams->ljc14.qi);
1933 gmx_fio_do_real(fio, iparams->ljc14.qj);
1934 gmx_fio_do_real(fio, iparams->ljc14.c6);
1935 gmx_fio_do_real(fio, iparams->ljc14.c12);
1937 case F_LJC_PAIRS_NB:
1938 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1939 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1940 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1941 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1947 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1948 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1949 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1951 /* Read the incorrectly stored multiplicity */
1952 gmx_fio_do_real(fio, iparams->harmonic.rB);
1953 gmx_fio_do_real(fio, iparams->harmonic.krB);
1954 iparams->pdihs.phiB = iparams->pdihs.phiA;
1955 iparams->pdihs.cpB = iparams->pdihs.cpA;
1959 gmx_fio_do_real(fio, iparams->pdihs.phiB);
1960 gmx_fio_do_real(fio, iparams->pdihs.cpB);
1961 gmx_fio_do_int(fio, iparams->pdihs.mult);
1965 gmx_fio_do_int(fio, iparams->disres.label);
1966 gmx_fio_do_int(fio, iparams->disres.type);
1967 gmx_fio_do_real(fio, iparams->disres.low);
1968 gmx_fio_do_real(fio, iparams->disres.up1);
1969 gmx_fio_do_real(fio, iparams->disres.up2);
1970 gmx_fio_do_real(fio, iparams->disres.kfac);
1973 gmx_fio_do_int(fio, iparams->orires.ex);
1974 gmx_fio_do_int(fio, iparams->orires.label);
1975 gmx_fio_do_int(fio, iparams->orires.power);
1976 gmx_fio_do_real(fio, iparams->orires.c);
1977 gmx_fio_do_real(fio, iparams->orires.obs);
1978 gmx_fio_do_real(fio, iparams->orires.kfac);
1981 if (file_version < 82)
1983 gmx_fio_do_int(fio, idum);
1984 gmx_fio_do_int(fio, idum);
1986 gmx_fio_do_real(fio, iparams->dihres.phiA);
1987 gmx_fio_do_real(fio, iparams->dihres.dphiA);
1988 gmx_fio_do_real(fio, iparams->dihres.kfacA);
1989 if (file_version >= 82)
1991 gmx_fio_do_real(fio, iparams->dihres.phiB);
1992 gmx_fio_do_real(fio, iparams->dihres.dphiB);
1993 gmx_fio_do_real(fio, iparams->dihres.kfacB);
1997 iparams->dihres.phiB = iparams->dihres.phiA;
1998 iparams->dihres.dphiB = iparams->dihres.dphiA;
1999 iparams->dihres.kfacB = iparams->dihres.kfacA;
2003 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
2004 gmx_fio_do_rvec(fio, iparams->posres.fcA);
2005 if (bRead && file_version < 27)
2007 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
2008 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
2012 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
2013 gmx_fio_do_rvec(fio, iparams->posres.fcB);
2017 gmx_fio_do_int(fio, iparams->fbposres.geom);
2018 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
2019 gmx_fio_do_real(fio, iparams->fbposres.r);
2020 gmx_fio_do_real(fio, iparams->fbposres.k);
2023 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2024 if (file_version >= 25)
2026 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2030 /* Fourier dihedrals are internally represented
2031 * as Ryckaert-Bellemans since those are faster to compute.
2033 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2034 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2038 gmx_fio_do_real(fio, iparams->constr.dA);
2039 gmx_fio_do_real(fio, iparams->constr.dB);
2042 gmx_fio_do_real(fio, iparams->settle.doh);
2043 gmx_fio_do_real(fio, iparams->settle.dhh);
2046 gmx_fio_do_real(fio, iparams->vsite.a);
2051 gmx_fio_do_real(fio, iparams->vsite.a);
2052 gmx_fio_do_real(fio, iparams->vsite.b);
2057 gmx_fio_do_real(fio, iparams->vsite.a);
2058 gmx_fio_do_real(fio, iparams->vsite.b);
2059 gmx_fio_do_real(fio, iparams->vsite.c);
2062 gmx_fio_do_int(fio, iparams->vsiten.n);
2063 gmx_fio_do_real(fio, iparams->vsiten.a);
2068 /* We got rid of some parameters in version 68 */
2069 if (bRead && file_version < 68)
2071 gmx_fio_do_real(fio, rdum);
2072 gmx_fio_do_real(fio, rdum);
2073 gmx_fio_do_real(fio, rdum);
2074 gmx_fio_do_real(fio, rdum);
2076 gmx_fio_do_real(fio, iparams->gb.sar);
2077 gmx_fio_do_real(fio, iparams->gb.st);
2078 gmx_fio_do_real(fio, iparams->gb.pi);
2079 gmx_fio_do_real(fio, iparams->gb.gbr);
2080 gmx_fio_do_real(fio, iparams->gb.bmlt);
2083 gmx_fio_do_int(fio, iparams->cmap.cmapA);
2084 gmx_fio_do_int(fio, iparams->cmap.cmapB);
2087 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
2088 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
2092 gmx_fio_unset_comment(fio);
2096 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
2103 gmx_fio_set_comment(fio, interaction_function[ftype].name);
2105 if (file_version < 44)
2107 for (i = 0; i < MAXNODES; i++)
2109 gmx_fio_do_int(fio, idum);
2112 gmx_fio_do_int(fio, ilist->nr);
2115 snew(ilist->iatoms, ilist->nr);
2117 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2120 gmx_fio_unset_comment(fio);
2124 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2125 gmx_bool bRead, int file_version)
2130 gmx_fio_do_int(fio, ffparams->atnr);
2131 if (file_version < 57)
2133 gmx_fio_do_int(fio, idum);
2135 gmx_fio_do_int(fio, ffparams->ntypes);
2138 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2139 ffparams->atnr, ffparams->ntypes);
2143 snew(ffparams->functype, ffparams->ntypes);
2144 snew(ffparams->iparams, ffparams->ntypes);
2146 /* Read/write all the function types */
2147 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2150 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2153 if (file_version >= 66)
2155 gmx_fio_do_double(fio, ffparams->reppow);
2159 ffparams->reppow = 12.0;
2162 if (file_version >= 57)
2164 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2167 /* Check whether all these function types are supported by the code.
2168 * In practice the code is backwards compatible, which means that the
2169 * numbering may have to be altered from old numbering to new numbering
2171 for (i = 0; (i < ffparams->ntypes); i++)
2175 /* Loop over file versions */
2176 for (k = 0; (k < NFTUPD); k++)
2178 /* Compare the read file_version to the update table */
2179 if ((file_version < ftupd[k].fvnr) &&
2180 (ffparams->functype[i] >= ftupd[k].ftype))
2182 ffparams->functype[i] += 1;
2185 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2186 i, ffparams->functype[i],
2187 interaction_function[ftupd[k].ftype].longname);
2194 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2198 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2203 static void add_settle_atoms(t_ilist *ilist)
2207 /* Settle used to only store the first atom: add the other two */
2208 srenew(ilist->iatoms, 2*ilist->nr);
2209 for (i = ilist->nr/2-1; i >= 0; i--)
2211 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2212 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2213 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2214 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2216 ilist->nr = 2*ilist->nr;
2219 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2222 int i, j, renum[F_NRE];
2226 for (j = 0; (j < F_NRE); j++)
2231 for (k = 0; k < NFTUPD; k++)
2233 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2242 ilist[j].iatoms = NULL;
2246 do_ilist(fio, &ilist[j], bRead, file_version, j);
2247 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2249 add_settle_atoms(&ilist[j]);
2253 if (bRead && gmx_debug_at)
2254 pr_ilist(debug,0,interaction_function[j].longname,
2255 functype,&ilist[j],TRUE);
2260 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2261 gmx_bool bRead, int file_version)
2263 do_ffparams(fio, ffparams, bRead, file_version);
2265 if (file_version >= 54)
2267 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2270 do_ilists(fio, molt->ilist, bRead, file_version);
2273 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2275 int i, idum, dum_nra, *dum_a;
2277 if (file_version < 44)
2279 for (i = 0; i < MAXNODES; i++)
2281 gmx_fio_do_int(fio, idum);
2284 gmx_fio_do_int(fio, block->nr);
2285 if (file_version < 51)
2287 gmx_fio_do_int(fio, dum_nra);
2291 if ((block->nalloc_index > 0) && (NULL != block->index))
2293 sfree(block->index);
2295 block->nalloc_index = block->nr+1;
2296 snew(block->index, block->nalloc_index);
2298 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2300 if (file_version < 51 && dum_nra > 0)
2302 snew(dum_a, dum_nra);
2303 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2308 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2313 if (file_version < 44)
2315 for (i = 0; i < MAXNODES; i++)
2317 gmx_fio_do_int(fio, idum);
2320 gmx_fio_do_int(fio, block->nr);
2321 gmx_fio_do_int(fio, block->nra);
2324 block->nalloc_index = block->nr+1;
2325 snew(block->index, block->nalloc_index);
2326 block->nalloc_a = block->nra;
2327 snew(block->a, block->nalloc_a);
2329 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2330 gmx_fio_ndo_int(fio, block->a, block->nra);
2333 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2334 int file_version, gmx_groups_t *groups, int atnr)
2338 gmx_fio_do_real(fio, atom->m);
2339 gmx_fio_do_real(fio, atom->q);
2340 gmx_fio_do_real(fio, atom->mB);
2341 gmx_fio_do_real(fio, atom->qB);
2342 gmx_fio_do_ushort(fio, atom->type);
2343 gmx_fio_do_ushort(fio, atom->typeB);
2344 gmx_fio_do_int(fio, atom->ptype);
2345 gmx_fio_do_int(fio, atom->resind);
2346 if (file_version >= 52)
2348 gmx_fio_do_int(fio, atom->atomnumber);
2352 atom->atomnumber = NOTSET;
2354 if (file_version < 23)
2358 else if (file_version < 39)
2367 if (file_version < 57)
2369 unsigned char uchar[egcNR];
2370 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2371 for (i = myngrp; (i < ngrp); i++)
2375 /* Copy the old data format to the groups struct */
2376 for (i = 0; i < ngrp; i++)
2378 groups->grpnr[i][atnr] = uchar[i];
2383 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2388 if (file_version < 23)
2392 else if (file_version < 39)
2401 for (j = 0; (j < ngrp); j++)
2405 gmx_fio_do_int(fio, grps[j].nr);
2408 snew(grps[j].nm_ind, grps[j].nr);
2410 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2415 snew(grps[j].nm_ind, grps[j].nr);
2420 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2426 gmx_fio_do_int(fio, ls);
2427 *nm = get_symtab_handle(symtab, ls);
2431 ls = lookup_symtab(symtab, *nm);
2432 gmx_fio_do_int(fio, ls);
2436 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2441 for (j = 0; (j < nstr); j++)
2443 do_symstr(fio, &(nm[j]), bRead, symtab);
2447 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2448 t_symtab *symtab, int file_version)
2452 for (j = 0; (j < n); j++)
2454 do_symstr(fio, &(ri[j].name), bRead, symtab);
2455 if (file_version >= 63)
2457 gmx_fio_do_int(fio, ri[j].nr);
2458 gmx_fio_do_uchar(fio, ri[j].ic);
2468 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2470 gmx_groups_t *groups)
2474 gmx_fio_do_int(fio, atoms->nr);
2475 gmx_fio_do_int(fio, atoms->nres);
2476 if (file_version < 57)
2478 gmx_fio_do_int(fio, groups->ngrpname);
2479 for (i = 0; i < egcNR; i++)
2481 groups->ngrpnr[i] = atoms->nr;
2482 snew(groups->grpnr[i], groups->ngrpnr[i]);
2487 snew(atoms->atom, atoms->nr);
2488 snew(atoms->atomname, atoms->nr);
2489 snew(atoms->atomtype, atoms->nr);
2490 snew(atoms->atomtypeB, atoms->nr);
2491 snew(atoms->resinfo, atoms->nres);
2492 if (file_version < 57)
2494 snew(groups->grpname, groups->ngrpname);
2496 atoms->pdbinfo = NULL;
2498 for (i = 0; (i < atoms->nr); i++)
2500 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2502 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2503 if (bRead && (file_version <= 20))
2505 for (i = 0; i < atoms->nr; i++)
2507 atoms->atomtype[i] = put_symtab(symtab, "?");
2508 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2513 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2514 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2516 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2518 if (file_version < 57)
2520 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2522 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2526 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2527 gmx_bool bRead, t_symtab *symtab,
2532 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2533 gmx_fio_do_int(fio, groups->ngrpname);
2536 snew(groups->grpname, groups->ngrpname);
2538 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2539 for (g = 0; g < egcNR; g++)
2541 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2542 if (groups->ngrpnr[g] == 0)
2546 groups->grpnr[g] = NULL;
2553 snew(groups->grpnr[g], groups->ngrpnr[g]);
2555 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2560 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2565 if (file_version > 25)
2567 gmx_fio_do_int(fio, atomtypes->nr);
2571 snew(atomtypes->radius, j);
2572 snew(atomtypes->vol, j);
2573 snew(atomtypes->surftens, j);
2574 snew(atomtypes->atomnumber, j);
2575 snew(atomtypes->gb_radius, j);
2576 snew(atomtypes->S_hct, j);
2578 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2579 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2580 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2581 if (file_version >= 40)
2583 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2585 if (file_version >= 60)
2587 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2588 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2593 /* File versions prior to 26 cannot do GBSA,
2594 * so they dont use this structure
2597 atomtypes->radius = NULL;
2598 atomtypes->vol = NULL;
2599 atomtypes->surftens = NULL;
2600 atomtypes->atomnumber = NULL;
2601 atomtypes->gb_radius = NULL;
2602 atomtypes->S_hct = NULL;
2606 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2612 gmx_fio_do_int(fio, symtab->nr);
2616 snew(symtab->symbuf, 1);
2617 symbuf = symtab->symbuf;
2618 symbuf->bufsize = nr;
2619 snew(symbuf->buf, nr);
2620 for (i = 0; (i < nr); i++)
2622 gmx_fio_do_string(fio, buf);
2623 symbuf->buf[i] = strdup(buf);
2628 symbuf = symtab->symbuf;
2629 while (symbuf != NULL)
2631 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2633 gmx_fio_do_string(fio, symbuf->buf[i]);
2636 symbuf = symbuf->next;
2640 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2645 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2647 int i, j, ngrid, gs, nelem;
2649 gmx_fio_do_int(fio, cmap_grid->ngrid);
2650 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2652 ngrid = cmap_grid->ngrid;
2653 gs = cmap_grid->grid_spacing;
2658 snew(cmap_grid->cmapdata, ngrid);
2660 for (i = 0; i < cmap_grid->ngrid; i++)
2662 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2666 for (i = 0; i < cmap_grid->ngrid; i++)
2668 for (j = 0; j < nelem; j++)
2670 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2671 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2672 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2673 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2679 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2681 int m, a, a0, a1, r;
2685 /* We always assign a new chain number, but save the chain id characters
2686 * for larger molecules.
2688 #define CHAIN_MIN_ATOMS 15
2692 for (m = 0; m < mols->nr; m++)
2694 a0 = mols->index[m];
2695 a1 = mols->index[m+1];
2696 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2705 for (a = a0; a < a1; a++)
2707 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2708 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2713 /* Blank out the chain id if there was only one chain */
2716 for (r = 0; r < atoms->nres; r++)
2718 atoms->resinfo[r].chainid = ' ';
2723 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2724 t_symtab *symtab, int file_version,
2725 gmx_groups_t *groups)
2729 if (file_version >= 57)
2731 do_symstr(fio, &(molt->name), bRead, symtab);
2734 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2736 if (bRead && gmx_debug_at)
2738 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2741 if (file_version >= 57)
2743 do_ilists(fio, molt->ilist, bRead, file_version);
2745 do_block(fio, &molt->cgs, bRead, file_version);
2746 if (bRead && gmx_debug_at)
2748 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2752 /* This used to be in the atoms struct */
2753 do_blocka(fio, &molt->excls, bRead, file_version);
2756 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2760 gmx_fio_do_int(fio, molb->type);
2761 gmx_fio_do_int(fio, molb->nmol);
2762 gmx_fio_do_int(fio, molb->natoms_mol);
2763 /* Position restraint coordinates */
2764 gmx_fio_do_int(fio, molb->nposres_xA);
2765 if (molb->nposres_xA > 0)
2769 snew(molb->posres_xA, molb->nposres_xA);
2771 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2773 gmx_fio_do_int(fio, molb->nposres_xB);
2774 if (molb->nposres_xB > 0)
2778 snew(molb->posres_xB, molb->nposres_xB);
2780 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2785 static t_block mtop_mols(gmx_mtop_t *mtop)
2791 for (mb = 0; mb < mtop->nmolblock; mb++)
2793 mols.nr += mtop->molblock[mb].nmol;
2795 mols.nalloc_index = mols.nr + 1;
2796 snew(mols.index, mols.nalloc_index);
2801 for (mb = 0; mb < mtop->nmolblock; mb++)
2803 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2805 a += mtop->molblock[mb].natoms_mol;
2814 static void add_posres_molblock(gmx_mtop_t *mtop)
2819 gmx_molblock_t *molb;
2822 /* posres reference positions are stored in ip->posres (if present) and
2823 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2824 posres.pos0A are identical to fbposres.pos0. */
2825 il = &mtop->moltype[0].ilist[F_POSRES];
2826 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2827 if (il->nr == 0 && ilfb->nr == 0)
2833 for (i = 0; i < il->nr; i += 2)
2835 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2836 am = max(am, il->iatoms[i+1]);
2837 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2838 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2839 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2844 /* This loop is required if we have only flat-bottomed posres:
2846 - bFE == FALSE (no B-state for flat-bottomed posres) */
2849 for (i = 0; i < ilfb->nr; i += 2)
2851 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2852 am = max(am, ilfb->iatoms[i+1]);
2855 /* Make the posres coordinate block end at a molecule end */
2857 while (am >= mtop->mols.index[mol+1])
2861 molb = &mtop->molblock[0];
2862 molb->nposres_xA = mtop->mols.index[mol+1];
2863 snew(molb->posres_xA, molb->nposres_xA);
2866 molb->nposres_xB = molb->nposres_xA;
2867 snew(molb->posres_xB, molb->nposres_xB);
2871 molb->nposres_xB = 0;
2873 for (i = 0; i < il->nr; i += 2)
2875 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2876 a = il->iatoms[i+1];
2877 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2878 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2879 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2882 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2883 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2884 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2889 /* If only flat-bottomed posres are present, take reference pos from them.
2890 Here: bFE == FALSE */
2891 for (i = 0; i < ilfb->nr; i += 2)
2893 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2894 a = ilfb->iatoms[i+1];
2895 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2896 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2897 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2902 static void set_disres_npair(gmx_mtop_t *mtop)
2909 ip = mtop->ffparams.iparams;
2911 for (mt = 0; mt < mtop->nmoltype; mt++)
2913 il = &mtop->moltype[mt].ilist[F_DISRES];
2918 for (i = 0; i < il->nr; i += 3)
2921 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2923 ip[a[i]].disres.npair = npair;
2931 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2941 do_symtab(fio, &(mtop->symtab), bRead);
2944 pr_symtab(debug, 0, "symtab", &mtop->symtab);
2947 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2949 if (file_version >= 57)
2951 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
2953 gmx_fio_do_int(fio, mtop->nmoltype);
2961 snew(mtop->moltype, mtop->nmoltype);
2962 if (file_version < 57)
2964 mtop->moltype[0].name = mtop->name;
2967 for (mt = 0; mt < mtop->nmoltype; mt++)
2969 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
2973 if (file_version >= 57)
2975 gmx_fio_do_int(fio, mtop->nmolblock);
2979 mtop->nmolblock = 1;
2983 snew(mtop->molblock, mtop->nmolblock);
2985 if (file_version >= 57)
2987 for (mb = 0; mb < mtop->nmolblock; mb++)
2989 do_molblock(fio, &mtop->molblock[mb], bRead);
2991 gmx_fio_do_int(fio, mtop->natoms);
2995 mtop->molblock[0].type = 0;
2996 mtop->molblock[0].nmol = 1;
2997 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
2998 mtop->molblock[0].nposres_xA = 0;
2999 mtop->molblock[0].nposres_xB = 0;
3002 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
3005 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
3008 if (file_version < 57)
3010 /* Debug statements are inside do_idef */
3011 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
3012 mtop->natoms = mtop->moltype[0].atoms.nr;
3015 if (file_version >= 65)
3017 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
3021 mtop->ffparams.cmap_grid.ngrid = 0;
3022 mtop->ffparams.cmap_grid.grid_spacing = 0;
3023 mtop->ffparams.cmap_grid.cmapdata = NULL;
3026 if (file_version >= 57)
3028 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
3031 if (file_version < 57)
3033 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
3034 if (bRead && gmx_debug_at)
3036 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
3038 do_block(fio, &mtop->mols, bRead, file_version);
3039 /* Add the posres coordinates to the molblock */
3040 add_posres_molblock(mtop);
3044 if (file_version >= 57)
3046 done_block(&mtop->mols);
3047 mtop->mols = mtop_mols(mtop);
3051 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
3055 if (file_version < 51)
3057 /* Here used to be the shake blocks */
3058 do_blocka(fio, &dumb, bRead, file_version);
3071 close_symtab(&(mtop->symtab));
3075 /* If TopOnlyOK is TRUE then we can read even future versions
3076 * of tpx files, provided the file_generation hasn't changed.
3077 * If it is FALSE, we need the inputrecord too, and bail out
3078 * if the file is newer than the program.
3080 * The version and generation if the topology (see top of this file)
3081 * are returned in the two last arguments.
3083 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
3085 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
3086 gmx_bool TopOnlyOK, int *file_version,
3087 int *file_generation)
3090 char file_tag[STRLEN];
3097 gmx_fio_checktype(fio);
3098 gmx_fio_setdebug(fio, bDebugMode());
3100 /* NEW! XDR tpb file */
3101 precision = sizeof(real);
3104 gmx_fio_do_string(fio, buf);
3105 if (strncmp(buf, "VERSION", 7))
3107 gmx_fatal(FARGS, "Can not read file %s,\n"
3108 " this file is from a Gromacs version which is older than 2.0\n"
3109 " Make a new one with grompp or use a gro or pdb file, if possible",
3110 gmx_fio_getname(fio));
3112 gmx_fio_do_int(fio, precision);
3113 bDouble = (precision == sizeof(double));
3114 if ((precision != sizeof(float)) && !bDouble)
3116 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3117 "instead of %d or %d",
3118 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3120 gmx_fio_setprecision(fio, bDouble);
3121 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3122 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3126 gmx_fio_write_string(fio, GromacsVersion());
3127 bDouble = (precision == sizeof(double));
3128 gmx_fio_setprecision(fio, bDouble);
3129 gmx_fio_do_int(fio, precision);
3131 sprintf(file_tag, "%s", tpx_tag);
3132 fgen = tpx_generation;
3135 /* Check versions! */
3136 gmx_fio_do_int(fio, fver);
3138 /* This is for backward compatibility with development versions 77-79
3139 * where the tag was, mistakenly, placed before the generation,
3140 * which would cause a segv instead of a proper error message
3141 * when reading the topology only from tpx with <77 code.
3143 if (fver >= 77 && fver <= 79)
3145 gmx_fio_do_string(fio, file_tag);
3150 gmx_fio_do_int(fio, fgen);
3159 gmx_fio_do_string(fio, file_tag);
3165 /* Versions before 77 don't have the tag, set it to release */
3166 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3169 if (strcmp(file_tag, tpx_tag) != 0)
3171 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3174 /* We only support reading tpx files with the same tag as the code
3175 * or tpx files with the release tag and with lower version number.
3177 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
3179 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3180 gmx_fio_getname(fio), fver, file_tag,
3181 tpx_version, tpx_tag);
3186 if (file_version != NULL)
3188 *file_version = fver;
3190 if (file_generation != NULL)
3192 *file_generation = fgen;
3196 if ((fver <= tpx_incompatible_version) ||
3197 ((fver > tpx_version) && !TopOnlyOK) ||
3198 (fgen > tpx_generation) ||
3199 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3201 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3202 gmx_fio_getname(fio), fver, tpx_version);
3205 do_section(fio, eitemHEADER, bRead);
3206 gmx_fio_do_int(fio, tpx->natoms);
3209 gmx_fio_do_int(fio, tpx->ngtc);
3217 gmx_fio_do_int(fio, idum);
3218 gmx_fio_do_real(fio, rdum);
3220 /*a better decision will eventually (5.0 or later) need to be made
3221 on how to treat the alchemical state of the system, which can now
3222 vary through a simulation, and cannot be completely described
3223 though a single lambda variable, or even a single state
3224 index. Eventually, should probably be a vector. MRS*/
3227 gmx_fio_do_int(fio, tpx->fep_state);
3229 gmx_fio_do_real(fio, tpx->lambda);
3230 gmx_fio_do_int(fio, tpx->bIr);
3231 gmx_fio_do_int(fio, tpx->bTop);
3232 gmx_fio_do_int(fio, tpx->bX);
3233 gmx_fio_do_int(fio, tpx->bV);
3234 gmx_fio_do_int(fio, tpx->bF);
3235 gmx_fio_do_int(fio, tpx->bBox);
3237 if ((fgen > tpx_generation))
3239 /* This can only happen if TopOnlyOK=TRUE */
3244 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3245 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3246 gmx_bool bXVallocated)
3252 int file_version, file_generation;
3256 gmx_bool bPeriodicMols;
3260 tpx.natoms = state->natoms;
3261 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3262 tpx.fep_state = state->fep_state;
3263 tpx.lambda = state->lambda[efptFEP];
3264 tpx.bIr = (ir != NULL);
3265 tpx.bTop = (mtop != NULL);
3266 tpx.bX = (state->x != NULL);
3267 tpx.bV = (state->v != NULL);
3268 tpx.bF = (f != NULL);
3272 TopOnlyOK = (ir == NULL);
3274 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3279 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3280 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3285 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3286 state->natoms = tpx.natoms;
3287 state->nalloc = tpx.natoms;
3293 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3297 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3299 do_test(fio, tpx.bBox, state->box);
3300 do_section(fio, eitemBOX, bRead);
3303 gmx_fio_ndo_rvec(fio, state->box, DIM);
3304 if (file_version >= 51)
3306 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3310 /* We initialize box_rel after reading the inputrec */
3311 clear_mat(state->box_rel);
3313 if (file_version >= 28)
3315 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3316 if (file_version < 56)
3319 gmx_fio_ndo_rvec(fio, mdum, DIM);
3324 if (state->ngtc > 0 && file_version >= 28)
3327 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3328 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3329 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3330 snew(dumv, state->ngtc);
3331 if (file_version < 69)
3333 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3335 /* These used to be the Berendsen tcoupl_lambda's */
3336 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3340 /* Prior to tpx version 26, the inputrec was here.
3341 * I moved it to enable partial forward-compatibility
3342 * for analysis/viewer programs.
3344 if (file_version < 26)
3346 do_test(fio, tpx.bIr, ir);
3347 do_section(fio, eitemIR, bRead);
3352 do_inputrec(fio, ir, bRead, file_version,
3353 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3356 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3361 do_inputrec(fio, &dum_ir, bRead, file_version,
3362 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3365 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3367 done_inputrec(&dum_ir);
3373 do_test(fio, tpx.bTop, mtop);
3374 do_section(fio, eitemTOP, bRead);
3379 do_mtop(fio, mtop, bRead, file_version);
3383 do_mtop(fio, &dum_top, bRead, file_version);
3384 done_mtop(&dum_top, TRUE);
3387 do_test(fio, tpx.bX, state->x);
3388 do_section(fio, eitemX, bRead);
3393 state->flags |= (1<<estX);
3395 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3398 do_test(fio, tpx.bV, state->v);
3399 do_section(fio, eitemV, bRead);
3404 state->flags |= (1<<estV);
3406 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3409 do_test(fio, tpx.bF, f);
3410 do_section(fio, eitemF, bRead);
3413 gmx_fio_ndo_rvec(fio, f, state->natoms);
3416 /* Starting with tpx version 26, we have the inputrec
3417 * at the end of the file, so we can ignore it
3418 * if the file is never than the software (but still the
3419 * same generation - see comments at the top of this file.
3424 bPeriodicMols = FALSE;
3425 if (file_version >= 26)
3427 do_test(fio, tpx.bIr, ir);
3428 do_section(fio, eitemIR, bRead);
3431 if (file_version >= 53)
3433 /* Removed the pbc info from do_inputrec, since we always want it */
3437 bPeriodicMols = ir->bPeriodicMols;
3439 gmx_fio_do_int(fio, ePBC);
3440 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3442 if (file_generation <= tpx_generation && ir)
3444 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3447 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3449 if (file_version < 51)
3451 set_box_rel(ir, state);
3453 if (file_version < 53)
3456 bPeriodicMols = ir->bPeriodicMols;
3459 if (bRead && ir && file_version >= 53)
3461 /* We need to do this after do_inputrec, since that initializes ir */
3463 ir->bPeriodicMols = bPeriodicMols;
3472 if (state->ngtc == 0)
3474 /* Reading old version without tcoupl state data: set it */
3475 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3477 if (tpx.bTop && mtop)
3479 if (file_version < 57)
3481 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3483 ir->eDisre = edrSimple;
3487 ir->eDisre = edrNone;
3490 set_disres_npair(mtop);
3494 if (tpx.bTop && mtop)
3496 gmx_mtop_finalize(mtop);
3499 if (file_version >= 57)
3503 env = getenv("GMX_NOCHARGEGROUPS");
3506 sscanf(env, "%d", &ienv);
3507 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3512 "Will make single atomic charge groups in non-solvent%s\n",
3513 ienv > 1 ? " and solvent" : "");
3514 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3516 fprintf(stderr, "\n");
3524 /************************************************************
3526 * The following routines are the exported ones
3528 ************************************************************/
3530 t_fileio *open_tpx(const char *fn, const char *mode)
3532 return gmx_fio_open(fn, mode);
3535 void close_tpx(t_fileio *fio)
3540 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3541 int *file_version, int *file_generation)
3545 fio = open_tpx(fn, "r");
3546 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3550 void write_tpx_state(const char *fn,
3551 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3555 fio = open_tpx(fn, "w");
3556 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3560 void read_tpx_state(const char *fn,
3561 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3565 fio = open_tpx(fn, "r");
3566 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3570 int read_tpx(const char *fn,
3571 t_inputrec *ir, matrix box, int *natoms,
3572 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3580 fio = open_tpx(fn, "r");
3581 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3583 *natoms = state.natoms;
3586 copy_mat(state.box, box);
3595 int read_tpx_top(const char *fn,
3596 t_inputrec *ir, matrix box, int *natoms,
3597 rvec *x, rvec *v, rvec *f, t_topology *top)
3603 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3605 *top = gmx_mtop_t_to_t_topology(&mtop);
3610 gmx_bool fn2bTPX(const char *file)
3612 switch (fn2ftp(file))
3623 static void done_gmx_groups_t(gmx_groups_t *g)
3627 for (i = 0; (i < egcNR); i++)
3629 if (NULL != g->grps[i].nm_ind)
3631 sfree(g->grps[i].nm_ind);
3632 g->grps[i].nm_ind = NULL;
3634 if (NULL != g->grpnr[i])
3640 /* The contents of this array is in symtab, don't free it here */
3644 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3645 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3648 int natoms, i, version, generation;
3649 gmx_bool bTop, bXNULL = FALSE;
3651 t_topology *topconv;
3654 bTop = fn2bTPX(infile);
3658 read_tpxheader(infile, &header, TRUE, &version, &generation);
3661 snew(*x, header.natoms);
3665 snew(*v, header.natoms);
3668 *ePBC = read_tpx(infile, NULL, box, &natoms,
3669 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3670 *top = gmx_mtop_t_to_t_topology(mtop);
3671 /* In this case we need to throw away the group data too */
3672 done_gmx_groups_t(&mtop->groups);
3674 strcpy(title, *top->name);
3675 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3679 get_stx_coordnum(infile, &natoms);
3680 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3691 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3699 aps = gmx_atomprop_init();
3700 for (i = 0; (i < natoms); i++)
3702 if (!gmx_atomprop_query(aps, epropMass,
3703 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3704 *top->atoms.atomname[i],
3705 &(top->atoms.atom[i].m)))
3709 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3710 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3711 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3712 *top->atoms.atomname[i]);
3716 gmx_atomprop_destroy(aps);
3718 top->idef.ntypes = -1;