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 /* This is the tag string which is stored in the tpx file.
65 * Change this if you want to change the tpx format in a feature branch.
66 * This ensures that there will not be different tpx formats around which
67 * can not be distinguished.
69 static const char *tpx_tag = TPX_TAG_RELEASE;
72 /* The tpx_version number should be increased whenever the file format changes!
74 * The following comment section helps to keep track of which feature has been
75 * added in which version.
77 * version feature added
78 * 96 support for ion/water position swaps (computational electrophysiology)
80 static const int tpx_version = 96;
83 /* This number should only be increased when you edit the TOPOLOGY section
84 * or the HEADER of the tpx format.
85 * This way we can maintain forward compatibility too for all analysis tools
86 * and/or external programs that only need to know the atom/residue names,
87 * charges, and bond connectivity.
89 * It first appeared in tpx version 26, when I also moved the inputrecord
90 * to the end of the tpx file, so we can just skip it if we only
93 static const int tpx_generation = 25;
95 /* This number should be the most recent backwards incompatible version
96 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
98 static const int tpx_incompatible_version = 9;
102 /* Struct used to maintain tpx compatibility when function types are added */
104 int fvnr; /* file version number in which the function type first appeared */
105 int ftype; /* function type */
109 * The entries should be ordered in:
110 * 1. ascending file version number
111 * 2. ascending function type number
113 /*static const t_ftupd ftupd[] = {
114 { 20, F_CUBICBONDS },
118 { 22, F_DISRESVIOL },
124 { 26, F_DIHRESVIOL },
125 { 30, F_CROSS_BOND_BONDS },
126 { 30, F_CROSS_BOND_ANGLES },
127 { 30, F_UREY_BRADLEY },
128 { 30, F_POLARIZATION },
132 * The entries should be ordered in:
133 * 1. ascending function type number
134 * 2. ascending file version number
136 /* question; what is the purpose of the commented code above? */
137 static const t_ftupd ftupd[] = {
138 { 20, F_CUBICBONDS },
143 { 43, F_TABBONDSNC },
144 { 70, F_RESTRBONDS },
145 { 76, F_LINEAR_ANGLES },
146 { 30, F_CROSS_BOND_BONDS },
147 { 30, F_CROSS_BOND_ANGLES },
148 { 30, F_UREY_BRADLEY },
149 { 34, F_QUARTIC_ANGLES },
159 { 72, F_NPSOLVATION },
161 { 41, F_LJC_PAIRS_NB },
164 { 32, F_COUL_RECIP },
167 { 30, F_POLARIZATION },
170 { 22, F_DISRESVIOL },
174 { 26, F_DIHRESVIOL },
179 { 46, F_ECONSERVED },
180 { 69, F_VTEMP_NOLONGERUSED},
182 { 54, F_DVDL_CONSTR },
183 { 76, F_ANHARM_POL },
186 { 79, F_DVDL_BONDED, },
187 { 79, F_DVDL_RESTRAINT },
188 { 79, F_DVDL_TEMPERATURE },
190 #define NFTUPD asize(ftupd)
192 /* Needed for backward compatibility */
195 static void _do_section(t_fileio *fio, int key, gmx_bool bRead, const char *src,
201 if (gmx_fio_getftp(fio) == efTPA)
205 gmx_fio_write_string(fio, itemstr[key]);
206 bDbg = gmx_fio_getdebug(fio);
207 gmx_fio_setdebug(fio, FALSE);
208 gmx_fio_write_string(fio, comment_str[key]);
209 gmx_fio_setdebug(fio, bDbg);
213 if (gmx_fio_getdebug(fio))
215 fprintf(stderr, "Looking for section %s (%s, %d)",
216 itemstr[key], src, line);
221 gmx_fio_do_string(fio, buf);
223 while ((gmx_strcasecmp(buf, itemstr[key]) != 0));
225 if (gmx_strcasecmp(buf, itemstr[key]) != 0)
227 gmx_fatal(FARGS, "\nCould not find section heading %s", itemstr[key]);
229 else if (gmx_fio_getdebug(fio))
231 fprintf(stderr, " and found it\n");
237 #define do_section(fio, key, bRead) _do_section(fio, key, bRead, __FILE__, __LINE__)
239 /**************************************************************
241 * Now the higer level routines that do io of the structures and arrays
243 **************************************************************/
244 static void do_pullgrp_tpx_pre95(t_fileio *fio,
253 gmx_fio_do_int(fio, pgrp->nat);
256 snew(pgrp->ind, pgrp->nat);
258 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
259 gmx_fio_do_int(fio, pgrp->nweight);
262 snew(pgrp->weight, pgrp->nweight);
264 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
265 gmx_fio_do_int(fio, pgrp->pbcatom);
266 gmx_fio_do_rvec(fio, pcrd->vec);
267 clear_rvec(pcrd->origin);
268 gmx_fio_do_rvec(fio, tmp);
270 gmx_fio_do_real(fio, pcrd->rate);
271 gmx_fio_do_real(fio, pcrd->k);
272 if (file_version >= 56)
274 gmx_fio_do_real(fio, pcrd->kB);
282 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
286 gmx_fio_do_int(fio, pgrp->nat);
289 snew(pgrp->ind, pgrp->nat);
291 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
292 gmx_fio_do_int(fio, pgrp->nweight);
295 snew(pgrp->weight, pgrp->nweight);
297 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
298 gmx_fio_do_int(fio, pgrp->pbcatom);
301 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd)
305 gmx_fio_do_int(fio, pcrd->group[0]);
306 gmx_fio_do_int(fio, pcrd->group[1]);
307 gmx_fio_do_rvec(fio, pcrd->origin);
308 gmx_fio_do_rvec(fio, pcrd->vec);
309 gmx_fio_do_real(fio, pcrd->init);
310 gmx_fio_do_real(fio, pcrd->rate);
311 gmx_fio_do_real(fio, pcrd->k);
312 gmx_fio_do_real(fio, pcrd->kB);
315 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
317 /* i is used in the ndo_double macro*/
321 int n_lambda = fepvals->n_lambda;
323 /* reset the lambda calculation window */
324 fepvals->lambda_start_n = 0;
325 fepvals->lambda_stop_n = n_lambda;
326 if (file_version >= 79)
332 snew(expand->init_lambda_weights, n_lambda);
334 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
335 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
338 gmx_fio_do_int(fio, expand->nstexpanded);
339 gmx_fio_do_int(fio, expand->elmcmove);
340 gmx_fio_do_int(fio, expand->elamstats);
341 gmx_fio_do_int(fio, expand->lmc_repeats);
342 gmx_fio_do_int(fio, expand->gibbsdeltalam);
343 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
344 gmx_fio_do_int(fio, expand->lmc_seed);
345 gmx_fio_do_real(fio, expand->mc_temp);
346 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
347 gmx_fio_do_int(fio, expand->nstTij);
348 gmx_fio_do_int(fio, expand->minvarmin);
349 gmx_fio_do_int(fio, expand->c_range);
350 gmx_fio_do_real(fio, expand->wl_scale);
351 gmx_fio_do_real(fio, expand->wl_ratio);
352 gmx_fio_do_real(fio, expand->init_wl_delta);
353 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
354 gmx_fio_do_int(fio, expand->elmceq);
355 gmx_fio_do_int(fio, expand->equil_steps);
356 gmx_fio_do_int(fio, expand->equil_samples);
357 gmx_fio_do_int(fio, expand->equil_n_at_lam);
358 gmx_fio_do_real(fio, expand->equil_wl_delta);
359 gmx_fio_do_real(fio, expand->equil_ratio);
363 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
366 if (file_version >= 79)
368 gmx_fio_do_int(fio, simtemp->eSimTempScale);
369 gmx_fio_do_real(fio, simtemp->simtemp_high);
370 gmx_fio_do_real(fio, simtemp->simtemp_low);
375 snew(simtemp->temperatures, n_lambda);
377 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
382 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
384 /* i is defined in the ndo_double macro; use g to iterate. */
389 /* free energy values */
391 if (file_version >= 79)
393 gmx_fio_do_int(fio, fepvals->init_fep_state);
394 gmx_fio_do_double(fio, fepvals->init_lambda);
395 gmx_fio_do_double(fio, fepvals->delta_lambda);
397 else if (file_version >= 59)
399 gmx_fio_do_double(fio, fepvals->init_lambda);
400 gmx_fio_do_double(fio, fepvals->delta_lambda);
404 gmx_fio_do_real(fio, rdum);
405 fepvals->init_lambda = rdum;
406 gmx_fio_do_real(fio, rdum);
407 fepvals->delta_lambda = rdum;
409 if (file_version >= 79)
411 gmx_fio_do_int(fio, fepvals->n_lambda);
414 snew(fepvals->all_lambda, efptNR);
416 for (g = 0; g < efptNR; g++)
418 if (fepvals->n_lambda > 0)
422 snew(fepvals->all_lambda[g], fepvals->n_lambda);
424 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
425 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
427 else if (fepvals->init_lambda >= 0)
429 fepvals->separate_dvdl[efptFEP] = TRUE;
433 else if (file_version >= 64)
435 gmx_fio_do_int(fio, fepvals->n_lambda);
440 snew(fepvals->all_lambda, efptNR);
441 /* still allocate the all_lambda array's contents. */
442 for (g = 0; g < efptNR; g++)
444 if (fepvals->n_lambda > 0)
446 snew(fepvals->all_lambda[g], fepvals->n_lambda);
450 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
452 if (fepvals->init_lambda >= 0)
456 fepvals->separate_dvdl[efptFEP] = TRUE;
460 /* copy the contents of the efptFEP lambda component to all
461 the other components */
462 for (g = 0; g < efptNR; g++)
464 for (h = 0; h < fepvals->n_lambda; h++)
468 fepvals->all_lambda[g][h] =
469 fepvals->all_lambda[efptFEP][h];
478 fepvals->n_lambda = 0;
479 fepvals->all_lambda = NULL;
480 if (fepvals->init_lambda >= 0)
482 fepvals->separate_dvdl[efptFEP] = TRUE;
485 if (file_version >= 13)
487 gmx_fio_do_real(fio, fepvals->sc_alpha);
491 fepvals->sc_alpha = 0;
493 if (file_version >= 38)
495 gmx_fio_do_int(fio, fepvals->sc_power);
499 fepvals->sc_power = 2;
501 if (file_version >= 79)
503 gmx_fio_do_real(fio, fepvals->sc_r_power);
507 fepvals->sc_r_power = 6.0;
509 if (file_version >= 15)
511 gmx_fio_do_real(fio, fepvals->sc_sigma);
515 fepvals->sc_sigma = 0.3;
519 if (file_version >= 71)
521 fepvals->sc_sigma_min = fepvals->sc_sigma;
525 fepvals->sc_sigma_min = 0;
528 if (file_version >= 79)
530 gmx_fio_do_int(fio, fepvals->bScCoul);
534 fepvals->bScCoul = TRUE;
536 if (file_version >= 64)
538 gmx_fio_do_int(fio, fepvals->nstdhdl);
542 fepvals->nstdhdl = 1;
545 if (file_version >= 73)
547 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
548 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
552 fepvals->separate_dhdl_file = esepdhdlfileYES;
553 fepvals->dhdl_derivatives = edhdlderivativesYES;
555 if (file_version >= 71)
557 gmx_fio_do_int(fio, fepvals->dh_hist_size);
558 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
562 fepvals->dh_hist_size = 0;
563 fepvals->dh_hist_spacing = 0.1;
565 if (file_version >= 79)
567 gmx_fio_do_int(fio, fepvals->bPrintEnergy);
571 fepvals->bPrintEnergy = FALSE;
574 /* handle lambda_neighbors */
575 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
577 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
578 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
579 (fepvals->init_lambda < 0) )
581 fepvals->lambda_start_n = (fepvals->init_fep_state -
582 fepvals->lambda_neighbors);
583 fepvals->lambda_stop_n = (fepvals->init_fep_state +
584 fepvals->lambda_neighbors + 1);
585 if (fepvals->lambda_start_n < 0)
587 fepvals->lambda_start_n = 0;;
589 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
591 fepvals->lambda_stop_n = fepvals->n_lambda;
596 fepvals->lambda_start_n = 0;
597 fepvals->lambda_stop_n = fepvals->n_lambda;
602 fepvals->lambda_start_n = 0;
603 fepvals->lambda_stop_n = fepvals->n_lambda;
607 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
611 if (file_version >= 95)
613 gmx_fio_do_int(fio, pull->ngroup);
615 gmx_fio_do_int(fio, pull->ncoord);
616 if (file_version < 95)
618 pull->ngroup = pull->ncoord + 1;
620 gmx_fio_do_int(fio, pull->eGeom);
621 gmx_fio_do_ivec(fio, pull->dim);
622 gmx_fio_do_real(fio, pull->cyl_r1);
623 gmx_fio_do_real(fio, pull->cyl_r0);
624 gmx_fio_do_real(fio, pull->constr_tol);
625 if (file_version >= 95)
627 gmx_fio_do_int(fio, pull->bPrintRef);
629 gmx_fio_do_int(fio, pull->nstxout);
630 gmx_fio_do_int(fio, pull->nstfout);
633 snew(pull->group, pull->ngroup);
634 snew(pull->coord, pull->ncoord);
636 if (file_version < 95)
638 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
639 if (pull->eGeom == epullgDIRPBC)
641 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
643 if (pull->eGeom > epullgDIRPBC)
648 for (g = 0; g < pull->ngroup; g++)
650 /* We read and ignore a pull coordinate for group 0 */
651 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[max(g-1, 0)],
652 bRead, file_version);
655 pull->coord[g-1].group[0] = 0;
656 pull->coord[g-1].group[1] = g;
660 pull->bPrintRef = (pull->group[0].nat > 0);
664 for (g = 0; g < pull->ngroup; g++)
666 do_pull_group(fio, &pull->group[g], bRead);
668 for (g = 0; g < pull->ncoord; g++)
670 do_pull_coord(fio, &pull->coord[g]);
676 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
680 gmx_fio_do_int(fio, rotg->eType);
681 gmx_fio_do_int(fio, rotg->bMassW);
682 gmx_fio_do_int(fio, rotg->nat);
685 snew(rotg->ind, rotg->nat);
687 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
690 snew(rotg->x_ref, rotg->nat);
692 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
693 gmx_fio_do_rvec(fio, rotg->vec);
694 gmx_fio_do_rvec(fio, rotg->pivot);
695 gmx_fio_do_real(fio, rotg->rate);
696 gmx_fio_do_real(fio, rotg->k);
697 gmx_fio_do_real(fio, rotg->slab_dist);
698 gmx_fio_do_real(fio, rotg->min_gaussian);
699 gmx_fio_do_real(fio, rotg->eps);
700 gmx_fio_do_int(fio, rotg->eFittype);
701 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
702 gmx_fio_do_real(fio, rotg->PotAngle_step);
705 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
709 gmx_fio_do_int(fio, rot->ngrp);
710 gmx_fio_do_int(fio, rot->nstrout);
711 gmx_fio_do_int(fio, rot->nstsout);
714 snew(rot->grp, rot->ngrp);
716 for (g = 0; g < rot->ngrp; g++)
718 do_rotgrp(fio, &rot->grp[g], bRead);
723 static void do_swapcoords(t_fileio *fio, t_swapcoords *swap, gmx_bool bRead)
728 gmx_fio_do_int(fio, swap->nat);
729 gmx_fio_do_int(fio, swap->nat_sol);
730 for (j = 0; j < 2; j++)
732 gmx_fio_do_int(fio, swap->nat_split[j]);
733 gmx_fio_do_int(fio, swap->massw_split[j]);
735 gmx_fio_do_int(fio, swap->nstswap);
736 gmx_fio_do_int(fio, swap->nAverage);
737 gmx_fio_do_real(fio, swap->threshold);
738 gmx_fio_do_real(fio, swap->cyl0r);
739 gmx_fio_do_real(fio, swap->cyl0u);
740 gmx_fio_do_real(fio, swap->cyl0l);
741 gmx_fio_do_real(fio, swap->cyl1r);
742 gmx_fio_do_real(fio, swap->cyl1u);
743 gmx_fio_do_real(fio, swap->cyl1l);
747 snew(swap->ind, swap->nat);
748 snew(swap->ind_sol, swap->nat_sol);
749 for (j = 0; j < 2; j++)
751 snew(swap->ind_split[j], swap->nat_split[j]);
755 gmx_fio_ndo_int(fio, swap->ind, swap->nat);
756 gmx_fio_ndo_int(fio, swap->ind_sol, swap->nat_sol);
757 for (j = 0; j < 2; j++)
759 gmx_fio_ndo_int(fio, swap->ind_split[j], swap->nat_split[j]);
762 for (j = 0; j < eCompNR; j++)
764 gmx_fio_do_int(fio, swap->nanions[j]);
765 gmx_fio_do_int(fio, swap->ncations[j]);
771 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
772 int file_version, real *fudgeQQ)
774 int i, j, k, *tmp, idum = 0;
778 real zerotemptime, finish_t, init_temp, finish_temp;
780 if (file_version != tpx_version)
782 /* Give a warning about features that are not accessible */
783 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
784 file_version, tpx_version);
792 if (file_version == 0)
797 /* Basic inputrec stuff */
798 gmx_fio_do_int(fio, ir->eI);
799 if (file_version >= 62)
801 gmx_fio_do_int64(fio, ir->nsteps);
805 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, ir->ndelta);
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 if (file_version > 25)
914 gmx_fio_do_int(fio, ir->nstcheckpoint);
918 ir->nstcheckpoint = 0;
921 gmx_fio_do_int(fio, ir->nstcgsteep);
923 if (file_version >= 30)
925 gmx_fio_do_int(fio, ir->nbfgscorr);
932 gmx_fio_do_int(fio, ir->nstlog);
933 gmx_fio_do_int(fio, ir->nstxout);
934 gmx_fio_do_int(fio, ir->nstvout);
935 gmx_fio_do_int(fio, ir->nstfout);
936 gmx_fio_do_int(fio, ir->nstenergy);
937 gmx_fio_do_int(fio, ir->nstxout_compressed);
938 if (file_version >= 59)
940 gmx_fio_do_double(fio, ir->init_t);
941 gmx_fio_do_double(fio, ir->delta_t);
945 gmx_fio_do_real(fio, rdum);
947 gmx_fio_do_real(fio, rdum);
950 gmx_fio_do_real(fio, ir->x_compression_precision);
951 if (file_version < 19)
953 gmx_fio_do_int(fio, idum);
954 gmx_fio_do_int(fio, idum);
956 if (file_version < 18)
958 gmx_fio_do_int(fio, idum);
960 if (file_version >= 81)
962 gmx_fio_do_real(fio, ir->verletbuf_tol);
966 ir->verletbuf_tol = 0;
968 gmx_fio_do_real(fio, ir->rlist);
969 if (file_version >= 67)
971 gmx_fio_do_real(fio, ir->rlistlong);
973 if (file_version >= 82 && file_version != 90)
975 gmx_fio_do_int(fio, ir->nstcalclr);
979 /* Calculate at NS steps */
980 ir->nstcalclr = ir->nstlist;
982 gmx_fio_do_int(fio, ir->coulombtype);
983 if (file_version < 32 && ir->coulombtype == eelRF)
985 ir->coulombtype = eelRF_NEC;
987 if (file_version >= 81)
989 gmx_fio_do_int(fio, ir->coulomb_modifier);
993 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
995 gmx_fio_do_real(fio, ir->rcoulomb_switch);
996 gmx_fio_do_real(fio, ir->rcoulomb);
997 gmx_fio_do_int(fio, ir->vdwtype);
998 if (file_version >= 81)
1000 gmx_fio_do_int(fio, ir->vdw_modifier);
1004 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
1006 gmx_fio_do_real(fio, ir->rvdw_switch);
1007 gmx_fio_do_real(fio, ir->rvdw);
1008 if (file_version < 67)
1010 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
1012 gmx_fio_do_int(fio, ir->eDispCorr);
1013 gmx_fio_do_real(fio, ir->epsilon_r);
1014 if (file_version >= 37)
1016 gmx_fio_do_real(fio, ir->epsilon_rf);
1020 if (EEL_RF(ir->coulombtype))
1022 ir->epsilon_rf = ir->epsilon_r;
1023 ir->epsilon_r = 1.0;
1027 ir->epsilon_rf = 1.0;
1030 if (file_version >= 29)
1032 gmx_fio_do_real(fio, ir->tabext);
1039 if (file_version > 25)
1041 gmx_fio_do_int(fio, ir->gb_algorithm);
1042 gmx_fio_do_int(fio, ir->nstgbradii);
1043 gmx_fio_do_real(fio, ir->rgbradii);
1044 gmx_fio_do_real(fio, ir->gb_saltconc);
1045 gmx_fio_do_int(fio, ir->implicit_solvent);
1049 ir->gb_algorithm = egbSTILL;
1052 ir->gb_saltconc = 0;
1053 ir->implicit_solvent = eisNO;
1055 if (file_version >= 55)
1057 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1058 gmx_fio_do_real(fio, ir->gb_obc_alpha);
1059 gmx_fio_do_real(fio, ir->gb_obc_beta);
1060 gmx_fio_do_real(fio, ir->gb_obc_gamma);
1061 if (file_version >= 60)
1063 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1064 gmx_fio_do_int(fio, ir->sa_algorithm);
1068 ir->gb_dielectric_offset = 0.009;
1069 ir->sa_algorithm = esaAPPROX;
1071 gmx_fio_do_real(fio, ir->sa_surface_tension);
1073 /* Override sa_surface_tension if it is not changed in the mpd-file */
1074 if (ir->sa_surface_tension < 0)
1076 if (ir->gb_algorithm == egbSTILL)
1078 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1080 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1082 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1089 /* Better use sensible values than insane (0.0) ones... */
1090 ir->gb_epsilon_solvent = 80;
1091 ir->gb_obc_alpha = 1.0;
1092 ir->gb_obc_beta = 0.8;
1093 ir->gb_obc_gamma = 4.85;
1094 ir->sa_surface_tension = 2.092;
1098 if (file_version >= 81)
1100 gmx_fio_do_real(fio, ir->fourier_spacing);
1104 ir->fourier_spacing = 0.0;
1106 gmx_fio_do_int(fio, ir->nkx);
1107 gmx_fio_do_int(fio, ir->nky);
1108 gmx_fio_do_int(fio, ir->nkz);
1109 gmx_fio_do_int(fio, ir->pme_order);
1110 gmx_fio_do_real(fio, ir->ewald_rtol);
1112 if (file_version >= 93)
1114 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1118 ir->ewald_rtol_lj = ir->ewald_rtol;
1121 if (file_version >= 24)
1123 gmx_fio_do_int(fio, ir->ewald_geometry);
1127 ir->ewald_geometry = eewg3D;
1130 if (file_version <= 17)
1132 ir->epsilon_surface = 0;
1133 if (file_version == 17)
1135 gmx_fio_do_int(fio, idum);
1140 gmx_fio_do_real(fio, ir->epsilon_surface);
1143 gmx_fio_do_gmx_bool(fio, ir->bOptFFT);
1145 if (file_version >= 93)
1147 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1149 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1150 gmx_fio_do_int(fio, ir->etc);
1151 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1152 * but the values 0 and 1 still mean no and
1153 * berendsen temperature coupling, respectively.
1155 if (file_version >= 79)
1157 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1159 if (file_version >= 71)
1161 gmx_fio_do_int(fio, ir->nsttcouple);
1165 ir->nsttcouple = ir->nstcalcenergy;
1167 if (file_version <= 15)
1169 gmx_fio_do_int(fio, idum);
1171 if (file_version <= 17)
1173 gmx_fio_do_int(fio, ir->epct);
1174 if (file_version <= 15)
1178 ir->epct = epctSURFACETENSION;
1180 gmx_fio_do_int(fio, idum);
1183 /* we have removed the NO alternative at the beginning */
1187 ir->epct = epctISOTROPIC;
1191 ir->epc = epcBERENDSEN;
1196 gmx_fio_do_int(fio, ir->epc);
1197 gmx_fio_do_int(fio, ir->epct);
1199 if (file_version >= 71)
1201 gmx_fio_do_int(fio, ir->nstpcouple);
1205 ir->nstpcouple = ir->nstcalcenergy;
1207 gmx_fio_do_real(fio, ir->tau_p);
1208 if (file_version <= 15)
1210 gmx_fio_do_rvec(fio, vdum);
1211 clear_mat(ir->ref_p);
1212 for (i = 0; i < DIM; i++)
1214 ir->ref_p[i][i] = vdum[i];
1219 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1220 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1221 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1223 if (file_version <= 15)
1225 gmx_fio_do_rvec(fio, vdum);
1226 clear_mat(ir->compress);
1227 for (i = 0; i < DIM; i++)
1229 ir->compress[i][i] = vdum[i];
1234 gmx_fio_do_rvec(fio, ir->compress[XX]);
1235 gmx_fio_do_rvec(fio, ir->compress[YY]);
1236 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1238 if (file_version >= 47)
1240 gmx_fio_do_int(fio, ir->refcoord_scaling);
1241 gmx_fio_do_rvec(fio, ir->posres_com);
1242 gmx_fio_do_rvec(fio, ir->posres_comB);
1246 ir->refcoord_scaling = erscNO;
1247 clear_rvec(ir->posres_com);
1248 clear_rvec(ir->posres_comB);
1250 if ((file_version > 25) && (file_version < 79))
1252 gmx_fio_do_int(fio, ir->andersen_seed);
1256 ir->andersen_seed = 0;
1258 if (file_version < 26)
1260 gmx_fio_do_gmx_bool(fio, bSimAnn);
1261 gmx_fio_do_real(fio, zerotemptime);
1264 if (file_version < 37)
1266 gmx_fio_do_real(fio, rdum);
1269 gmx_fio_do_real(fio, ir->shake_tol);
1270 if (file_version < 54)
1272 gmx_fio_do_real(fio, *fudgeQQ);
1275 gmx_fio_do_int(fio, ir->efep);
1276 if (file_version <= 14 && ir->efep != efepNO)
1280 do_fepvals(fio, ir->fepvals, bRead, file_version);
1282 if (file_version >= 79)
1284 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1287 ir->bSimTemp = TRUE;
1292 ir->bSimTemp = FALSE;
1296 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1299 if (file_version >= 79)
1301 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1304 ir->bExpanded = TRUE;
1308 ir->bExpanded = FALSE;
1313 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1315 if (file_version >= 57)
1317 gmx_fio_do_int(fio, ir->eDisre);
1319 gmx_fio_do_int(fio, ir->eDisreWeighting);
1320 if (file_version < 22)
1322 if (ir->eDisreWeighting == 0)
1324 ir->eDisreWeighting = edrwEqual;
1328 ir->eDisreWeighting = edrwConservative;
1331 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1332 gmx_fio_do_real(fio, ir->dr_fc);
1333 gmx_fio_do_real(fio, ir->dr_tau);
1334 gmx_fio_do_int(fio, ir->nstdisreout);
1335 if (file_version >= 22)
1337 gmx_fio_do_real(fio, ir->orires_fc);
1338 gmx_fio_do_real(fio, ir->orires_tau);
1339 gmx_fio_do_int(fio, ir->nstorireout);
1345 ir->nstorireout = 0;
1347 if (file_version >= 26 && file_version < 79)
1349 gmx_fio_do_real(fio, ir->dihre_fc);
1350 if (file_version < 56)
1352 gmx_fio_do_real(fio, rdum);
1353 gmx_fio_do_int(fio, idum);
1361 gmx_fio_do_real(fio, ir->em_stepsize);
1362 gmx_fio_do_real(fio, ir->em_tol);
1363 if (file_version >= 22)
1365 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1369 ir->bShakeSOR = TRUE;
1371 if (file_version >= 11)
1373 gmx_fio_do_int(fio, ir->niter);
1378 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1381 if (file_version >= 21)
1383 gmx_fio_do_real(fio, ir->fc_stepsize);
1387 ir->fc_stepsize = 0;
1389 gmx_fio_do_int(fio, ir->eConstrAlg);
1390 gmx_fio_do_int(fio, ir->nProjOrder);
1391 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1392 if (file_version <= 14)
1394 gmx_fio_do_int(fio, idum);
1396 if (file_version >= 26)
1398 gmx_fio_do_int(fio, ir->nLincsIter);
1403 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1406 if (file_version < 33)
1408 gmx_fio_do_real(fio, bd_temp);
1410 gmx_fio_do_real(fio, ir->bd_fric);
1411 gmx_fio_do_int(fio, ir->ld_seed);
1412 if (file_version >= 33)
1414 for (i = 0; i < DIM; i++)
1416 gmx_fio_do_rvec(fio, ir->deform[i]);
1421 for (i = 0; i < DIM; i++)
1423 clear_rvec(ir->deform[i]);
1426 if (file_version >= 14)
1428 gmx_fio_do_real(fio, ir->cos_accel);
1434 gmx_fio_do_int(fio, ir->userint1);
1435 gmx_fio_do_int(fio, ir->userint2);
1436 gmx_fio_do_int(fio, ir->userint3);
1437 gmx_fio_do_int(fio, ir->userint4);
1438 gmx_fio_do_real(fio, ir->userreal1);
1439 gmx_fio_do_real(fio, ir->userreal2);
1440 gmx_fio_do_real(fio, ir->userreal3);
1441 gmx_fio_do_real(fio, ir->userreal4);
1444 if (file_version >= 77)
1446 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1451 snew(ir->adress, 1);
1453 gmx_fio_do_int(fio, ir->adress->type);
1454 gmx_fio_do_real(fio, ir->adress->const_wf);
1455 gmx_fio_do_real(fio, ir->adress->ex_width);
1456 gmx_fio_do_real(fio, ir->adress->hy_width);
1457 gmx_fio_do_int(fio, ir->adress->icor);
1458 gmx_fio_do_int(fio, ir->adress->site);
1459 gmx_fio_do_rvec(fio, ir->adress->refs);
1460 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1461 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1462 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1463 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1467 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1469 if (ir->adress->n_tf_grps > 0)
1471 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1475 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1477 if (ir->adress->n_energy_grps > 0)
1479 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1485 ir->bAdress = FALSE;
1489 if (file_version >= 48)
1491 gmx_fio_do_int(fio, ir->ePull);
1492 if (ir->ePull != epullNO)
1498 do_pull(fio, ir->pull, bRead, file_version);
1503 ir->ePull = epullNO;
1506 /* Enforced rotation */
1507 if (file_version >= 74)
1509 gmx_fio_do_int(fio, ir->bRot);
1510 if (ir->bRot == TRUE)
1516 do_rot(fio, ir->rot, bRead);
1525 gmx_fio_do_int(fio, ir->opts.ngtc);
1526 if (file_version >= 69)
1528 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1532 ir->opts.nhchainlength = 1;
1534 gmx_fio_do_int(fio, ir->opts.ngacc);
1535 gmx_fio_do_int(fio, ir->opts.ngfrz);
1536 gmx_fio_do_int(fio, ir->opts.ngener);
1540 snew(ir->opts.nrdf, ir->opts.ngtc);
1541 snew(ir->opts.ref_t, ir->opts.ngtc);
1542 snew(ir->opts.annealing, ir->opts.ngtc);
1543 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1544 snew(ir->opts.anneal_time, ir->opts.ngtc);
1545 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1546 snew(ir->opts.tau_t, ir->opts.ngtc);
1547 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1548 snew(ir->opts.acc, ir->opts.ngacc);
1549 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1551 if (ir->opts.ngtc > 0)
1553 if (bRead && file_version < 13)
1555 snew(tmp, ir->opts.ngtc);
1556 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1557 for (i = 0; i < ir->opts.ngtc; i++)
1559 ir->opts.nrdf[i] = tmp[i];
1565 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1567 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1568 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1569 if (file_version < 33 && ir->eI == eiBD)
1571 for (i = 0; i < ir->opts.ngtc; i++)
1573 ir->opts.tau_t[i] = bd_temp;
1577 if (ir->opts.ngfrz > 0)
1579 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1581 if (ir->opts.ngacc > 0)
1583 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1585 if (file_version >= 12)
1587 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1588 ir->opts.ngener*ir->opts.ngener);
1591 if (bRead && file_version < 26)
1593 for (i = 0; i < ir->opts.ngtc; i++)
1597 ir->opts.annealing[i] = eannSINGLE;
1598 ir->opts.anneal_npoints[i] = 2;
1599 snew(ir->opts.anneal_time[i], 2);
1600 snew(ir->opts.anneal_temp[i], 2);
1601 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1602 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1603 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1604 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1605 ir->opts.anneal_time[i][0] = ir->init_t;
1606 ir->opts.anneal_time[i][1] = finish_t;
1607 ir->opts.anneal_temp[i][0] = init_temp;
1608 ir->opts.anneal_temp[i][1] = finish_temp;
1612 ir->opts.annealing[i] = eannNO;
1613 ir->opts.anneal_npoints[i] = 0;
1619 /* file version 26 or later */
1620 /* First read the lists with annealing and npoints for each group */
1621 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1622 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1623 for (j = 0; j < (ir->opts.ngtc); j++)
1625 k = ir->opts.anneal_npoints[j];
1628 snew(ir->opts.anneal_time[j], k);
1629 snew(ir->opts.anneal_temp[j], k);
1631 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1632 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1636 if (file_version >= 45)
1638 gmx_fio_do_int(fio, ir->nwall);
1639 gmx_fio_do_int(fio, ir->wall_type);
1640 if (file_version >= 50)
1642 gmx_fio_do_real(fio, ir->wall_r_linpot);
1646 ir->wall_r_linpot = -1;
1648 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1649 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1650 gmx_fio_do_real(fio, ir->wall_density[0]);
1651 gmx_fio_do_real(fio, ir->wall_density[1]);
1652 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1658 ir->wall_atomtype[0] = -1;
1659 ir->wall_atomtype[1] = -1;
1660 ir->wall_density[0] = 0;
1661 ir->wall_density[1] = 0;
1662 ir->wall_ewald_zfac = 3;
1664 /* Cosine stuff for electric fields */
1665 for (j = 0; (j < DIM); j++)
1667 gmx_fio_do_int(fio, ir->ex[j].n);
1668 gmx_fio_do_int(fio, ir->et[j].n);
1671 snew(ir->ex[j].a, ir->ex[j].n);
1672 snew(ir->ex[j].phi, ir->ex[j].n);
1673 snew(ir->et[j].a, ir->et[j].n);
1674 snew(ir->et[j].phi, ir->et[j].n);
1676 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1677 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1678 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1679 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1683 if (file_version >= 96)
1685 gmx_fio_do_int(fio, ir->eSwapCoords);
1686 if (ir->eSwapCoords != eswapNO)
1692 do_swapcoords(fio, ir->swap, bRead);
1697 if (file_version >= 39)
1699 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1700 gmx_fio_do_int(fio, ir->QMMMscheme);
1701 gmx_fio_do_real(fio, ir->scalefactor);
1702 gmx_fio_do_int(fio, ir->opts.ngQM);
1705 snew(ir->opts.QMmethod, ir->opts.ngQM);
1706 snew(ir->opts.QMbasis, ir->opts.ngQM);
1707 snew(ir->opts.QMcharge, ir->opts.ngQM);
1708 snew(ir->opts.QMmult, ir->opts.ngQM);
1709 snew(ir->opts.bSH, ir->opts.ngQM);
1710 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1711 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1712 snew(ir->opts.SAon, ir->opts.ngQM);
1713 snew(ir->opts.SAoff, ir->opts.ngQM);
1714 snew(ir->opts.SAsteps, ir->opts.ngQM);
1715 snew(ir->opts.bOPT, ir->opts.ngQM);
1716 snew(ir->opts.bTS, ir->opts.ngQM);
1718 if (ir->opts.ngQM > 0)
1720 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1721 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1722 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1723 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1724 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1725 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1726 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1727 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1728 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1729 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1730 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1731 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1733 /* end of QMMM stuff */
1738 static void do_harm(t_fileio *fio, t_iparams *iparams)
1740 gmx_fio_do_real(fio, iparams->harmonic.rA);
1741 gmx_fio_do_real(fio, iparams->harmonic.krA);
1742 gmx_fio_do_real(fio, iparams->harmonic.rB);
1743 gmx_fio_do_real(fio, iparams->harmonic.krB);
1746 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1747 gmx_bool bRead, int file_version)
1754 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1764 do_harm(fio, iparams);
1765 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1767 /* Correct incorrect storage of parameters */
1768 iparams->pdihs.phiB = iparams->pdihs.phiA;
1769 iparams->pdihs.cpB = iparams->pdihs.cpA;
1772 case F_LINEAR_ANGLES:
1773 gmx_fio_do_real(fio, iparams->linangle.klinA);
1774 gmx_fio_do_real(fio, iparams->linangle.aA);
1775 gmx_fio_do_real(fio, iparams->linangle.klinB);
1776 gmx_fio_do_real(fio, iparams->linangle.aB);
1779 gmx_fio_do_real(fio, iparams->fene.bm);
1780 gmx_fio_do_real(fio, iparams->fene.kb);
1783 gmx_fio_do_real(fio, iparams->restraint.lowA);
1784 gmx_fio_do_real(fio, iparams->restraint.up1A);
1785 gmx_fio_do_real(fio, iparams->restraint.up2A);
1786 gmx_fio_do_real(fio, iparams->restraint.kA);
1787 gmx_fio_do_real(fio, iparams->restraint.lowB);
1788 gmx_fio_do_real(fio, iparams->restraint.up1B);
1789 gmx_fio_do_real(fio, iparams->restraint.up2B);
1790 gmx_fio_do_real(fio, iparams->restraint.kB);
1796 gmx_fio_do_real(fio, iparams->tab.kA);
1797 gmx_fio_do_int(fio, iparams->tab.table);
1798 gmx_fio_do_real(fio, iparams->tab.kB);
1800 case F_CROSS_BOND_BONDS:
1801 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1802 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1803 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1805 case F_CROSS_BOND_ANGLES:
1806 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1807 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1808 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1809 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1811 case F_UREY_BRADLEY:
1812 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1813 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1814 gmx_fio_do_real(fio, iparams->u_b.r13A);
1815 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1816 if (file_version >= 79)
1818 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1819 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1820 gmx_fio_do_real(fio, iparams->u_b.r13B);
1821 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1825 iparams->u_b.thetaB = iparams->u_b.thetaA;
1826 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1827 iparams->u_b.r13B = iparams->u_b.r13A;
1828 iparams->u_b.kUBB = iparams->u_b.kUBA;
1831 case F_QUARTIC_ANGLES:
1832 gmx_fio_do_real(fio, iparams->qangle.theta);
1833 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1836 gmx_fio_do_real(fio, iparams->bham.a);
1837 gmx_fio_do_real(fio, iparams->bham.b);
1838 gmx_fio_do_real(fio, iparams->bham.c);
1841 gmx_fio_do_real(fio, iparams->morse.b0A);
1842 gmx_fio_do_real(fio, iparams->morse.cbA);
1843 gmx_fio_do_real(fio, iparams->morse.betaA);
1844 if (file_version >= 79)
1846 gmx_fio_do_real(fio, iparams->morse.b0B);
1847 gmx_fio_do_real(fio, iparams->morse.cbB);
1848 gmx_fio_do_real(fio, iparams->morse.betaB);
1852 iparams->morse.b0B = iparams->morse.b0A;
1853 iparams->morse.cbB = iparams->morse.cbA;
1854 iparams->morse.betaB = iparams->morse.betaA;
1858 gmx_fio_do_real(fio, iparams->cubic.b0);
1859 gmx_fio_do_real(fio, iparams->cubic.kb);
1860 gmx_fio_do_real(fio, iparams->cubic.kcub);
1864 case F_POLARIZATION:
1865 gmx_fio_do_real(fio, iparams->polarize.alpha);
1868 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1869 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1870 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1873 if (file_version < 31)
1875 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1877 gmx_fio_do_real(fio, iparams->wpol.al_x);
1878 gmx_fio_do_real(fio, iparams->wpol.al_y);
1879 gmx_fio_do_real(fio, iparams->wpol.al_z);
1880 gmx_fio_do_real(fio, iparams->wpol.rOH);
1881 gmx_fio_do_real(fio, iparams->wpol.rHH);
1882 gmx_fio_do_real(fio, iparams->wpol.rOD);
1885 gmx_fio_do_real(fio, iparams->thole.a);
1886 gmx_fio_do_real(fio, iparams->thole.alpha1);
1887 gmx_fio_do_real(fio, iparams->thole.alpha2);
1888 gmx_fio_do_real(fio, iparams->thole.rfac);
1891 gmx_fio_do_real(fio, iparams->lj.c6);
1892 gmx_fio_do_real(fio, iparams->lj.c12);
1895 gmx_fio_do_real(fio, iparams->lj14.c6A);
1896 gmx_fio_do_real(fio, iparams->lj14.c12A);
1897 gmx_fio_do_real(fio, iparams->lj14.c6B);
1898 gmx_fio_do_real(fio, iparams->lj14.c12B);
1901 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1902 gmx_fio_do_real(fio, iparams->ljc14.qi);
1903 gmx_fio_do_real(fio, iparams->ljc14.qj);
1904 gmx_fio_do_real(fio, iparams->ljc14.c6);
1905 gmx_fio_do_real(fio, iparams->ljc14.c12);
1907 case F_LJC_PAIRS_NB:
1908 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1909 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1910 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1911 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1917 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1918 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1919 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1921 /* Read the incorrectly stored multiplicity */
1922 gmx_fio_do_real(fio, iparams->harmonic.rB);
1923 gmx_fio_do_real(fio, iparams->harmonic.krB);
1924 iparams->pdihs.phiB = iparams->pdihs.phiA;
1925 iparams->pdihs.cpB = iparams->pdihs.cpA;
1929 gmx_fio_do_real(fio, iparams->pdihs.phiB);
1930 gmx_fio_do_real(fio, iparams->pdihs.cpB);
1931 gmx_fio_do_int(fio, iparams->pdihs.mult);
1935 gmx_fio_do_int(fio, iparams->disres.label);
1936 gmx_fio_do_int(fio, iparams->disres.type);
1937 gmx_fio_do_real(fio, iparams->disres.low);
1938 gmx_fio_do_real(fio, iparams->disres.up1);
1939 gmx_fio_do_real(fio, iparams->disres.up2);
1940 gmx_fio_do_real(fio, iparams->disres.kfac);
1943 gmx_fio_do_int(fio, iparams->orires.ex);
1944 gmx_fio_do_int(fio, iparams->orires.label);
1945 gmx_fio_do_int(fio, iparams->orires.power);
1946 gmx_fio_do_real(fio, iparams->orires.c);
1947 gmx_fio_do_real(fio, iparams->orires.obs);
1948 gmx_fio_do_real(fio, iparams->orires.kfac);
1951 if (file_version < 82)
1953 gmx_fio_do_int(fio, idum);
1954 gmx_fio_do_int(fio, idum);
1956 gmx_fio_do_real(fio, iparams->dihres.phiA);
1957 gmx_fio_do_real(fio, iparams->dihres.dphiA);
1958 gmx_fio_do_real(fio, iparams->dihres.kfacA);
1959 if (file_version >= 82)
1961 gmx_fio_do_real(fio, iparams->dihres.phiB);
1962 gmx_fio_do_real(fio, iparams->dihres.dphiB);
1963 gmx_fio_do_real(fio, iparams->dihres.kfacB);
1967 iparams->dihres.phiB = iparams->dihres.phiA;
1968 iparams->dihres.dphiB = iparams->dihres.dphiA;
1969 iparams->dihres.kfacB = iparams->dihres.kfacA;
1973 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
1974 gmx_fio_do_rvec(fio, iparams->posres.fcA);
1975 if (bRead && file_version < 27)
1977 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
1978 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
1982 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
1983 gmx_fio_do_rvec(fio, iparams->posres.fcB);
1987 gmx_fio_do_int(fio, iparams->fbposres.geom);
1988 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
1989 gmx_fio_do_real(fio, iparams->fbposres.r);
1990 gmx_fio_do_real(fio, iparams->fbposres.k);
1993 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1994 if (file_version >= 25)
1996 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2000 /* Fourier dihedrals are internally represented
2001 * as Ryckaert-Bellemans since those are faster to compute.
2003 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
2004 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
2008 gmx_fio_do_real(fio, iparams->constr.dA);
2009 gmx_fio_do_real(fio, iparams->constr.dB);
2012 gmx_fio_do_real(fio, iparams->settle.doh);
2013 gmx_fio_do_real(fio, iparams->settle.dhh);
2016 gmx_fio_do_real(fio, iparams->vsite.a);
2021 gmx_fio_do_real(fio, iparams->vsite.a);
2022 gmx_fio_do_real(fio, iparams->vsite.b);
2027 gmx_fio_do_real(fio, iparams->vsite.a);
2028 gmx_fio_do_real(fio, iparams->vsite.b);
2029 gmx_fio_do_real(fio, iparams->vsite.c);
2032 gmx_fio_do_int(fio, iparams->vsiten.n);
2033 gmx_fio_do_real(fio, iparams->vsiten.a);
2038 /* We got rid of some parameters in version 68 */
2039 if (bRead && file_version < 68)
2041 gmx_fio_do_real(fio, rdum);
2042 gmx_fio_do_real(fio, rdum);
2043 gmx_fio_do_real(fio, rdum);
2044 gmx_fio_do_real(fio, rdum);
2046 gmx_fio_do_real(fio, iparams->gb.sar);
2047 gmx_fio_do_real(fio, iparams->gb.st);
2048 gmx_fio_do_real(fio, iparams->gb.pi);
2049 gmx_fio_do_real(fio, iparams->gb.gbr);
2050 gmx_fio_do_real(fio, iparams->gb.bmlt);
2053 gmx_fio_do_int(fio, iparams->cmap.cmapA);
2054 gmx_fio_do_int(fio, iparams->cmap.cmapB);
2057 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
2058 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
2062 gmx_fio_unset_comment(fio);
2066 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
2073 gmx_fio_set_comment(fio, interaction_function[ftype].name);
2075 if (file_version < 44)
2077 for (i = 0; i < MAXNODES; i++)
2079 gmx_fio_do_int(fio, idum);
2082 gmx_fio_do_int(fio, ilist->nr);
2085 snew(ilist->iatoms, ilist->nr);
2087 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2090 gmx_fio_unset_comment(fio);
2094 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2095 gmx_bool bRead, int file_version)
2100 gmx_fio_do_int(fio, ffparams->atnr);
2101 if (file_version < 57)
2103 gmx_fio_do_int(fio, idum);
2105 gmx_fio_do_int(fio, ffparams->ntypes);
2108 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2109 ffparams->atnr, ffparams->ntypes);
2113 snew(ffparams->functype, ffparams->ntypes);
2114 snew(ffparams->iparams, ffparams->ntypes);
2116 /* Read/write all the function types */
2117 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2120 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2123 if (file_version >= 66)
2125 gmx_fio_do_double(fio, ffparams->reppow);
2129 ffparams->reppow = 12.0;
2132 if (file_version >= 57)
2134 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2137 /* Check whether all these function types are supported by the code.
2138 * In practice the code is backwards compatible, which means that the
2139 * numbering may have to be altered from old numbering to new numbering
2141 for (i = 0; (i < ffparams->ntypes); i++)
2145 /* Loop over file versions */
2146 for (k = 0; (k < NFTUPD); k++)
2148 /* Compare the read file_version to the update table */
2149 if ((file_version < ftupd[k].fvnr) &&
2150 (ffparams->functype[i] >= ftupd[k].ftype))
2152 ffparams->functype[i] += 1;
2155 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2156 i, ffparams->functype[i],
2157 interaction_function[ftupd[k].ftype].longname);
2164 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2168 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2173 static void add_settle_atoms(t_ilist *ilist)
2177 /* Settle used to only store the first atom: add the other two */
2178 srenew(ilist->iatoms, 2*ilist->nr);
2179 for (i = ilist->nr/2-1; i >= 0; i--)
2181 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2182 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2183 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2184 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2186 ilist->nr = 2*ilist->nr;
2189 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2192 int i, j, renum[F_NRE];
2196 for (j = 0; (j < F_NRE); j++)
2201 for (k = 0; k < NFTUPD; k++)
2203 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2212 ilist[j].iatoms = NULL;
2216 do_ilist(fio, &ilist[j], bRead, file_version, j);
2217 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2219 add_settle_atoms(&ilist[j]);
2223 if (bRead && gmx_debug_at)
2224 pr_ilist(debug,0,interaction_function[j].longname,
2225 functype,&ilist[j],TRUE);
2230 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2231 gmx_bool bRead, int file_version)
2233 do_ffparams(fio, ffparams, bRead, file_version);
2235 if (file_version >= 54)
2237 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2240 do_ilists(fio, molt->ilist, bRead, file_version);
2243 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2245 int i, idum, dum_nra, *dum_a;
2247 if (file_version < 44)
2249 for (i = 0; i < MAXNODES; i++)
2251 gmx_fio_do_int(fio, idum);
2254 gmx_fio_do_int(fio, block->nr);
2255 if (file_version < 51)
2257 gmx_fio_do_int(fio, dum_nra);
2261 if ((block->nalloc_index > 0) && (NULL != block->index))
2263 sfree(block->index);
2265 block->nalloc_index = block->nr+1;
2266 snew(block->index, block->nalloc_index);
2268 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2270 if (file_version < 51 && dum_nra > 0)
2272 snew(dum_a, dum_nra);
2273 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2278 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2283 if (file_version < 44)
2285 for (i = 0; i < MAXNODES; i++)
2287 gmx_fio_do_int(fio, idum);
2290 gmx_fio_do_int(fio, block->nr);
2291 gmx_fio_do_int(fio, block->nra);
2294 block->nalloc_index = block->nr+1;
2295 snew(block->index, block->nalloc_index);
2296 block->nalloc_a = block->nra;
2297 snew(block->a, block->nalloc_a);
2299 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2300 gmx_fio_ndo_int(fio, block->a, block->nra);
2303 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2304 int file_version, gmx_groups_t *groups, int atnr)
2308 gmx_fio_do_real(fio, atom->m);
2309 gmx_fio_do_real(fio, atom->q);
2310 gmx_fio_do_real(fio, atom->mB);
2311 gmx_fio_do_real(fio, atom->qB);
2312 gmx_fio_do_ushort(fio, atom->type);
2313 gmx_fio_do_ushort(fio, atom->typeB);
2314 gmx_fio_do_int(fio, atom->ptype);
2315 gmx_fio_do_int(fio, atom->resind);
2316 if (file_version >= 52)
2318 gmx_fio_do_int(fio, atom->atomnumber);
2322 atom->atomnumber = NOTSET;
2324 if (file_version < 23)
2328 else if (file_version < 39)
2337 if (file_version < 57)
2339 unsigned char uchar[egcNR];
2340 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2341 for (i = myngrp; (i < ngrp); i++)
2345 /* Copy the old data format to the groups struct */
2346 for (i = 0; i < ngrp; i++)
2348 groups->grpnr[i][atnr] = uchar[i];
2353 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2358 if (file_version < 23)
2362 else if (file_version < 39)
2371 for (j = 0; (j < ngrp); j++)
2375 gmx_fio_do_int(fio, grps[j].nr);
2378 snew(grps[j].nm_ind, grps[j].nr);
2380 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2385 snew(grps[j].nm_ind, grps[j].nr);
2390 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2396 gmx_fio_do_int(fio, ls);
2397 *nm = get_symtab_handle(symtab, ls);
2401 ls = lookup_symtab(symtab, *nm);
2402 gmx_fio_do_int(fio, ls);
2406 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2411 for (j = 0; (j < nstr); j++)
2413 do_symstr(fio, &(nm[j]), bRead, symtab);
2417 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2418 t_symtab *symtab, int file_version)
2422 for (j = 0; (j < n); j++)
2424 do_symstr(fio, &(ri[j].name), bRead, symtab);
2425 if (file_version >= 63)
2427 gmx_fio_do_int(fio, ri[j].nr);
2428 gmx_fio_do_uchar(fio, ri[j].ic);
2438 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2440 gmx_groups_t *groups)
2444 gmx_fio_do_int(fio, atoms->nr);
2445 gmx_fio_do_int(fio, atoms->nres);
2446 if (file_version < 57)
2448 gmx_fio_do_int(fio, groups->ngrpname);
2449 for (i = 0; i < egcNR; i++)
2451 groups->ngrpnr[i] = atoms->nr;
2452 snew(groups->grpnr[i], groups->ngrpnr[i]);
2457 snew(atoms->atom, atoms->nr);
2458 snew(atoms->atomname, atoms->nr);
2459 snew(atoms->atomtype, atoms->nr);
2460 snew(atoms->atomtypeB, atoms->nr);
2461 snew(atoms->resinfo, atoms->nres);
2462 if (file_version < 57)
2464 snew(groups->grpname, groups->ngrpname);
2466 atoms->pdbinfo = NULL;
2468 for (i = 0; (i < atoms->nr); i++)
2470 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2472 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2473 if (bRead && (file_version <= 20))
2475 for (i = 0; i < atoms->nr; i++)
2477 atoms->atomtype[i] = put_symtab(symtab, "?");
2478 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2483 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2484 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2486 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2488 if (file_version < 57)
2490 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2492 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2496 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2497 gmx_bool bRead, t_symtab *symtab,
2502 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2503 gmx_fio_do_int(fio, groups->ngrpname);
2506 snew(groups->grpname, groups->ngrpname);
2508 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2509 for (g = 0; g < egcNR; g++)
2511 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2512 if (groups->ngrpnr[g] == 0)
2516 groups->grpnr[g] = NULL;
2523 snew(groups->grpnr[g], groups->ngrpnr[g]);
2525 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2530 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2535 if (file_version > 25)
2537 gmx_fio_do_int(fio, atomtypes->nr);
2541 snew(atomtypes->radius, j);
2542 snew(atomtypes->vol, j);
2543 snew(atomtypes->surftens, j);
2544 snew(atomtypes->atomnumber, j);
2545 snew(atomtypes->gb_radius, j);
2546 snew(atomtypes->S_hct, j);
2548 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2549 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2550 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2551 if (file_version >= 40)
2553 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2555 if (file_version >= 60)
2557 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2558 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2563 /* File versions prior to 26 cannot do GBSA,
2564 * so they dont use this structure
2567 atomtypes->radius = NULL;
2568 atomtypes->vol = NULL;
2569 atomtypes->surftens = NULL;
2570 atomtypes->atomnumber = NULL;
2571 atomtypes->gb_radius = NULL;
2572 atomtypes->S_hct = NULL;
2576 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2582 gmx_fio_do_int(fio, symtab->nr);
2586 snew(symtab->symbuf, 1);
2587 symbuf = symtab->symbuf;
2588 symbuf->bufsize = nr;
2589 snew(symbuf->buf, nr);
2590 for (i = 0; (i < nr); i++)
2592 gmx_fio_do_string(fio, buf);
2593 symbuf->buf[i] = strdup(buf);
2598 symbuf = symtab->symbuf;
2599 while (symbuf != NULL)
2601 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2603 gmx_fio_do_string(fio, symbuf->buf[i]);
2606 symbuf = symbuf->next;
2610 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2615 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2617 int i, j, ngrid, gs, nelem;
2619 gmx_fio_do_int(fio, cmap_grid->ngrid);
2620 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2622 ngrid = cmap_grid->ngrid;
2623 gs = cmap_grid->grid_spacing;
2628 snew(cmap_grid->cmapdata, ngrid);
2630 for (i = 0; i < cmap_grid->ngrid; i++)
2632 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2636 for (i = 0; i < cmap_grid->ngrid; i++)
2638 for (j = 0; j < nelem; j++)
2640 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2641 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2642 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2643 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2649 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2651 int m, a, a0, a1, r;
2655 /* We always assign a new chain number, but save the chain id characters
2656 * for larger molecules.
2658 #define CHAIN_MIN_ATOMS 15
2662 for (m = 0; m < mols->nr; m++)
2664 a0 = mols->index[m];
2665 a1 = mols->index[m+1];
2666 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2675 for (a = a0; a < a1; a++)
2677 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2678 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2683 /* Blank out the chain id if there was only one chain */
2686 for (r = 0; r < atoms->nres; r++)
2688 atoms->resinfo[r].chainid = ' ';
2693 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2694 t_symtab *symtab, int file_version,
2695 gmx_groups_t *groups)
2699 if (file_version >= 57)
2701 do_symstr(fio, &(molt->name), bRead, symtab);
2704 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2706 if (bRead && gmx_debug_at)
2708 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2711 if (file_version >= 57)
2713 do_ilists(fio, molt->ilist, bRead, file_version);
2715 do_block(fio, &molt->cgs, bRead, file_version);
2716 if (bRead && gmx_debug_at)
2718 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2722 /* This used to be in the atoms struct */
2723 do_blocka(fio, &molt->excls, bRead, file_version);
2726 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2730 gmx_fio_do_int(fio, molb->type);
2731 gmx_fio_do_int(fio, molb->nmol);
2732 gmx_fio_do_int(fio, molb->natoms_mol);
2733 /* Position restraint coordinates */
2734 gmx_fio_do_int(fio, molb->nposres_xA);
2735 if (molb->nposres_xA > 0)
2739 snew(molb->posres_xA, molb->nposres_xA);
2741 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2743 gmx_fio_do_int(fio, molb->nposres_xB);
2744 if (molb->nposres_xB > 0)
2748 snew(molb->posres_xB, molb->nposres_xB);
2750 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2755 static t_block mtop_mols(gmx_mtop_t *mtop)
2761 for (mb = 0; mb < mtop->nmolblock; mb++)
2763 mols.nr += mtop->molblock[mb].nmol;
2765 mols.nalloc_index = mols.nr + 1;
2766 snew(mols.index, mols.nalloc_index);
2771 for (mb = 0; mb < mtop->nmolblock; mb++)
2773 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2775 a += mtop->molblock[mb].natoms_mol;
2784 static void add_posres_molblock(gmx_mtop_t *mtop)
2789 gmx_molblock_t *molb;
2792 /* posres reference positions are stored in ip->posres (if present) and
2793 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2794 posres.pos0A are identical to fbposres.pos0. */
2795 il = &mtop->moltype[0].ilist[F_POSRES];
2796 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2797 if (il->nr == 0 && ilfb->nr == 0)
2803 for (i = 0; i < il->nr; i += 2)
2805 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2806 am = max(am, il->iatoms[i+1]);
2807 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2808 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2809 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2814 /* This loop is required if we have only flat-bottomed posres:
2816 - bFE == FALSE (no B-state for flat-bottomed posres) */
2819 for (i = 0; i < ilfb->nr; i += 2)
2821 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2822 am = max(am, ilfb->iatoms[i+1]);
2825 /* Make the posres coordinate block end at a molecule end */
2827 while (am >= mtop->mols.index[mol+1])
2831 molb = &mtop->molblock[0];
2832 molb->nposres_xA = mtop->mols.index[mol+1];
2833 snew(molb->posres_xA, molb->nposres_xA);
2836 molb->nposres_xB = molb->nposres_xA;
2837 snew(molb->posres_xB, molb->nposres_xB);
2841 molb->nposres_xB = 0;
2843 for (i = 0; i < il->nr; i += 2)
2845 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2846 a = il->iatoms[i+1];
2847 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2848 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2849 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2852 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2853 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2854 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2859 /* If only flat-bottomed posres are present, take reference pos from them.
2860 Here: bFE == FALSE */
2861 for (i = 0; i < ilfb->nr; i += 2)
2863 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2864 a = ilfb->iatoms[i+1];
2865 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2866 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2867 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2872 static void set_disres_npair(gmx_mtop_t *mtop)
2879 ip = mtop->ffparams.iparams;
2881 for (mt = 0; mt < mtop->nmoltype; mt++)
2883 il = &mtop->moltype[mt].ilist[F_DISRES];
2888 for (i = 0; i < il->nr; i += 3)
2891 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2893 ip[a[i]].disres.npair = npair;
2901 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2911 do_symtab(fio, &(mtop->symtab), bRead);
2914 pr_symtab(debug, 0, "symtab", &mtop->symtab);
2917 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2919 if (file_version >= 57)
2921 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
2923 gmx_fio_do_int(fio, mtop->nmoltype);
2931 snew(mtop->moltype, mtop->nmoltype);
2932 if (file_version < 57)
2934 mtop->moltype[0].name = mtop->name;
2937 for (mt = 0; mt < mtop->nmoltype; mt++)
2939 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
2943 if (file_version >= 57)
2945 gmx_fio_do_int(fio, mtop->nmolblock);
2949 mtop->nmolblock = 1;
2953 snew(mtop->molblock, mtop->nmolblock);
2955 if (file_version >= 57)
2957 for (mb = 0; mb < mtop->nmolblock; mb++)
2959 do_molblock(fio, &mtop->molblock[mb], bRead);
2961 gmx_fio_do_int(fio, mtop->natoms);
2965 mtop->molblock[0].type = 0;
2966 mtop->molblock[0].nmol = 1;
2967 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
2968 mtop->molblock[0].nposres_xA = 0;
2969 mtop->molblock[0].nposres_xB = 0;
2972 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
2975 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
2978 if (file_version < 57)
2980 /* Debug statements are inside do_idef */
2981 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
2982 mtop->natoms = mtop->moltype[0].atoms.nr;
2985 if (file_version >= 65)
2987 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
2991 mtop->ffparams.cmap_grid.ngrid = 0;
2992 mtop->ffparams.cmap_grid.grid_spacing = 0;
2993 mtop->ffparams.cmap_grid.cmapdata = NULL;
2996 if (file_version >= 57)
2998 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
3001 if (file_version < 57)
3003 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
3004 if (bRead && gmx_debug_at)
3006 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
3008 do_block(fio, &mtop->mols, bRead, file_version);
3009 /* Add the posres coordinates to the molblock */
3010 add_posres_molblock(mtop);
3014 if (file_version >= 57)
3016 done_block(&mtop->mols);
3017 mtop->mols = mtop_mols(mtop);
3021 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
3025 if (file_version < 51)
3027 /* Here used to be the shake blocks */
3028 do_blocka(fio, &dumb, bRead, file_version);
3041 close_symtab(&(mtop->symtab));
3045 /* If TopOnlyOK is TRUE then we can read even future versions
3046 * of tpx files, provided the file_generation hasn't changed.
3047 * If it is FALSE, we need the inputrecord too, and bail out
3048 * if the file is newer than the program.
3050 * The version and generation if the topology (see top of this file)
3051 * are returned in the two last arguments.
3053 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
3055 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
3056 gmx_bool TopOnlyOK, int *file_version,
3057 int *file_generation)
3060 char file_tag[STRLEN];
3067 gmx_fio_checktype(fio);
3068 gmx_fio_setdebug(fio, bDebugMode());
3070 /* NEW! XDR tpb file */
3071 precision = sizeof(real);
3074 gmx_fio_do_string(fio, buf);
3075 if (strncmp(buf, "VERSION", 7))
3077 gmx_fatal(FARGS, "Can not read file %s,\n"
3078 " this file is from a Gromacs version which is older than 2.0\n"
3079 " Make a new one with grompp or use a gro or pdb file, if possible",
3080 gmx_fio_getname(fio));
3082 gmx_fio_do_int(fio, precision);
3083 bDouble = (precision == sizeof(double));
3084 if ((precision != sizeof(float)) && !bDouble)
3086 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3087 "instead of %d or %d",
3088 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3090 gmx_fio_setprecision(fio, bDouble);
3091 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3092 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3096 gmx_fio_write_string(fio, GromacsVersion());
3097 bDouble = (precision == sizeof(double));
3098 gmx_fio_setprecision(fio, bDouble);
3099 gmx_fio_do_int(fio, precision);
3101 sprintf(file_tag, "%s", tpx_tag);
3102 fgen = tpx_generation;
3105 /* Check versions! */
3106 gmx_fio_do_int(fio, fver);
3108 /* This is for backward compatibility with development versions 77-79
3109 * where the tag was, mistakenly, placed before the generation,
3110 * which would cause a segv instead of a proper error message
3111 * when reading the topology only from tpx with <77 code.
3113 if (fver >= 77 && fver <= 79)
3115 gmx_fio_do_string(fio, file_tag);
3120 gmx_fio_do_int(fio, fgen);
3129 gmx_fio_do_string(fio, file_tag);
3135 /* Versions before 77 don't have the tag, set it to release */
3136 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3139 if (strcmp(file_tag, tpx_tag) != 0)
3141 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3144 /* We only support reading tpx files with the same tag as the code
3145 * or tpx files with the release tag and with lower version number.
3147 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
3149 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3150 gmx_fio_getname(fio), fver, file_tag,
3151 tpx_version, tpx_tag);
3156 if (file_version != NULL)
3158 *file_version = fver;
3160 if (file_generation != NULL)
3162 *file_generation = fgen;
3166 if ((fver <= tpx_incompatible_version) ||
3167 ((fver > tpx_version) && !TopOnlyOK) ||
3168 (fgen > tpx_generation) ||
3169 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3171 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3172 gmx_fio_getname(fio), fver, tpx_version);
3175 do_section(fio, eitemHEADER, bRead);
3176 gmx_fio_do_int(fio, tpx->natoms);
3179 gmx_fio_do_int(fio, tpx->ngtc);
3187 gmx_fio_do_int(fio, idum);
3188 gmx_fio_do_real(fio, rdum);
3190 /*a better decision will eventually (5.0 or later) need to be made
3191 on how to treat the alchemical state of the system, which can now
3192 vary through a simulation, and cannot be completely described
3193 though a single lambda variable, or even a single state
3194 index. Eventually, should probably be a vector. MRS*/
3197 gmx_fio_do_int(fio, tpx->fep_state);
3199 gmx_fio_do_real(fio, tpx->lambda);
3200 gmx_fio_do_int(fio, tpx->bIr);
3201 gmx_fio_do_int(fio, tpx->bTop);
3202 gmx_fio_do_int(fio, tpx->bX);
3203 gmx_fio_do_int(fio, tpx->bV);
3204 gmx_fio_do_int(fio, tpx->bF);
3205 gmx_fio_do_int(fio, tpx->bBox);
3207 if ((fgen > tpx_generation))
3209 /* This can only happen if TopOnlyOK=TRUE */
3214 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3215 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3216 gmx_bool bXVallocated)
3222 int file_version, file_generation;
3226 gmx_bool bPeriodicMols;
3230 tpx.natoms = state->natoms;
3231 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3232 tpx.fep_state = state->fep_state;
3233 tpx.lambda = state->lambda[efptFEP];
3234 tpx.bIr = (ir != NULL);
3235 tpx.bTop = (mtop != NULL);
3236 tpx.bX = (state->x != NULL);
3237 tpx.bV = (state->v != NULL);
3238 tpx.bF = (f != NULL);
3242 TopOnlyOK = (ir == NULL);
3244 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3249 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3250 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3255 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3256 state->natoms = tpx.natoms;
3257 state->nalloc = tpx.natoms;
3263 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3267 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3269 do_test(fio, tpx.bBox, state->box);
3270 do_section(fio, eitemBOX, bRead);
3273 gmx_fio_ndo_rvec(fio, state->box, DIM);
3274 if (file_version >= 51)
3276 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3280 /* We initialize box_rel after reading the inputrec */
3281 clear_mat(state->box_rel);
3283 if (file_version >= 28)
3285 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3286 if (file_version < 56)
3289 gmx_fio_ndo_rvec(fio, mdum, DIM);
3294 if (state->ngtc > 0 && file_version >= 28)
3297 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3298 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3299 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3300 snew(dumv, state->ngtc);
3301 if (file_version < 69)
3303 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3305 /* These used to be the Berendsen tcoupl_lambda's */
3306 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3310 /* Prior to tpx version 26, the inputrec was here.
3311 * I moved it to enable partial forward-compatibility
3312 * for analysis/viewer programs.
3314 if (file_version < 26)
3316 do_test(fio, tpx.bIr, ir);
3317 do_section(fio, eitemIR, bRead);
3322 do_inputrec(fio, ir, bRead, file_version,
3323 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3326 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3331 do_inputrec(fio, &dum_ir, bRead, file_version,
3332 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3335 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3337 done_inputrec(&dum_ir);
3343 do_test(fio, tpx.bTop, mtop);
3344 do_section(fio, eitemTOP, bRead);
3349 do_mtop(fio, mtop, bRead, file_version);
3353 do_mtop(fio, &dum_top, bRead, file_version);
3354 done_mtop(&dum_top, TRUE);
3357 do_test(fio, tpx.bX, state->x);
3358 do_section(fio, eitemX, bRead);
3363 state->flags |= (1<<estX);
3365 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3368 do_test(fio, tpx.bV, state->v);
3369 do_section(fio, eitemV, bRead);
3374 state->flags |= (1<<estV);
3376 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3379 do_test(fio, tpx.bF, f);
3380 do_section(fio, eitemF, bRead);
3383 gmx_fio_ndo_rvec(fio, f, state->natoms);
3386 /* Starting with tpx version 26, we have the inputrec
3387 * at the end of the file, so we can ignore it
3388 * if the file is never than the software (but still the
3389 * same generation - see comments at the top of this file.
3394 bPeriodicMols = FALSE;
3395 if (file_version >= 26)
3397 do_test(fio, tpx.bIr, ir);
3398 do_section(fio, eitemIR, bRead);
3401 if (file_version >= 53)
3403 /* Removed the pbc info from do_inputrec, since we always want it */
3407 bPeriodicMols = ir->bPeriodicMols;
3409 gmx_fio_do_int(fio, ePBC);
3410 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3412 if (file_generation <= tpx_generation && ir)
3414 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3417 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3419 if (file_version < 51)
3421 set_box_rel(ir, state);
3423 if (file_version < 53)
3426 bPeriodicMols = ir->bPeriodicMols;
3429 if (bRead && ir && file_version >= 53)
3431 /* We need to do this after do_inputrec, since that initializes ir */
3433 ir->bPeriodicMols = bPeriodicMols;
3442 if (state->ngtc == 0)
3444 /* Reading old version without tcoupl state data: set it */
3445 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3447 if (tpx.bTop && mtop)
3449 if (file_version < 57)
3451 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3453 ir->eDisre = edrSimple;
3457 ir->eDisre = edrNone;
3460 set_disres_npair(mtop);
3464 if (tpx.bTop && mtop)
3466 gmx_mtop_finalize(mtop);
3469 if (file_version >= 57)
3473 env = getenv("GMX_NOCHARGEGROUPS");
3476 sscanf(env, "%d", &ienv);
3477 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3482 "Will make single atomic charge groups in non-solvent%s\n",
3483 ienv > 1 ? " and solvent" : "");
3484 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3486 fprintf(stderr, "\n");
3494 /************************************************************
3496 * The following routines are the exported ones
3498 ************************************************************/
3500 t_fileio *open_tpx(const char *fn, const char *mode)
3502 return gmx_fio_open(fn, mode);
3505 void close_tpx(t_fileio *fio)
3510 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3511 int *file_version, int *file_generation)
3515 fio = open_tpx(fn, "r");
3516 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3520 void write_tpx_state(const char *fn,
3521 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3525 fio = open_tpx(fn, "w");
3526 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3530 void read_tpx_state(const char *fn,
3531 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3535 fio = open_tpx(fn, "r");
3536 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3540 int read_tpx(const char *fn,
3541 t_inputrec *ir, matrix box, int *natoms,
3542 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3550 fio = open_tpx(fn, "r");
3551 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3553 *natoms = state.natoms;
3556 copy_mat(state.box, box);
3565 int read_tpx_top(const char *fn,
3566 t_inputrec *ir, matrix box, int *natoms,
3567 rvec *x, rvec *v, rvec *f, t_topology *top)
3573 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3575 *top = gmx_mtop_t_to_t_topology(&mtop);
3580 gmx_bool fn2bTPX(const char *file)
3582 switch (fn2ftp(file))
3593 static void done_gmx_groups_t(gmx_groups_t *g)
3597 for (i = 0; (i < egcNR); i++)
3599 if (NULL != g->grps[i].nm_ind)
3601 sfree(g->grps[i].nm_ind);
3602 g->grps[i].nm_ind = NULL;
3604 if (NULL != g->grpnr[i])
3610 /* The contents of this array is in symtab, don't free it here */
3614 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3615 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3618 int natoms, i, version, generation;
3619 gmx_bool bTop, bXNULL = FALSE;
3621 t_topology *topconv;
3624 bTop = fn2bTPX(infile);
3628 read_tpxheader(infile, &header, TRUE, &version, &generation);
3631 snew(*x, header.natoms);
3635 snew(*v, header.natoms);
3638 *ePBC = read_tpx(infile, NULL, box, &natoms,
3639 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3640 *top = gmx_mtop_t_to_t_topology(mtop);
3641 /* In this case we need to throw away the group data too */
3642 done_gmx_groups_t(&mtop->groups);
3644 strcpy(title, *top->name);
3645 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3649 get_stx_coordnum(infile, &natoms);
3650 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3661 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3669 aps = gmx_atomprop_init();
3670 for (i = 0; (i < natoms); i++)
3672 if (!gmx_atomprop_query(aps, epropMass,
3673 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3674 *top->atoms.atomname[i],
3675 &(top->atoms.atom[i].m)))
3679 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3680 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3681 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3682 *top->atoms.atomname[i]);
3686 gmx_atomprop_destroy(aps);
3688 top->idef.ntypes = -1;