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, 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"
61 #include "mtop_util.h"
63 #define TPX_TAG_RELEASE "release"
65 /* This is the tag string which is stored in the tpx file.
66 * Change this if you want to change the tpx format in a feature branch.
67 * This ensures that there will not be different tpx formats around which
68 * can not be distinguished.
70 static const char *tpx_tag = TPX_TAG_RELEASE;
72 /* This number should be increased whenever the file format changes! */
73 static const int tpx_version = 95;
75 /* This number should only be increased when you edit the TOPOLOGY section
76 * or the HEADER of the tpx format.
77 * This way we can maintain forward compatibility too for all analysis tools
78 * and/or external programs that only need to know the atom/residue names,
79 * charges, and bond connectivity.
81 * It first appeared in tpx version 26, when I also moved the inputrecord
82 * to the end of the tpx file, so we can just skip it if we only
85 static const int tpx_generation = 25;
87 /* This number should be the most recent backwards incompatible version
88 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
90 static const int tpx_incompatible_version = 9;
94 /* Struct used to maintain tpx compatibility when function types are added */
96 int fvnr; /* file version number in which the function type first appeared */
97 int ftype; /* function type */
101 * The entries should be ordered in:
102 * 1. ascending file version number
103 * 2. ascending function type number
105 /*static const t_ftupd ftupd[] = {
106 { 20, F_CUBICBONDS },
110 { 22, F_DISRESVIOL },
116 { 26, F_DIHRESVIOL },
117 { 30, F_CROSS_BOND_BONDS },
118 { 30, F_CROSS_BOND_ANGLES },
119 { 30, F_UREY_BRADLEY },
120 { 30, F_POLARIZATION },
124 * The entries should be ordered in:
125 * 1. ascending function type number
126 * 2. ascending file version number
128 /* question; what is the purpose of the commented code above? */
129 static const t_ftupd ftupd[] = {
130 { 20, F_CUBICBONDS },
135 { 43, F_TABBONDSNC },
136 { 70, F_RESTRBONDS },
137 { 76, F_LINEAR_ANGLES },
138 { 30, F_CROSS_BOND_BONDS },
139 { 30, F_CROSS_BOND_ANGLES },
140 { 30, F_UREY_BRADLEY },
141 { 34, F_QUARTIC_ANGLES },
151 { 72, F_NPSOLVATION },
153 { 41, F_LJC_PAIRS_NB },
156 { 32, F_COUL_RECIP },
159 { 30, F_POLARIZATION },
162 { 22, F_DISRESVIOL },
166 { 26, F_DIHRESVIOL },
171 { 46, F_ECONSERVED },
172 { 69, F_VTEMP_NOLONGERUSED},
174 { 54, F_DVDL_CONSTR },
175 { 76, F_ANHARM_POL },
178 { 79, F_DVDL_BONDED, },
179 { 79, F_DVDL_RESTRAINT },
180 { 79, F_DVDL_TEMPERATURE },
182 #define NFTUPD asize(ftupd)
184 /* Needed for backward compatibility */
187 static void _do_section(t_fileio *fio, int key, gmx_bool bRead, const char *src,
193 if (gmx_fio_getftp(fio) == efTPA)
197 gmx_fio_write_string(fio, itemstr[key]);
198 bDbg = gmx_fio_getdebug(fio);
199 gmx_fio_setdebug(fio, FALSE);
200 gmx_fio_write_string(fio, comment_str[key]);
201 gmx_fio_setdebug(fio, bDbg);
205 if (gmx_fio_getdebug(fio))
207 fprintf(stderr, "Looking for section %s (%s, %d)",
208 itemstr[key], src, line);
213 gmx_fio_do_string(fio, buf);
215 while ((gmx_strcasecmp(buf, itemstr[key]) != 0));
217 if (gmx_strcasecmp(buf, itemstr[key]) != 0)
219 gmx_fatal(FARGS, "\nCould not find section heading %s", itemstr[key]);
221 else if (gmx_fio_getdebug(fio))
223 fprintf(stderr, " and found it\n");
229 #define do_section(fio, key, bRead) _do_section(fio, key, bRead, __FILE__, __LINE__)
231 /**************************************************************
233 * Now the higer level routines that do io of the structures and arrays
235 **************************************************************/
236 static void do_pullgrp_tpx_pre95(t_fileio *fio,
245 gmx_fio_do_int(fio, pgrp->nat);
248 snew(pgrp->ind, pgrp->nat);
250 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
251 gmx_fio_do_int(fio, pgrp->nweight);
254 snew(pgrp->weight, pgrp->nweight);
256 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
257 gmx_fio_do_int(fio, pgrp->pbcatom);
258 gmx_fio_do_rvec(fio, pcrd->vec);
259 clear_rvec(pcrd->origin);
260 gmx_fio_do_rvec(fio, tmp);
262 gmx_fio_do_real(fio, pcrd->rate);
263 gmx_fio_do_real(fio, pcrd->k);
264 if (file_version >= 56)
266 gmx_fio_do_real(fio, pcrd->kB);
274 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
278 gmx_fio_do_int(fio, pgrp->nat);
281 snew(pgrp->ind, pgrp->nat);
283 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
284 gmx_fio_do_int(fio, pgrp->nweight);
287 snew(pgrp->weight, pgrp->nweight);
289 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
290 gmx_fio_do_int(fio, pgrp->pbcatom);
293 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd)
297 gmx_fio_do_int(fio, pcrd->group[0]);
298 gmx_fio_do_int(fio, pcrd->group[1]);
299 gmx_fio_do_rvec(fio, pcrd->origin);
300 gmx_fio_do_rvec(fio, pcrd->vec);
301 gmx_fio_do_real(fio, pcrd->init);
302 gmx_fio_do_real(fio, pcrd->rate);
303 gmx_fio_do_real(fio, pcrd->k);
304 gmx_fio_do_real(fio, pcrd->kB);
307 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
309 /* i is used in the ndo_double macro*/
313 int n_lambda = fepvals->n_lambda;
315 /* reset the lambda calculation window */
316 fepvals->lambda_start_n = 0;
317 fepvals->lambda_stop_n = n_lambda;
318 if (file_version >= 79)
324 snew(expand->init_lambda_weights, n_lambda);
326 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
327 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
330 gmx_fio_do_int(fio, expand->nstexpanded);
331 gmx_fio_do_int(fio, expand->elmcmove);
332 gmx_fio_do_int(fio, expand->elamstats);
333 gmx_fio_do_int(fio, expand->lmc_repeats);
334 gmx_fio_do_int(fio, expand->gibbsdeltalam);
335 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
336 gmx_fio_do_int(fio, expand->lmc_seed);
337 gmx_fio_do_real(fio, expand->mc_temp);
338 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
339 gmx_fio_do_int(fio, expand->nstTij);
340 gmx_fio_do_int(fio, expand->minvarmin);
341 gmx_fio_do_int(fio, expand->c_range);
342 gmx_fio_do_real(fio, expand->wl_scale);
343 gmx_fio_do_real(fio, expand->wl_ratio);
344 gmx_fio_do_real(fio, expand->init_wl_delta);
345 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
346 gmx_fio_do_int(fio, expand->elmceq);
347 gmx_fio_do_int(fio, expand->equil_steps);
348 gmx_fio_do_int(fio, expand->equil_samples);
349 gmx_fio_do_int(fio, expand->equil_n_at_lam);
350 gmx_fio_do_real(fio, expand->equil_wl_delta);
351 gmx_fio_do_real(fio, expand->equil_ratio);
355 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
358 if (file_version >= 79)
360 gmx_fio_do_int(fio, simtemp->eSimTempScale);
361 gmx_fio_do_real(fio, simtemp->simtemp_high);
362 gmx_fio_do_real(fio, simtemp->simtemp_low);
367 snew(simtemp->temperatures, n_lambda);
369 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
374 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
376 /* i is defined in the ndo_double macro; use g to iterate. */
381 /* free energy values */
383 if (file_version >= 79)
385 gmx_fio_do_int(fio, fepvals->init_fep_state);
386 gmx_fio_do_double(fio, fepvals->init_lambda);
387 gmx_fio_do_double(fio, fepvals->delta_lambda);
389 else if (file_version >= 59)
391 gmx_fio_do_double(fio, fepvals->init_lambda);
392 gmx_fio_do_double(fio, fepvals->delta_lambda);
396 gmx_fio_do_real(fio, rdum);
397 fepvals->init_lambda = rdum;
398 gmx_fio_do_real(fio, rdum);
399 fepvals->delta_lambda = rdum;
401 if (file_version >= 79)
403 gmx_fio_do_int(fio, fepvals->n_lambda);
406 snew(fepvals->all_lambda, efptNR);
408 for (g = 0; g < efptNR; g++)
410 if (fepvals->n_lambda > 0)
414 snew(fepvals->all_lambda[g], fepvals->n_lambda);
416 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
417 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
419 else if (fepvals->init_lambda >= 0)
421 fepvals->separate_dvdl[efptFEP] = TRUE;
425 else if (file_version >= 64)
427 gmx_fio_do_int(fio, fepvals->n_lambda);
432 snew(fepvals->all_lambda, efptNR);
433 /* still allocate the all_lambda array's contents. */
434 for (g = 0; g < efptNR; g++)
436 if (fepvals->n_lambda > 0)
438 snew(fepvals->all_lambda[g], fepvals->n_lambda);
442 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
444 if (fepvals->init_lambda >= 0)
448 fepvals->separate_dvdl[efptFEP] = TRUE;
452 /* copy the contents of the efptFEP lambda component to all
453 the other components */
454 for (g = 0; g < efptNR; g++)
456 for (h = 0; h < fepvals->n_lambda; h++)
460 fepvals->all_lambda[g][h] =
461 fepvals->all_lambda[efptFEP][h];
470 fepvals->n_lambda = 0;
471 fepvals->all_lambda = NULL;
472 if (fepvals->init_lambda >= 0)
474 fepvals->separate_dvdl[efptFEP] = TRUE;
477 if (file_version >= 13)
479 gmx_fio_do_real(fio, fepvals->sc_alpha);
483 fepvals->sc_alpha = 0;
485 if (file_version >= 38)
487 gmx_fio_do_int(fio, fepvals->sc_power);
491 fepvals->sc_power = 2;
493 if (file_version >= 79)
495 gmx_fio_do_real(fio, fepvals->sc_r_power);
499 fepvals->sc_r_power = 6.0;
501 if (file_version >= 15)
503 gmx_fio_do_real(fio, fepvals->sc_sigma);
507 fepvals->sc_sigma = 0.3;
511 if (file_version >= 71)
513 fepvals->sc_sigma_min = fepvals->sc_sigma;
517 fepvals->sc_sigma_min = 0;
520 if (file_version >= 79)
522 gmx_fio_do_int(fio, fepvals->bScCoul);
526 fepvals->bScCoul = TRUE;
528 if (file_version >= 64)
530 gmx_fio_do_int(fio, fepvals->nstdhdl);
534 fepvals->nstdhdl = 1;
537 if (file_version >= 73)
539 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
540 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
544 fepvals->separate_dhdl_file = esepdhdlfileYES;
545 fepvals->dhdl_derivatives = edhdlderivativesYES;
547 if (file_version >= 71)
549 gmx_fio_do_int(fio, fepvals->dh_hist_size);
550 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
554 fepvals->dh_hist_size = 0;
555 fepvals->dh_hist_spacing = 0.1;
557 if (file_version >= 79)
559 gmx_fio_do_int(fio, fepvals->bPrintEnergy);
563 fepvals->bPrintEnergy = FALSE;
566 /* handle lambda_neighbors */
567 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
569 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
570 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
571 (fepvals->init_lambda < 0) )
573 fepvals->lambda_start_n = (fepvals->init_fep_state -
574 fepvals->lambda_neighbors);
575 fepvals->lambda_stop_n = (fepvals->init_fep_state +
576 fepvals->lambda_neighbors + 1);
577 if (fepvals->lambda_start_n < 0)
579 fepvals->lambda_start_n = 0;;
581 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
583 fepvals->lambda_stop_n = fepvals->n_lambda;
588 fepvals->lambda_start_n = 0;
589 fepvals->lambda_stop_n = fepvals->n_lambda;
594 fepvals->lambda_start_n = 0;
595 fepvals->lambda_stop_n = fepvals->n_lambda;
599 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
603 if (file_version >= 95)
605 gmx_fio_do_int(fio, pull->ngroup);
607 gmx_fio_do_int(fio, pull->ncoord);
608 if (file_version < 95)
610 pull->ngroup = pull->ncoord + 1;
612 gmx_fio_do_int(fio, pull->eGeom);
613 gmx_fio_do_ivec(fio, pull->dim);
614 gmx_fio_do_real(fio, pull->cyl_r1);
615 gmx_fio_do_real(fio, pull->cyl_r0);
616 gmx_fio_do_real(fio, pull->constr_tol);
617 if (file_version >= 95)
619 gmx_fio_do_int(fio, pull->bPrintRef);
621 gmx_fio_do_int(fio, pull->nstxout);
622 gmx_fio_do_int(fio, pull->nstfout);
625 snew(pull->group, pull->ngroup);
626 snew(pull->coord, pull->ncoord);
628 if (file_version < 95)
630 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
631 if (pull->eGeom == epullgDIRPBC)
633 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
635 if (pull->eGeom > epullgDIRPBC)
640 for (g = 0; g < pull->ngroup; g++)
642 /* We read and ignore a pull coordinate for group 0 */
643 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[max(g-1,0)],
644 bRead, file_version);
647 pull->coord[g-1].group[0] = 0;
648 pull->coord[g-1].group[1] = g;
652 pull->bPrintRef = (pull->group[0].nat > 0);
656 for (g = 0; g < pull->ngroup; g++)
658 do_pull_group(fio, &pull->group[g], bRead);
660 for (g = 0; g < pull->ncoord; g++)
662 do_pull_coord(fio, &pull->coord[g]);
668 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
672 gmx_fio_do_int(fio, rotg->eType);
673 gmx_fio_do_int(fio, rotg->bMassW);
674 gmx_fio_do_int(fio, rotg->nat);
677 snew(rotg->ind, rotg->nat);
679 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
682 snew(rotg->x_ref, rotg->nat);
684 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
685 gmx_fio_do_rvec(fio, rotg->vec);
686 gmx_fio_do_rvec(fio, rotg->pivot);
687 gmx_fio_do_real(fio, rotg->rate);
688 gmx_fio_do_real(fio, rotg->k);
689 gmx_fio_do_real(fio, rotg->slab_dist);
690 gmx_fio_do_real(fio, rotg->min_gaussian);
691 gmx_fio_do_real(fio, rotg->eps);
692 gmx_fio_do_int(fio, rotg->eFittype);
693 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
694 gmx_fio_do_real(fio, rotg->PotAngle_step);
697 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
701 gmx_fio_do_int(fio, rot->ngrp);
702 gmx_fio_do_int(fio, rot->nstrout);
703 gmx_fio_do_int(fio, rot->nstsout);
706 snew(rot->grp, rot->ngrp);
708 for (g = 0; g < rot->ngrp; g++)
710 do_rotgrp(fio, &rot->grp[g], bRead);
715 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
716 int file_version, real *fudgeQQ)
718 int i, j, k, *tmp, idum = 0;
722 real zerotemptime, finish_t, init_temp, finish_temp;
724 if (file_version != tpx_version)
726 /* Give a warning about features that are not accessible */
727 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
728 file_version, tpx_version);
736 if (file_version == 0)
741 /* Basic inputrec stuff */
742 gmx_fio_do_int(fio, ir->eI);
743 if (file_version >= 62)
745 gmx_fio_do_int64(fio, ir->nsteps);
749 gmx_fio_do_int(fio, idum);
752 if (file_version > 25)
754 if (file_version >= 62)
756 gmx_fio_do_int64(fio, ir->init_step);
760 gmx_fio_do_int(fio, idum);
761 ir->init_step = idum;
769 if (file_version >= 58)
771 gmx_fio_do_int(fio, ir->simulation_part);
775 ir->simulation_part = 1;
778 if (file_version >= 67)
780 gmx_fio_do_int(fio, ir->nstcalcenergy);
784 ir->nstcalcenergy = 1;
786 if (file_version < 53)
788 /* The pbc info has been moved out of do_inputrec,
789 * since we always want it, also without reading the inputrec.
791 gmx_fio_do_int(fio, ir->ePBC);
792 if ((file_version <= 15) && (ir->ePBC == 2))
796 if (file_version >= 45)
798 gmx_fio_do_int(fio, ir->bPeriodicMols);
805 ir->bPeriodicMols = TRUE;
809 ir->bPeriodicMols = FALSE;
813 if (file_version >= 81)
815 gmx_fio_do_int(fio, ir->cutoff_scheme);
816 if (file_version < 94)
818 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
823 ir->cutoff_scheme = ecutsGROUP;
825 gmx_fio_do_int(fio, ir->ns_type);
826 gmx_fio_do_int(fio, ir->nstlist);
827 gmx_fio_do_int(fio, ir->ndelta);
828 if (file_version < 41)
830 gmx_fio_do_int(fio, idum);
831 gmx_fio_do_int(fio, idum);
833 if (file_version >= 45)
835 gmx_fio_do_real(fio, ir->rtpi);
841 gmx_fio_do_int(fio, ir->nstcomm);
842 if (file_version > 34)
844 gmx_fio_do_int(fio, ir->comm_mode);
846 else if (ir->nstcomm < 0)
848 ir->comm_mode = ecmANGULAR;
852 ir->comm_mode = ecmLINEAR;
854 ir->nstcomm = abs(ir->nstcomm);
856 if (file_version > 25)
858 gmx_fio_do_int(fio, ir->nstcheckpoint);
862 ir->nstcheckpoint = 0;
865 gmx_fio_do_int(fio, ir->nstcgsteep);
867 if (file_version >= 30)
869 gmx_fio_do_int(fio, ir->nbfgscorr);
876 gmx_fio_do_int(fio, ir->nstlog);
877 gmx_fio_do_int(fio, ir->nstxout);
878 gmx_fio_do_int(fio, ir->nstvout);
879 gmx_fio_do_int(fio, ir->nstfout);
880 gmx_fio_do_int(fio, ir->nstenergy);
881 gmx_fio_do_int(fio, ir->nstxtcout);
882 if (file_version >= 59)
884 gmx_fio_do_double(fio, ir->init_t);
885 gmx_fio_do_double(fio, ir->delta_t);
889 gmx_fio_do_real(fio, rdum);
891 gmx_fio_do_real(fio, rdum);
894 gmx_fio_do_real(fio, ir->xtcprec);
895 if (file_version < 19)
897 gmx_fio_do_int(fio, idum);
898 gmx_fio_do_int(fio, idum);
900 if (file_version < 18)
902 gmx_fio_do_int(fio, idum);
904 if (file_version >= 81)
906 gmx_fio_do_real(fio, ir->verletbuf_tol);
910 ir->verletbuf_tol = 0;
912 gmx_fio_do_real(fio, ir->rlist);
913 if (file_version >= 67)
915 gmx_fio_do_real(fio, ir->rlistlong);
917 if (file_version >= 82 && file_version != 90)
919 gmx_fio_do_int(fio, ir->nstcalclr);
923 /* Calculate at NS steps */
924 ir->nstcalclr = ir->nstlist;
926 gmx_fio_do_int(fio, ir->coulombtype);
927 if (file_version < 32 && ir->coulombtype == eelRF)
929 ir->coulombtype = eelRF_NEC;
931 if (file_version >= 81)
933 gmx_fio_do_int(fio, ir->coulomb_modifier);
937 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
939 gmx_fio_do_real(fio, ir->rcoulomb_switch);
940 gmx_fio_do_real(fio, ir->rcoulomb);
941 gmx_fio_do_int(fio, ir->vdwtype);
942 if (file_version >= 81)
944 gmx_fio_do_int(fio, ir->vdw_modifier);
948 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
950 gmx_fio_do_real(fio, ir->rvdw_switch);
951 gmx_fio_do_real(fio, ir->rvdw);
952 if (file_version < 67)
954 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
956 gmx_fio_do_int(fio, ir->eDispCorr);
957 gmx_fio_do_real(fio, ir->epsilon_r);
958 if (file_version >= 37)
960 gmx_fio_do_real(fio, ir->epsilon_rf);
964 if (EEL_RF(ir->coulombtype))
966 ir->epsilon_rf = ir->epsilon_r;
971 ir->epsilon_rf = 1.0;
974 if (file_version >= 29)
976 gmx_fio_do_real(fio, ir->tabext);
983 if (file_version > 25)
985 gmx_fio_do_int(fio, ir->gb_algorithm);
986 gmx_fio_do_int(fio, ir->nstgbradii);
987 gmx_fio_do_real(fio, ir->rgbradii);
988 gmx_fio_do_real(fio, ir->gb_saltconc);
989 gmx_fio_do_int(fio, ir->implicit_solvent);
993 ir->gb_algorithm = egbSTILL;
997 ir->implicit_solvent = eisNO;
999 if (file_version >= 55)
1001 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1002 gmx_fio_do_real(fio, ir->gb_obc_alpha);
1003 gmx_fio_do_real(fio, ir->gb_obc_beta);
1004 gmx_fio_do_real(fio, ir->gb_obc_gamma);
1005 if (file_version >= 60)
1007 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1008 gmx_fio_do_int(fio, ir->sa_algorithm);
1012 ir->gb_dielectric_offset = 0.009;
1013 ir->sa_algorithm = esaAPPROX;
1015 gmx_fio_do_real(fio, ir->sa_surface_tension);
1017 /* Override sa_surface_tension if it is not changed in the mpd-file */
1018 if (ir->sa_surface_tension < 0)
1020 if (ir->gb_algorithm == egbSTILL)
1022 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1024 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1026 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1033 /* Better use sensible values than insane (0.0) ones... */
1034 ir->gb_epsilon_solvent = 80;
1035 ir->gb_obc_alpha = 1.0;
1036 ir->gb_obc_beta = 0.8;
1037 ir->gb_obc_gamma = 4.85;
1038 ir->sa_surface_tension = 2.092;
1042 if (file_version >= 81)
1044 gmx_fio_do_real(fio, ir->fourier_spacing);
1048 ir->fourier_spacing = 0.0;
1050 gmx_fio_do_int(fio, ir->nkx);
1051 gmx_fio_do_int(fio, ir->nky);
1052 gmx_fio_do_int(fio, ir->nkz);
1053 gmx_fio_do_int(fio, ir->pme_order);
1054 gmx_fio_do_real(fio, ir->ewald_rtol);
1056 if (file_version >= 93)
1058 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1062 ir->ewald_rtol_lj = ir->ewald_rtol;
1065 if (file_version >= 24)
1067 gmx_fio_do_int(fio, ir->ewald_geometry);
1071 ir->ewald_geometry = eewg3D;
1074 if (file_version <= 17)
1076 ir->epsilon_surface = 0;
1077 if (file_version == 17)
1079 gmx_fio_do_int(fio, idum);
1084 gmx_fio_do_real(fio, ir->epsilon_surface);
1087 gmx_fio_do_gmx_bool(fio, ir->bOptFFT);
1089 if (file_version >= 93)
1091 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1093 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1094 gmx_fio_do_int(fio, ir->etc);
1095 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1096 * but the values 0 and 1 still mean no and
1097 * berendsen temperature coupling, respectively.
1099 if (file_version >= 79)
1101 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1103 if (file_version >= 71)
1105 gmx_fio_do_int(fio, ir->nsttcouple);
1109 ir->nsttcouple = ir->nstcalcenergy;
1111 if (file_version <= 15)
1113 gmx_fio_do_int(fio, idum);
1115 if (file_version <= 17)
1117 gmx_fio_do_int(fio, ir->epct);
1118 if (file_version <= 15)
1122 ir->epct = epctSURFACETENSION;
1124 gmx_fio_do_int(fio, idum);
1127 /* we have removed the NO alternative at the beginning */
1131 ir->epct = epctISOTROPIC;
1135 ir->epc = epcBERENDSEN;
1140 gmx_fio_do_int(fio, ir->epc);
1141 gmx_fio_do_int(fio, ir->epct);
1143 if (file_version >= 71)
1145 gmx_fio_do_int(fio, ir->nstpcouple);
1149 ir->nstpcouple = ir->nstcalcenergy;
1151 gmx_fio_do_real(fio, ir->tau_p);
1152 if (file_version <= 15)
1154 gmx_fio_do_rvec(fio, vdum);
1155 clear_mat(ir->ref_p);
1156 for (i = 0; i < DIM; i++)
1158 ir->ref_p[i][i] = vdum[i];
1163 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1164 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1165 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1167 if (file_version <= 15)
1169 gmx_fio_do_rvec(fio, vdum);
1170 clear_mat(ir->compress);
1171 for (i = 0; i < DIM; i++)
1173 ir->compress[i][i] = vdum[i];
1178 gmx_fio_do_rvec(fio, ir->compress[XX]);
1179 gmx_fio_do_rvec(fio, ir->compress[YY]);
1180 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1182 if (file_version >= 47)
1184 gmx_fio_do_int(fio, ir->refcoord_scaling);
1185 gmx_fio_do_rvec(fio, ir->posres_com);
1186 gmx_fio_do_rvec(fio, ir->posres_comB);
1190 ir->refcoord_scaling = erscNO;
1191 clear_rvec(ir->posres_com);
1192 clear_rvec(ir->posres_comB);
1194 if ((file_version > 25) && (file_version < 79))
1196 gmx_fio_do_int(fio, ir->andersen_seed);
1200 ir->andersen_seed = 0;
1202 if (file_version < 26)
1204 gmx_fio_do_gmx_bool(fio, bSimAnn);
1205 gmx_fio_do_real(fio, zerotemptime);
1208 if (file_version < 37)
1210 gmx_fio_do_real(fio, rdum);
1213 gmx_fio_do_real(fio, ir->shake_tol);
1214 if (file_version < 54)
1216 gmx_fio_do_real(fio, *fudgeQQ);
1219 gmx_fio_do_int(fio, ir->efep);
1220 if (file_version <= 14 && ir->efep != efepNO)
1224 do_fepvals(fio, ir->fepvals, bRead, file_version);
1226 if (file_version >= 79)
1228 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1231 ir->bSimTemp = TRUE;
1236 ir->bSimTemp = FALSE;
1240 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1243 if (file_version >= 79)
1245 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1248 ir->bExpanded = TRUE;
1252 ir->bExpanded = FALSE;
1257 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1259 if (file_version >= 57)
1261 gmx_fio_do_int(fio, ir->eDisre);
1263 gmx_fio_do_int(fio, ir->eDisreWeighting);
1264 if (file_version < 22)
1266 if (ir->eDisreWeighting == 0)
1268 ir->eDisreWeighting = edrwEqual;
1272 ir->eDisreWeighting = edrwConservative;
1275 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1276 gmx_fio_do_real(fio, ir->dr_fc);
1277 gmx_fio_do_real(fio, ir->dr_tau);
1278 gmx_fio_do_int(fio, ir->nstdisreout);
1279 if (file_version >= 22)
1281 gmx_fio_do_real(fio, ir->orires_fc);
1282 gmx_fio_do_real(fio, ir->orires_tau);
1283 gmx_fio_do_int(fio, ir->nstorireout);
1289 ir->nstorireout = 0;
1291 if (file_version >= 26 && file_version < 79)
1293 gmx_fio_do_real(fio, ir->dihre_fc);
1294 if (file_version < 56)
1296 gmx_fio_do_real(fio, rdum);
1297 gmx_fio_do_int(fio, idum);
1305 gmx_fio_do_real(fio, ir->em_stepsize);
1306 gmx_fio_do_real(fio, ir->em_tol);
1307 if (file_version >= 22)
1309 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1313 ir->bShakeSOR = TRUE;
1315 if (file_version >= 11)
1317 gmx_fio_do_int(fio, ir->niter);
1322 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1325 if (file_version >= 21)
1327 gmx_fio_do_real(fio, ir->fc_stepsize);
1331 ir->fc_stepsize = 0;
1333 gmx_fio_do_int(fio, ir->eConstrAlg);
1334 gmx_fio_do_int(fio, ir->nProjOrder);
1335 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1336 if (file_version <= 14)
1338 gmx_fio_do_int(fio, idum);
1340 if (file_version >= 26)
1342 gmx_fio_do_int(fio, ir->nLincsIter);
1347 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1350 if (file_version < 33)
1352 gmx_fio_do_real(fio, bd_temp);
1354 gmx_fio_do_real(fio, ir->bd_fric);
1355 gmx_fio_do_int(fio, ir->ld_seed);
1356 if (file_version >= 33)
1358 for (i = 0; i < DIM; i++)
1360 gmx_fio_do_rvec(fio, ir->deform[i]);
1365 for (i = 0; i < DIM; i++)
1367 clear_rvec(ir->deform[i]);
1370 if (file_version >= 14)
1372 gmx_fio_do_real(fio, ir->cos_accel);
1378 gmx_fio_do_int(fio, ir->userint1);
1379 gmx_fio_do_int(fio, ir->userint2);
1380 gmx_fio_do_int(fio, ir->userint3);
1381 gmx_fio_do_int(fio, ir->userint4);
1382 gmx_fio_do_real(fio, ir->userreal1);
1383 gmx_fio_do_real(fio, ir->userreal2);
1384 gmx_fio_do_real(fio, ir->userreal3);
1385 gmx_fio_do_real(fio, ir->userreal4);
1388 if (file_version >= 77)
1390 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1395 snew(ir->adress, 1);
1397 gmx_fio_do_int(fio, ir->adress->type);
1398 gmx_fio_do_real(fio, ir->adress->const_wf);
1399 gmx_fio_do_real(fio, ir->adress->ex_width);
1400 gmx_fio_do_real(fio, ir->adress->hy_width);
1401 gmx_fio_do_int(fio, ir->adress->icor);
1402 gmx_fio_do_int(fio, ir->adress->site);
1403 gmx_fio_do_rvec(fio, ir->adress->refs);
1404 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1405 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1406 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1407 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1411 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1413 if (ir->adress->n_tf_grps > 0)
1415 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1419 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1421 if (ir->adress->n_energy_grps > 0)
1423 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1429 ir->bAdress = FALSE;
1433 if (file_version >= 48)
1435 gmx_fio_do_int(fio, ir->ePull);
1436 if (ir->ePull != epullNO)
1442 do_pull(fio, ir->pull, bRead, file_version);
1447 ir->ePull = epullNO;
1450 /* Enforced rotation */
1451 if (file_version >= 74)
1453 gmx_fio_do_int(fio, ir->bRot);
1454 if (ir->bRot == TRUE)
1460 do_rot(fio, ir->rot, bRead);
1469 gmx_fio_do_int(fio, ir->opts.ngtc);
1470 if (file_version >= 69)
1472 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1476 ir->opts.nhchainlength = 1;
1478 gmx_fio_do_int(fio, ir->opts.ngacc);
1479 gmx_fio_do_int(fio, ir->opts.ngfrz);
1480 gmx_fio_do_int(fio, ir->opts.ngener);
1484 snew(ir->opts.nrdf, ir->opts.ngtc);
1485 snew(ir->opts.ref_t, ir->opts.ngtc);
1486 snew(ir->opts.annealing, ir->opts.ngtc);
1487 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1488 snew(ir->opts.anneal_time, ir->opts.ngtc);
1489 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1490 snew(ir->opts.tau_t, ir->opts.ngtc);
1491 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1492 snew(ir->opts.acc, ir->opts.ngacc);
1493 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1495 if (ir->opts.ngtc > 0)
1497 if (bRead && file_version < 13)
1499 snew(tmp, ir->opts.ngtc);
1500 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1501 for (i = 0; i < ir->opts.ngtc; i++)
1503 ir->opts.nrdf[i] = tmp[i];
1509 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1511 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1512 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1513 if (file_version < 33 && ir->eI == eiBD)
1515 for (i = 0; i < ir->opts.ngtc; i++)
1517 ir->opts.tau_t[i] = bd_temp;
1521 if (ir->opts.ngfrz > 0)
1523 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1525 if (ir->opts.ngacc > 0)
1527 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1529 if (file_version >= 12)
1531 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1532 ir->opts.ngener*ir->opts.ngener);
1535 if (bRead && file_version < 26)
1537 for (i = 0; i < ir->opts.ngtc; i++)
1541 ir->opts.annealing[i] = eannSINGLE;
1542 ir->opts.anneal_npoints[i] = 2;
1543 snew(ir->opts.anneal_time[i], 2);
1544 snew(ir->opts.anneal_temp[i], 2);
1545 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1546 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1547 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1548 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1549 ir->opts.anneal_time[i][0] = ir->init_t;
1550 ir->opts.anneal_time[i][1] = finish_t;
1551 ir->opts.anneal_temp[i][0] = init_temp;
1552 ir->opts.anneal_temp[i][1] = finish_temp;
1556 ir->opts.annealing[i] = eannNO;
1557 ir->opts.anneal_npoints[i] = 0;
1563 /* file version 26 or later */
1564 /* First read the lists with annealing and npoints for each group */
1565 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1566 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1567 for (j = 0; j < (ir->opts.ngtc); j++)
1569 k = ir->opts.anneal_npoints[j];
1572 snew(ir->opts.anneal_time[j], k);
1573 snew(ir->opts.anneal_temp[j], k);
1575 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1576 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1580 if (file_version >= 45)
1582 gmx_fio_do_int(fio, ir->nwall);
1583 gmx_fio_do_int(fio, ir->wall_type);
1584 if (file_version >= 50)
1586 gmx_fio_do_real(fio, ir->wall_r_linpot);
1590 ir->wall_r_linpot = -1;
1592 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1593 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1594 gmx_fio_do_real(fio, ir->wall_density[0]);
1595 gmx_fio_do_real(fio, ir->wall_density[1]);
1596 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1602 ir->wall_atomtype[0] = -1;
1603 ir->wall_atomtype[1] = -1;
1604 ir->wall_density[0] = 0;
1605 ir->wall_density[1] = 0;
1606 ir->wall_ewald_zfac = 3;
1608 /* Cosine stuff for electric fields */
1609 for (j = 0; (j < DIM); j++)
1611 gmx_fio_do_int(fio, ir->ex[j].n);
1612 gmx_fio_do_int(fio, ir->et[j].n);
1615 snew(ir->ex[j].a, ir->ex[j].n);
1616 snew(ir->ex[j].phi, ir->ex[j].n);
1617 snew(ir->et[j].a, ir->et[j].n);
1618 snew(ir->et[j].phi, ir->et[j].n);
1620 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1621 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1622 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1623 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1627 if (file_version >= 39)
1629 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1630 gmx_fio_do_int(fio, ir->QMMMscheme);
1631 gmx_fio_do_real(fio, ir->scalefactor);
1632 gmx_fio_do_int(fio, ir->opts.ngQM);
1635 snew(ir->opts.QMmethod, ir->opts.ngQM);
1636 snew(ir->opts.QMbasis, ir->opts.ngQM);
1637 snew(ir->opts.QMcharge, ir->opts.ngQM);
1638 snew(ir->opts.QMmult, ir->opts.ngQM);
1639 snew(ir->opts.bSH, ir->opts.ngQM);
1640 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1641 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1642 snew(ir->opts.SAon, ir->opts.ngQM);
1643 snew(ir->opts.SAoff, ir->opts.ngQM);
1644 snew(ir->opts.SAsteps, ir->opts.ngQM);
1645 snew(ir->opts.bOPT, ir->opts.ngQM);
1646 snew(ir->opts.bTS, ir->opts.ngQM);
1648 if (ir->opts.ngQM > 0)
1650 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1651 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1652 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1653 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1654 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1655 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1656 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1657 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1658 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1659 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1660 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1661 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1663 /* end of QMMM stuff */
1668 static void do_harm(t_fileio *fio, t_iparams *iparams)
1670 gmx_fio_do_real(fio, iparams->harmonic.rA);
1671 gmx_fio_do_real(fio, iparams->harmonic.krA);
1672 gmx_fio_do_real(fio, iparams->harmonic.rB);
1673 gmx_fio_do_real(fio, iparams->harmonic.krB);
1676 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1677 gmx_bool bRead, int file_version)
1684 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1694 do_harm(fio, iparams);
1695 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1697 /* Correct incorrect storage of parameters */
1698 iparams->pdihs.phiB = iparams->pdihs.phiA;
1699 iparams->pdihs.cpB = iparams->pdihs.cpA;
1702 case F_LINEAR_ANGLES:
1703 gmx_fio_do_real(fio, iparams->linangle.klinA);
1704 gmx_fio_do_real(fio, iparams->linangle.aA);
1705 gmx_fio_do_real(fio, iparams->linangle.klinB);
1706 gmx_fio_do_real(fio, iparams->linangle.aB);
1709 gmx_fio_do_real(fio, iparams->fene.bm);
1710 gmx_fio_do_real(fio, iparams->fene.kb);
1713 gmx_fio_do_real(fio, iparams->restraint.lowA);
1714 gmx_fio_do_real(fio, iparams->restraint.up1A);
1715 gmx_fio_do_real(fio, iparams->restraint.up2A);
1716 gmx_fio_do_real(fio, iparams->restraint.kA);
1717 gmx_fio_do_real(fio, iparams->restraint.lowB);
1718 gmx_fio_do_real(fio, iparams->restraint.up1B);
1719 gmx_fio_do_real(fio, iparams->restraint.up2B);
1720 gmx_fio_do_real(fio, iparams->restraint.kB);
1726 gmx_fio_do_real(fio, iparams->tab.kA);
1727 gmx_fio_do_int(fio, iparams->tab.table);
1728 gmx_fio_do_real(fio, iparams->tab.kB);
1730 case F_CROSS_BOND_BONDS:
1731 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1732 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1733 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1735 case F_CROSS_BOND_ANGLES:
1736 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1737 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1738 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1739 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1741 case F_UREY_BRADLEY:
1742 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1743 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1744 gmx_fio_do_real(fio, iparams->u_b.r13A);
1745 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1746 if (file_version >= 79)
1748 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1749 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1750 gmx_fio_do_real(fio, iparams->u_b.r13B);
1751 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1755 iparams->u_b.thetaB = iparams->u_b.thetaA;
1756 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1757 iparams->u_b.r13B = iparams->u_b.r13A;
1758 iparams->u_b.kUBB = iparams->u_b.kUBA;
1761 case F_QUARTIC_ANGLES:
1762 gmx_fio_do_real(fio, iparams->qangle.theta);
1763 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1766 gmx_fio_do_real(fio, iparams->bham.a);
1767 gmx_fio_do_real(fio, iparams->bham.b);
1768 gmx_fio_do_real(fio, iparams->bham.c);
1771 gmx_fio_do_real(fio, iparams->morse.b0A);
1772 gmx_fio_do_real(fio, iparams->morse.cbA);
1773 gmx_fio_do_real(fio, iparams->morse.betaA);
1774 if (file_version >= 79)
1776 gmx_fio_do_real(fio, iparams->morse.b0B);
1777 gmx_fio_do_real(fio, iparams->morse.cbB);
1778 gmx_fio_do_real(fio, iparams->morse.betaB);
1782 iparams->morse.b0B = iparams->morse.b0A;
1783 iparams->morse.cbB = iparams->morse.cbA;
1784 iparams->morse.betaB = iparams->morse.betaA;
1788 gmx_fio_do_real(fio, iparams->cubic.b0);
1789 gmx_fio_do_real(fio, iparams->cubic.kb);
1790 gmx_fio_do_real(fio, iparams->cubic.kcub);
1794 case F_POLARIZATION:
1795 gmx_fio_do_real(fio, iparams->polarize.alpha);
1798 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1799 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1800 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1803 if (file_version < 31)
1805 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1807 gmx_fio_do_real(fio, iparams->wpol.al_x);
1808 gmx_fio_do_real(fio, iparams->wpol.al_y);
1809 gmx_fio_do_real(fio, iparams->wpol.al_z);
1810 gmx_fio_do_real(fio, iparams->wpol.rOH);
1811 gmx_fio_do_real(fio, iparams->wpol.rHH);
1812 gmx_fio_do_real(fio, iparams->wpol.rOD);
1815 gmx_fio_do_real(fio, iparams->thole.a);
1816 gmx_fio_do_real(fio, iparams->thole.alpha1);
1817 gmx_fio_do_real(fio, iparams->thole.alpha2);
1818 gmx_fio_do_real(fio, iparams->thole.rfac);
1821 gmx_fio_do_real(fio, iparams->lj.c6);
1822 gmx_fio_do_real(fio, iparams->lj.c12);
1825 gmx_fio_do_real(fio, iparams->lj14.c6A);
1826 gmx_fio_do_real(fio, iparams->lj14.c12A);
1827 gmx_fio_do_real(fio, iparams->lj14.c6B);
1828 gmx_fio_do_real(fio, iparams->lj14.c12B);
1831 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1832 gmx_fio_do_real(fio, iparams->ljc14.qi);
1833 gmx_fio_do_real(fio, iparams->ljc14.qj);
1834 gmx_fio_do_real(fio, iparams->ljc14.c6);
1835 gmx_fio_do_real(fio, iparams->ljc14.c12);
1837 case F_LJC_PAIRS_NB:
1838 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1839 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1840 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1841 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1847 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1848 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1849 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1851 /* Read the incorrectly stored multiplicity */
1852 gmx_fio_do_real(fio, iparams->harmonic.rB);
1853 gmx_fio_do_real(fio, iparams->harmonic.krB);
1854 iparams->pdihs.phiB = iparams->pdihs.phiA;
1855 iparams->pdihs.cpB = iparams->pdihs.cpA;
1859 gmx_fio_do_real(fio, iparams->pdihs.phiB);
1860 gmx_fio_do_real(fio, iparams->pdihs.cpB);
1861 gmx_fio_do_int(fio, iparams->pdihs.mult);
1865 gmx_fio_do_int(fio, iparams->disres.label);
1866 gmx_fio_do_int(fio, iparams->disres.type);
1867 gmx_fio_do_real(fio, iparams->disres.low);
1868 gmx_fio_do_real(fio, iparams->disres.up1);
1869 gmx_fio_do_real(fio, iparams->disres.up2);
1870 gmx_fio_do_real(fio, iparams->disres.kfac);
1873 gmx_fio_do_int(fio, iparams->orires.ex);
1874 gmx_fio_do_int(fio, iparams->orires.label);
1875 gmx_fio_do_int(fio, iparams->orires.power);
1876 gmx_fio_do_real(fio, iparams->orires.c);
1877 gmx_fio_do_real(fio, iparams->orires.obs);
1878 gmx_fio_do_real(fio, iparams->orires.kfac);
1881 if (file_version < 82)
1883 gmx_fio_do_int(fio, idum);
1884 gmx_fio_do_int(fio, idum);
1886 gmx_fio_do_real(fio, iparams->dihres.phiA);
1887 gmx_fio_do_real(fio, iparams->dihres.dphiA);
1888 gmx_fio_do_real(fio, iparams->dihres.kfacA);
1889 if (file_version >= 82)
1891 gmx_fio_do_real(fio, iparams->dihres.phiB);
1892 gmx_fio_do_real(fio, iparams->dihres.dphiB);
1893 gmx_fio_do_real(fio, iparams->dihres.kfacB);
1897 iparams->dihres.phiB = iparams->dihres.phiA;
1898 iparams->dihres.dphiB = iparams->dihres.dphiA;
1899 iparams->dihres.kfacB = iparams->dihres.kfacA;
1903 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
1904 gmx_fio_do_rvec(fio, iparams->posres.fcA);
1905 if (bRead && file_version < 27)
1907 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
1908 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
1912 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
1913 gmx_fio_do_rvec(fio, iparams->posres.fcB);
1917 gmx_fio_do_int(fio, iparams->fbposres.geom);
1918 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
1919 gmx_fio_do_real(fio, iparams->fbposres.r);
1920 gmx_fio_do_real(fio, iparams->fbposres.k);
1923 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1924 if (file_version >= 25)
1926 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1930 /* Fourier dihedrals are internally represented
1931 * as Ryckaert-Bellemans since those are faster to compute.
1933 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1934 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1938 gmx_fio_do_real(fio, iparams->constr.dA);
1939 gmx_fio_do_real(fio, iparams->constr.dB);
1942 gmx_fio_do_real(fio, iparams->settle.doh);
1943 gmx_fio_do_real(fio, iparams->settle.dhh);
1946 gmx_fio_do_real(fio, iparams->vsite.a);
1951 gmx_fio_do_real(fio, iparams->vsite.a);
1952 gmx_fio_do_real(fio, iparams->vsite.b);
1957 gmx_fio_do_real(fio, iparams->vsite.a);
1958 gmx_fio_do_real(fio, iparams->vsite.b);
1959 gmx_fio_do_real(fio, iparams->vsite.c);
1962 gmx_fio_do_int(fio, iparams->vsiten.n);
1963 gmx_fio_do_real(fio, iparams->vsiten.a);
1968 /* We got rid of some parameters in version 68 */
1969 if (bRead && file_version < 68)
1971 gmx_fio_do_real(fio, rdum);
1972 gmx_fio_do_real(fio, rdum);
1973 gmx_fio_do_real(fio, rdum);
1974 gmx_fio_do_real(fio, rdum);
1976 gmx_fio_do_real(fio, iparams->gb.sar);
1977 gmx_fio_do_real(fio, iparams->gb.st);
1978 gmx_fio_do_real(fio, iparams->gb.pi);
1979 gmx_fio_do_real(fio, iparams->gb.gbr);
1980 gmx_fio_do_real(fio, iparams->gb.bmlt);
1983 gmx_fio_do_int(fio, iparams->cmap.cmapA);
1984 gmx_fio_do_int(fio, iparams->cmap.cmapB);
1987 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1988 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1992 gmx_fio_unset_comment(fio);
1996 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
2003 gmx_fio_set_comment(fio, interaction_function[ftype].name);
2005 if (file_version < 44)
2007 for (i = 0; i < MAXNODES; i++)
2009 gmx_fio_do_int(fio, idum);
2012 gmx_fio_do_int(fio, ilist->nr);
2015 snew(ilist->iatoms, ilist->nr);
2017 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2020 gmx_fio_unset_comment(fio);
2024 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2025 gmx_bool bRead, int file_version)
2030 gmx_fio_do_int(fio, ffparams->atnr);
2031 if (file_version < 57)
2033 gmx_fio_do_int(fio, idum);
2035 gmx_fio_do_int(fio, ffparams->ntypes);
2038 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2039 ffparams->atnr, ffparams->ntypes);
2043 snew(ffparams->functype, ffparams->ntypes);
2044 snew(ffparams->iparams, ffparams->ntypes);
2046 /* Read/write all the function types */
2047 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2050 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2053 if (file_version >= 66)
2055 gmx_fio_do_double(fio, ffparams->reppow);
2059 ffparams->reppow = 12.0;
2062 if (file_version >= 57)
2064 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2067 /* Check whether all these function types are supported by the code.
2068 * In practice the code is backwards compatible, which means that the
2069 * numbering may have to be altered from old numbering to new numbering
2071 for (i = 0; (i < ffparams->ntypes); i++)
2075 /* Loop over file versions */
2076 for (k = 0; (k < NFTUPD); k++)
2078 /* Compare the read file_version to the update table */
2079 if ((file_version < ftupd[k].fvnr) &&
2080 (ffparams->functype[i] >= ftupd[k].ftype))
2082 ffparams->functype[i] += 1;
2085 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2086 i, ffparams->functype[i],
2087 interaction_function[ftupd[k].ftype].longname);
2094 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2098 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2103 static void add_settle_atoms(t_ilist *ilist)
2107 /* Settle used to only store the first atom: add the other two */
2108 srenew(ilist->iatoms, 2*ilist->nr);
2109 for (i = ilist->nr/2-1; i >= 0; i--)
2111 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2112 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2113 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2114 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2116 ilist->nr = 2*ilist->nr;
2119 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2122 int i, j, renum[F_NRE];
2126 for (j = 0; (j < F_NRE); j++)
2131 for (k = 0; k < NFTUPD; k++)
2133 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2142 ilist[j].iatoms = NULL;
2146 do_ilist(fio, &ilist[j], bRead, file_version, j);
2147 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2149 add_settle_atoms(&ilist[j]);
2153 if (bRead && gmx_debug_at)
2154 pr_ilist(debug,0,interaction_function[j].longname,
2155 functype,&ilist[j],TRUE);
2160 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2161 gmx_bool bRead, int file_version)
2163 do_ffparams(fio, ffparams, bRead, file_version);
2165 if (file_version >= 54)
2167 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2170 do_ilists(fio, molt->ilist, bRead, file_version);
2173 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2175 int i, idum, dum_nra, *dum_a;
2177 if (file_version < 44)
2179 for (i = 0; i < MAXNODES; i++)
2181 gmx_fio_do_int(fio, idum);
2184 gmx_fio_do_int(fio, block->nr);
2185 if (file_version < 51)
2187 gmx_fio_do_int(fio, dum_nra);
2191 if ((block->nalloc_index > 0) && (NULL != block->index))
2193 sfree(block->index);
2195 block->nalloc_index = block->nr+1;
2196 snew(block->index, block->nalloc_index);
2198 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2200 if (file_version < 51 && dum_nra > 0)
2202 snew(dum_a, dum_nra);
2203 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2208 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2213 if (file_version < 44)
2215 for (i = 0; i < MAXNODES; i++)
2217 gmx_fio_do_int(fio, idum);
2220 gmx_fio_do_int(fio, block->nr);
2221 gmx_fio_do_int(fio, block->nra);
2224 block->nalloc_index = block->nr+1;
2225 snew(block->index, block->nalloc_index);
2226 block->nalloc_a = block->nra;
2227 snew(block->a, block->nalloc_a);
2229 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2230 gmx_fio_ndo_int(fio, block->a, block->nra);
2233 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2234 int file_version, gmx_groups_t *groups, int atnr)
2238 gmx_fio_do_real(fio, atom->m);
2239 gmx_fio_do_real(fio, atom->q);
2240 gmx_fio_do_real(fio, atom->mB);
2241 gmx_fio_do_real(fio, atom->qB);
2242 gmx_fio_do_ushort(fio, atom->type);
2243 gmx_fio_do_ushort(fio, atom->typeB);
2244 gmx_fio_do_int(fio, atom->ptype);
2245 gmx_fio_do_int(fio, atom->resind);
2246 if (file_version >= 52)
2248 gmx_fio_do_int(fio, atom->atomnumber);
2252 atom->atomnumber = NOTSET;
2254 if (file_version < 23)
2258 else if (file_version < 39)
2267 if (file_version < 57)
2269 unsigned char uchar[egcNR];
2270 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2271 for (i = myngrp; (i < ngrp); i++)
2275 /* Copy the old data format to the groups struct */
2276 for (i = 0; i < ngrp; i++)
2278 groups->grpnr[i][atnr] = uchar[i];
2283 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2288 if (file_version < 23)
2292 else if (file_version < 39)
2301 for (j = 0; (j < ngrp); j++)
2305 gmx_fio_do_int(fio, grps[j].nr);
2308 snew(grps[j].nm_ind, grps[j].nr);
2310 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2315 snew(grps[j].nm_ind, grps[j].nr);
2320 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2326 gmx_fio_do_int(fio, ls);
2327 *nm = get_symtab_handle(symtab, ls);
2331 ls = lookup_symtab(symtab, *nm);
2332 gmx_fio_do_int(fio, ls);
2336 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2341 for (j = 0; (j < nstr); j++)
2343 do_symstr(fio, &(nm[j]), bRead, symtab);
2347 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2348 t_symtab *symtab, int file_version)
2352 for (j = 0; (j < n); j++)
2354 do_symstr(fio, &(ri[j].name), bRead, symtab);
2355 if (file_version >= 63)
2357 gmx_fio_do_int(fio, ri[j].nr);
2358 gmx_fio_do_uchar(fio, ri[j].ic);
2368 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2370 gmx_groups_t *groups)
2374 gmx_fio_do_int(fio, atoms->nr);
2375 gmx_fio_do_int(fio, atoms->nres);
2376 if (file_version < 57)
2378 gmx_fio_do_int(fio, groups->ngrpname);
2379 for (i = 0; i < egcNR; i++)
2381 groups->ngrpnr[i] = atoms->nr;
2382 snew(groups->grpnr[i], groups->ngrpnr[i]);
2387 snew(atoms->atom, atoms->nr);
2388 snew(atoms->atomname, atoms->nr);
2389 snew(atoms->atomtype, atoms->nr);
2390 snew(atoms->atomtypeB, atoms->nr);
2391 snew(atoms->resinfo, atoms->nres);
2392 if (file_version < 57)
2394 snew(groups->grpname, groups->ngrpname);
2396 atoms->pdbinfo = NULL;
2398 for (i = 0; (i < atoms->nr); i++)
2400 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2402 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2403 if (bRead && (file_version <= 20))
2405 for (i = 0; i < atoms->nr; i++)
2407 atoms->atomtype[i] = put_symtab(symtab, "?");
2408 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2413 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2414 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2416 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2418 if (file_version < 57)
2420 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2422 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2426 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2427 gmx_bool bRead, t_symtab *symtab,
2432 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2433 gmx_fio_do_int(fio, groups->ngrpname);
2436 snew(groups->grpname, groups->ngrpname);
2438 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2439 for (g = 0; g < egcNR; g++)
2441 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2442 if (groups->ngrpnr[g] == 0)
2446 groups->grpnr[g] = NULL;
2453 snew(groups->grpnr[g], groups->ngrpnr[g]);
2455 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2460 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2465 if (file_version > 25)
2467 gmx_fio_do_int(fio, atomtypes->nr);
2471 snew(atomtypes->radius, j);
2472 snew(atomtypes->vol, j);
2473 snew(atomtypes->surftens, j);
2474 snew(atomtypes->atomnumber, j);
2475 snew(atomtypes->gb_radius, j);
2476 snew(atomtypes->S_hct, j);
2478 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2479 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2480 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2481 if (file_version >= 40)
2483 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2485 if (file_version >= 60)
2487 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2488 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2493 /* File versions prior to 26 cannot do GBSA,
2494 * so they dont use this structure
2497 atomtypes->radius = NULL;
2498 atomtypes->vol = NULL;
2499 atomtypes->surftens = NULL;
2500 atomtypes->atomnumber = NULL;
2501 atomtypes->gb_radius = NULL;
2502 atomtypes->S_hct = NULL;
2506 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2512 gmx_fio_do_int(fio, symtab->nr);
2516 snew(symtab->symbuf, 1);
2517 symbuf = symtab->symbuf;
2518 symbuf->bufsize = nr;
2519 snew(symbuf->buf, nr);
2520 for (i = 0; (i < nr); i++)
2522 gmx_fio_do_string(fio, buf);
2523 symbuf->buf[i] = strdup(buf);
2528 symbuf = symtab->symbuf;
2529 while (symbuf != NULL)
2531 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2533 gmx_fio_do_string(fio, symbuf->buf[i]);
2536 symbuf = symbuf->next;
2540 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2545 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2547 int i, j, ngrid, gs, nelem;
2549 gmx_fio_do_int(fio, cmap_grid->ngrid);
2550 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2552 ngrid = cmap_grid->ngrid;
2553 gs = cmap_grid->grid_spacing;
2558 snew(cmap_grid->cmapdata, ngrid);
2560 for (i = 0; i < cmap_grid->ngrid; i++)
2562 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2566 for (i = 0; i < cmap_grid->ngrid; i++)
2568 for (j = 0; j < nelem; j++)
2570 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2571 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2572 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2573 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2579 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2581 int m, a, a0, a1, r;
2585 /* We always assign a new chain number, but save the chain id characters
2586 * for larger molecules.
2588 #define CHAIN_MIN_ATOMS 15
2592 for (m = 0; m < mols->nr; m++)
2594 a0 = mols->index[m];
2595 a1 = mols->index[m+1];
2596 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2605 for (a = a0; a < a1; a++)
2607 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2608 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2613 /* Blank out the chain id if there was only one chain */
2616 for (r = 0; r < atoms->nres; r++)
2618 atoms->resinfo[r].chainid = ' ';
2623 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2624 t_symtab *symtab, int file_version,
2625 gmx_groups_t *groups)
2629 if (file_version >= 57)
2631 do_symstr(fio, &(molt->name), bRead, symtab);
2634 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2636 if (bRead && gmx_debug_at)
2638 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2641 if (file_version >= 57)
2643 do_ilists(fio, molt->ilist, bRead, file_version);
2645 do_block(fio, &molt->cgs, bRead, file_version);
2646 if (bRead && gmx_debug_at)
2648 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2652 /* This used to be in the atoms struct */
2653 do_blocka(fio, &molt->excls, bRead, file_version);
2656 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2660 gmx_fio_do_int(fio, molb->type);
2661 gmx_fio_do_int(fio, molb->nmol);
2662 gmx_fio_do_int(fio, molb->natoms_mol);
2663 /* Position restraint coordinates */
2664 gmx_fio_do_int(fio, molb->nposres_xA);
2665 if (molb->nposres_xA > 0)
2669 snew(molb->posres_xA, molb->nposres_xA);
2671 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2673 gmx_fio_do_int(fio, molb->nposres_xB);
2674 if (molb->nposres_xB > 0)
2678 snew(molb->posres_xB, molb->nposres_xB);
2680 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2685 static t_block mtop_mols(gmx_mtop_t *mtop)
2691 for (mb = 0; mb < mtop->nmolblock; mb++)
2693 mols.nr += mtop->molblock[mb].nmol;
2695 mols.nalloc_index = mols.nr + 1;
2696 snew(mols.index, mols.nalloc_index);
2701 for (mb = 0; mb < mtop->nmolblock; mb++)
2703 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2705 a += mtop->molblock[mb].natoms_mol;
2714 static void add_posres_molblock(gmx_mtop_t *mtop)
2719 gmx_molblock_t *molb;
2722 /* posres reference positions are stored in ip->posres (if present) and
2723 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2724 posres.pos0A are identical to fbposres.pos0. */
2725 il = &mtop->moltype[0].ilist[F_POSRES];
2726 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2727 if (il->nr == 0 && ilfb->nr == 0)
2733 for (i = 0; i < il->nr; i += 2)
2735 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2736 am = max(am, il->iatoms[i+1]);
2737 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2738 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2739 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2744 /* This loop is required if we have only flat-bottomed posres:
2746 - bFE == FALSE (no B-state for flat-bottomed posres) */
2749 for (i = 0; i < ilfb->nr; i += 2)
2751 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2752 am = max(am, ilfb->iatoms[i+1]);
2755 /* Make the posres coordinate block end at a molecule end */
2757 while (am >= mtop->mols.index[mol+1])
2761 molb = &mtop->molblock[0];
2762 molb->nposres_xA = mtop->mols.index[mol+1];
2763 snew(molb->posres_xA, molb->nposres_xA);
2766 molb->nposres_xB = molb->nposres_xA;
2767 snew(molb->posres_xB, molb->nposres_xB);
2771 molb->nposres_xB = 0;
2773 for (i = 0; i < il->nr; i += 2)
2775 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2776 a = il->iatoms[i+1];
2777 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2778 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2779 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2782 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2783 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2784 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2789 /* If only flat-bottomed posres are present, take reference pos from them.
2790 Here: bFE == FALSE */
2791 for (i = 0; i < ilfb->nr; i += 2)
2793 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2794 a = ilfb->iatoms[i+1];
2795 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2796 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2797 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2802 static void set_disres_npair(gmx_mtop_t *mtop)
2809 ip = mtop->ffparams.iparams;
2811 for (mt = 0; mt < mtop->nmoltype; mt++)
2813 il = &mtop->moltype[mt].ilist[F_DISRES];
2818 for (i = 0; i < il->nr; i += 3)
2821 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2823 ip[a[i]].disres.npair = npair;
2831 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2841 do_symtab(fio, &(mtop->symtab), bRead);
2844 pr_symtab(debug, 0, "symtab", &mtop->symtab);
2847 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2849 if (file_version >= 57)
2851 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
2853 gmx_fio_do_int(fio, mtop->nmoltype);
2861 snew(mtop->moltype, mtop->nmoltype);
2862 if (file_version < 57)
2864 mtop->moltype[0].name = mtop->name;
2867 for (mt = 0; mt < mtop->nmoltype; mt++)
2869 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
2873 if (file_version >= 57)
2875 gmx_fio_do_int(fio, mtop->nmolblock);
2879 mtop->nmolblock = 1;
2883 snew(mtop->molblock, mtop->nmolblock);
2885 if (file_version >= 57)
2887 for (mb = 0; mb < mtop->nmolblock; mb++)
2889 do_molblock(fio, &mtop->molblock[mb], bRead);
2891 gmx_fio_do_int(fio, mtop->natoms);
2895 mtop->molblock[0].type = 0;
2896 mtop->molblock[0].nmol = 1;
2897 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
2898 mtop->molblock[0].nposres_xA = 0;
2899 mtop->molblock[0].nposres_xB = 0;
2902 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
2905 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
2908 if (file_version < 57)
2910 /* Debug statements are inside do_idef */
2911 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
2912 mtop->natoms = mtop->moltype[0].atoms.nr;
2915 if (file_version >= 65)
2917 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
2921 mtop->ffparams.cmap_grid.ngrid = 0;
2922 mtop->ffparams.cmap_grid.grid_spacing = 0;
2923 mtop->ffparams.cmap_grid.cmapdata = NULL;
2926 if (file_version >= 57)
2928 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
2931 if (file_version < 57)
2933 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
2934 if (bRead && gmx_debug_at)
2936 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
2938 do_block(fio, &mtop->mols, bRead, file_version);
2939 /* Add the posres coordinates to the molblock */
2940 add_posres_molblock(mtop);
2944 if (file_version >= 57)
2946 done_block(&mtop->mols);
2947 mtop->mols = mtop_mols(mtop);
2951 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
2955 if (file_version < 51)
2957 /* Here used to be the shake blocks */
2958 do_blocka(fio, &dumb, bRead, file_version);
2971 close_symtab(&(mtop->symtab));
2975 /* If TopOnlyOK is TRUE then we can read even future versions
2976 * of tpx files, provided the file_generation hasn't changed.
2977 * If it is FALSE, we need the inputrecord too, and bail out
2978 * if the file is newer than the program.
2980 * The version and generation if the topology (see top of this file)
2981 * are returned in the two last arguments.
2983 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2985 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
2986 gmx_bool TopOnlyOK, int *file_version,
2987 int *file_generation)
2990 char file_tag[STRLEN];
2997 gmx_fio_checktype(fio);
2998 gmx_fio_setdebug(fio, bDebugMode());
3000 /* NEW! XDR tpb file */
3001 precision = sizeof(real);
3004 gmx_fio_do_string(fio, buf);
3005 if (strncmp(buf, "VERSION", 7))
3007 gmx_fatal(FARGS, "Can not read file %s,\n"
3008 " this file is from a Gromacs version which is older than 2.0\n"
3009 " Make a new one with grompp or use a gro or pdb file, if possible",
3010 gmx_fio_getname(fio));
3012 gmx_fio_do_int(fio, precision);
3013 bDouble = (precision == sizeof(double));
3014 if ((precision != sizeof(float)) && !bDouble)
3016 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3017 "instead of %d or %d",
3018 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3020 gmx_fio_setprecision(fio, bDouble);
3021 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3022 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3026 gmx_fio_write_string(fio, GromacsVersion());
3027 bDouble = (precision == sizeof(double));
3028 gmx_fio_setprecision(fio, bDouble);
3029 gmx_fio_do_int(fio, precision);
3031 sprintf(file_tag, "%s", tpx_tag);
3032 fgen = tpx_generation;
3035 /* Check versions! */
3036 gmx_fio_do_int(fio, fver);
3038 /* This is for backward compatibility with development versions 77-79
3039 * where the tag was, mistakenly, placed before the generation,
3040 * which would cause a segv instead of a proper error message
3041 * when reading the topology only from tpx with <77 code.
3043 if (fver >= 77 && fver <= 79)
3045 gmx_fio_do_string(fio, file_tag);
3050 gmx_fio_do_int(fio, fgen);
3059 gmx_fio_do_string(fio, file_tag);
3065 /* Versions before 77 don't have the tag, set it to release */
3066 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3069 if (strcmp(file_tag, tpx_tag) != 0)
3071 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3074 /* We only support reading tpx files with the same tag as the code
3075 * or tpx files with the release tag and with lower version number.
3077 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
3079 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3080 gmx_fio_getname(fio), fver, file_tag,
3081 tpx_version, tpx_tag);
3086 if (file_version != NULL)
3088 *file_version = fver;
3090 if (file_generation != NULL)
3092 *file_generation = fgen;
3096 if ((fver <= tpx_incompatible_version) ||
3097 ((fver > tpx_version) && !TopOnlyOK) ||
3098 (fgen > tpx_generation) ||
3099 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3101 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3102 gmx_fio_getname(fio), fver, tpx_version);
3105 do_section(fio, eitemHEADER, bRead);
3106 gmx_fio_do_int(fio, tpx->natoms);
3109 gmx_fio_do_int(fio, tpx->ngtc);
3117 gmx_fio_do_int(fio, idum);
3118 gmx_fio_do_real(fio, rdum);
3120 /*a better decision will eventually (5.0 or later) need to be made
3121 on how to treat the alchemical state of the system, which can now
3122 vary through a simulation, and cannot be completely described
3123 though a single lambda variable, or even a single state
3124 index. Eventually, should probably be a vector. MRS*/
3127 gmx_fio_do_int(fio, tpx->fep_state);
3129 gmx_fio_do_real(fio, tpx->lambda);
3130 gmx_fio_do_int(fio, tpx->bIr);
3131 gmx_fio_do_int(fio, tpx->bTop);
3132 gmx_fio_do_int(fio, tpx->bX);
3133 gmx_fio_do_int(fio, tpx->bV);
3134 gmx_fio_do_int(fio, tpx->bF);
3135 gmx_fio_do_int(fio, tpx->bBox);
3137 if ((fgen > tpx_generation))
3139 /* This can only happen if TopOnlyOK=TRUE */
3144 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3145 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3146 gmx_bool bXVallocated)
3152 int file_version, file_generation;
3156 gmx_bool bPeriodicMols;
3160 tpx.natoms = state->natoms;
3161 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3162 tpx.fep_state = state->fep_state;
3163 tpx.lambda = state->lambda[efptFEP];
3164 tpx.bIr = (ir != NULL);
3165 tpx.bTop = (mtop != NULL);
3166 tpx.bX = (state->x != NULL);
3167 tpx.bV = (state->v != NULL);
3168 tpx.bF = (f != NULL);
3172 TopOnlyOK = (ir == NULL);
3174 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3179 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3180 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3185 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3186 state->natoms = tpx.natoms;
3187 state->nalloc = tpx.natoms;
3193 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3197 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3199 do_test(fio, tpx.bBox, state->box);
3200 do_section(fio, eitemBOX, bRead);
3203 gmx_fio_ndo_rvec(fio, state->box, DIM);
3204 if (file_version >= 51)
3206 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3210 /* We initialize box_rel after reading the inputrec */
3211 clear_mat(state->box_rel);
3213 if (file_version >= 28)
3215 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3216 if (file_version < 56)
3219 gmx_fio_ndo_rvec(fio, mdum, DIM);
3224 if (state->ngtc > 0 && file_version >= 28)
3227 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3228 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3229 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3230 snew(dumv, state->ngtc);
3231 if (file_version < 69)
3233 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3235 /* These used to be the Berendsen tcoupl_lambda's */
3236 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3240 /* Prior to tpx version 26, the inputrec was here.
3241 * I moved it to enable partial forward-compatibility
3242 * for analysis/viewer programs.
3244 if (file_version < 26)
3246 do_test(fio, tpx.bIr, ir);
3247 do_section(fio, eitemIR, bRead);
3252 do_inputrec(fio, ir, bRead, file_version,
3253 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3256 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3261 do_inputrec(fio, &dum_ir, bRead, file_version,
3262 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3265 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3267 done_inputrec(&dum_ir);
3273 do_test(fio, tpx.bTop, mtop);
3274 do_section(fio, eitemTOP, bRead);
3279 do_mtop(fio, mtop, bRead, file_version);
3283 do_mtop(fio, &dum_top, bRead, file_version);
3284 done_mtop(&dum_top, TRUE);
3287 do_test(fio, tpx.bX, state->x);
3288 do_section(fio, eitemX, bRead);
3293 state->flags |= (1<<estX);
3295 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3298 do_test(fio, tpx.bV, state->v);
3299 do_section(fio, eitemV, bRead);
3304 state->flags |= (1<<estV);
3306 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3309 do_test(fio, tpx.bF, f);
3310 do_section(fio, eitemF, bRead);
3313 gmx_fio_ndo_rvec(fio, f, state->natoms);
3316 /* Starting with tpx version 26, we have the inputrec
3317 * at the end of the file, so we can ignore it
3318 * if the file is never than the software (but still the
3319 * same generation - see comments at the top of this file.
3324 bPeriodicMols = FALSE;
3325 if (file_version >= 26)
3327 do_test(fio, tpx.bIr, ir);
3328 do_section(fio, eitemIR, bRead);
3331 if (file_version >= 53)
3333 /* Removed the pbc info from do_inputrec, since we always want it */
3337 bPeriodicMols = ir->bPeriodicMols;
3339 gmx_fio_do_int(fio, ePBC);
3340 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3342 if (file_generation <= tpx_generation && ir)
3344 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3347 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3349 if (file_version < 51)
3351 set_box_rel(ir, state);
3353 if (file_version < 53)
3356 bPeriodicMols = ir->bPeriodicMols;
3359 if (bRead && ir && file_version >= 53)
3361 /* We need to do this after do_inputrec, since that initializes ir */
3363 ir->bPeriodicMols = bPeriodicMols;
3372 if (state->ngtc == 0)
3374 /* Reading old version without tcoupl state data: set it */
3375 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3377 if (tpx.bTop && mtop)
3379 if (file_version < 57)
3381 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3383 ir->eDisre = edrSimple;
3387 ir->eDisre = edrNone;
3390 set_disres_npair(mtop);
3394 if (tpx.bTop && mtop)
3396 gmx_mtop_finalize(mtop);
3399 if (file_version >= 57)
3403 env = getenv("GMX_NOCHARGEGROUPS");
3406 sscanf(env, "%d", &ienv);
3407 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3412 "Will make single atomic charge groups in non-solvent%s\n",
3413 ienv > 1 ? " and solvent" : "");
3414 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3416 fprintf(stderr, "\n");
3424 /************************************************************
3426 * The following routines are the exported ones
3428 ************************************************************/
3430 t_fileio *open_tpx(const char *fn, const char *mode)
3432 return gmx_fio_open(fn, mode);
3435 void close_tpx(t_fileio *fio)
3440 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3441 int *file_version, int *file_generation)
3445 fio = open_tpx(fn, "r");
3446 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3450 void write_tpx_state(const char *fn,
3451 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3455 fio = open_tpx(fn, "w");
3456 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3460 void read_tpx_state(const char *fn,
3461 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3465 fio = open_tpx(fn, "r");
3466 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3470 int read_tpx(const char *fn,
3471 t_inputrec *ir, matrix box, int *natoms,
3472 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3480 fio = open_tpx(fn, "r");
3481 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3483 *natoms = state.natoms;
3486 copy_mat(state.box, box);
3495 int read_tpx_top(const char *fn,
3496 t_inputrec *ir, matrix box, int *natoms,
3497 rvec *x, rvec *v, rvec *f, t_topology *top)
3503 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3505 *top = gmx_mtop_t_to_t_topology(&mtop);
3510 gmx_bool fn2bTPX(const char *file)
3512 switch (fn2ftp(file))
3523 static void done_gmx_groups_t(gmx_groups_t *g)
3527 for (i = 0; (i < egcNR); i++)
3529 if (NULL != g->grps[i].nm_ind)
3531 sfree(g->grps[i].nm_ind);
3532 g->grps[i].nm_ind = NULL;
3534 if (NULL != g->grpnr[i])
3540 /* The contents of this array is in symtab, don't free it here */
3544 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3545 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3548 int natoms, i, version, generation;
3549 gmx_bool bTop, bXNULL = FALSE;
3551 t_topology *topconv;
3554 bTop = fn2bTPX(infile);
3558 read_tpxheader(infile, &header, TRUE, &version, &generation);
3561 snew(*x, header.natoms);
3565 snew(*v, header.natoms);
3568 *ePBC = read_tpx(infile, NULL, box, &natoms,
3569 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3570 *top = gmx_mtop_t_to_t_topology(mtop);
3571 /* In this case we need to throw away the group data too */
3572 done_gmx_groups_t(&mtop->groups);
3574 strcpy(title, *top->name);
3575 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3579 get_stx_coordnum(infile, &natoms);
3580 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3591 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3599 aps = gmx_atomprop_init();
3600 for (i = 0; (i < natoms); i++)
3602 if (!gmx_atomprop_query(aps, epropMass,
3603 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3604 *top->atoms.atomname[i],
3605 &(top->atoms.atom[i].m)))
3609 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3610 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3611 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3612 *top->atoms.atomname[i]);
3616 gmx_atomprop_destroy(aps);
3618 top->idef.ntypes = -1;