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"
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;
71 /* This number should be increased whenever the file format changes! */
72 static const int tpx_version = 95;
74 /* This number should only be increased when you edit the TOPOLOGY section
75 * or the HEADER of the tpx format.
76 * This way we can maintain forward compatibility too for all analysis tools
77 * and/or external programs that only need to know the atom/residue names,
78 * charges, and bond connectivity.
80 * It first appeared in tpx version 26, when I also moved the inputrecord
81 * to the end of the tpx file, so we can just skip it if we only
84 static const int tpx_generation = 25;
86 /* This number should be the most recent backwards incompatible version
87 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
89 static const int tpx_incompatible_version = 9;
93 /* Struct used to maintain tpx compatibility when function types are added */
95 int fvnr; /* file version number in which the function type first appeared */
96 int ftype; /* function type */
100 * The entries should be ordered in:
101 * 1. ascending file version number
102 * 2. ascending function type number
104 /*static const t_ftupd ftupd[] = {
105 { 20, F_CUBICBONDS },
109 { 22, F_DISRESVIOL },
115 { 26, F_DIHRESVIOL },
116 { 30, F_CROSS_BOND_BONDS },
117 { 30, F_CROSS_BOND_ANGLES },
118 { 30, F_UREY_BRADLEY },
119 { 30, F_POLARIZATION },
123 * The entries should be ordered in:
124 * 1. ascending function type number
125 * 2. ascending file version number
127 /* question; what is the purpose of the commented code above? */
128 static const t_ftupd ftupd[] = {
129 { 20, F_CUBICBONDS },
134 { 43, F_TABBONDSNC },
135 { 70, F_RESTRBONDS },
136 { 76, F_LINEAR_ANGLES },
137 { 30, F_CROSS_BOND_BONDS },
138 { 30, F_CROSS_BOND_ANGLES },
139 { 30, F_UREY_BRADLEY },
140 { 34, F_QUARTIC_ANGLES },
150 { 72, F_NPSOLVATION },
152 { 41, F_LJC_PAIRS_NB },
155 { 32, F_COUL_RECIP },
158 { 30, F_POLARIZATION },
161 { 22, F_DISRESVIOL },
165 { 26, F_DIHRESVIOL },
170 { 46, F_ECONSERVED },
171 { 69, F_VTEMP_NOLONGERUSED},
173 { 54, F_DVDL_CONSTR },
174 { 76, F_ANHARM_POL },
177 { 79, F_DVDL_BONDED, },
178 { 79, F_DVDL_RESTRAINT },
179 { 79, F_DVDL_TEMPERATURE },
181 #define NFTUPD asize(ftupd)
183 /* Needed for backward compatibility */
186 static void _do_section(t_fileio *fio, int key, gmx_bool bRead, const char *src,
192 if (gmx_fio_getftp(fio) == efTPA)
196 gmx_fio_write_string(fio, itemstr[key]);
197 bDbg = gmx_fio_getdebug(fio);
198 gmx_fio_setdebug(fio, FALSE);
199 gmx_fio_write_string(fio, comment_str[key]);
200 gmx_fio_setdebug(fio, bDbg);
204 if (gmx_fio_getdebug(fio))
206 fprintf(stderr, "Looking for section %s (%s, %d)",
207 itemstr[key], src, line);
212 gmx_fio_do_string(fio, buf);
214 while ((gmx_strcasecmp(buf, itemstr[key]) != 0));
216 if (gmx_strcasecmp(buf, itemstr[key]) != 0)
218 gmx_fatal(FARGS, "\nCould not find section heading %s", itemstr[key]);
220 else if (gmx_fio_getdebug(fio))
222 fprintf(stderr, " and found it\n");
228 #define do_section(fio, key, bRead) _do_section(fio, key, bRead, __FILE__, __LINE__)
230 /**************************************************************
232 * Now the higer level routines that do io of the structures and arrays
234 **************************************************************/
235 static void do_pullgrp_tpx_pre95(t_fileio *fio,
244 gmx_fio_do_int(fio, pgrp->nat);
247 snew(pgrp->ind, pgrp->nat);
249 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
250 gmx_fio_do_int(fio, pgrp->nweight);
253 snew(pgrp->weight, pgrp->nweight);
255 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
256 gmx_fio_do_int(fio, pgrp->pbcatom);
257 gmx_fio_do_rvec(fio, pcrd->vec);
258 clear_rvec(pcrd->origin);
259 gmx_fio_do_rvec(fio, tmp);
261 gmx_fio_do_real(fio, pcrd->rate);
262 gmx_fio_do_real(fio, pcrd->k);
263 if (file_version >= 56)
265 gmx_fio_do_real(fio, pcrd->kB);
273 static void do_pull_group(t_fileio *fio, t_pull_group *pgrp, gmx_bool bRead)
277 gmx_fio_do_int(fio, pgrp->nat);
280 snew(pgrp->ind, pgrp->nat);
282 gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
283 gmx_fio_do_int(fio, pgrp->nweight);
286 snew(pgrp->weight, pgrp->nweight);
288 gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
289 gmx_fio_do_int(fio, pgrp->pbcatom);
292 static void do_pull_coord(t_fileio *fio, t_pull_coord *pcrd)
296 gmx_fio_do_int(fio, pcrd->group[0]);
297 gmx_fio_do_int(fio, pcrd->group[1]);
298 gmx_fio_do_rvec(fio, pcrd->origin);
299 gmx_fio_do_rvec(fio, pcrd->vec);
300 gmx_fio_do_real(fio, pcrd->init);
301 gmx_fio_do_real(fio, pcrd->rate);
302 gmx_fio_do_real(fio, pcrd->k);
303 gmx_fio_do_real(fio, pcrd->kB);
306 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
308 /* i is used in the ndo_double macro*/
312 int n_lambda = fepvals->n_lambda;
314 /* reset the lambda calculation window */
315 fepvals->lambda_start_n = 0;
316 fepvals->lambda_stop_n = n_lambda;
317 if (file_version >= 79)
323 snew(expand->init_lambda_weights, n_lambda);
325 gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
326 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
329 gmx_fio_do_int(fio, expand->nstexpanded);
330 gmx_fio_do_int(fio, expand->elmcmove);
331 gmx_fio_do_int(fio, expand->elamstats);
332 gmx_fio_do_int(fio, expand->lmc_repeats);
333 gmx_fio_do_int(fio, expand->gibbsdeltalam);
334 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
335 gmx_fio_do_int(fio, expand->lmc_seed);
336 gmx_fio_do_real(fio, expand->mc_temp);
337 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
338 gmx_fio_do_int(fio, expand->nstTij);
339 gmx_fio_do_int(fio, expand->minvarmin);
340 gmx_fio_do_int(fio, expand->c_range);
341 gmx_fio_do_real(fio, expand->wl_scale);
342 gmx_fio_do_real(fio, expand->wl_ratio);
343 gmx_fio_do_real(fio, expand->init_wl_delta);
344 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
345 gmx_fio_do_int(fio, expand->elmceq);
346 gmx_fio_do_int(fio, expand->equil_steps);
347 gmx_fio_do_int(fio, expand->equil_samples);
348 gmx_fio_do_int(fio, expand->equil_n_at_lam);
349 gmx_fio_do_real(fio, expand->equil_wl_delta);
350 gmx_fio_do_real(fio, expand->equil_ratio);
354 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
357 if (file_version >= 79)
359 gmx_fio_do_int(fio, simtemp->eSimTempScale);
360 gmx_fio_do_real(fio, simtemp->simtemp_high);
361 gmx_fio_do_real(fio, simtemp->simtemp_low);
366 snew(simtemp->temperatures, n_lambda);
368 gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
373 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
375 /* i is defined in the ndo_double macro; use g to iterate. */
380 /* free energy values */
382 if (file_version >= 79)
384 gmx_fio_do_int(fio, fepvals->init_fep_state);
385 gmx_fio_do_double(fio, fepvals->init_lambda);
386 gmx_fio_do_double(fio, fepvals->delta_lambda);
388 else if (file_version >= 59)
390 gmx_fio_do_double(fio, fepvals->init_lambda);
391 gmx_fio_do_double(fio, fepvals->delta_lambda);
395 gmx_fio_do_real(fio, rdum);
396 fepvals->init_lambda = rdum;
397 gmx_fio_do_real(fio, rdum);
398 fepvals->delta_lambda = rdum;
400 if (file_version >= 79)
402 gmx_fio_do_int(fio, fepvals->n_lambda);
405 snew(fepvals->all_lambda, efptNR);
407 for (g = 0; g < efptNR; g++)
409 if (fepvals->n_lambda > 0)
413 snew(fepvals->all_lambda[g], fepvals->n_lambda);
415 gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
416 gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
418 else if (fepvals->init_lambda >= 0)
420 fepvals->separate_dvdl[efptFEP] = TRUE;
424 else if (file_version >= 64)
426 gmx_fio_do_int(fio, fepvals->n_lambda);
431 snew(fepvals->all_lambda, efptNR);
432 /* still allocate the all_lambda array's contents. */
433 for (g = 0; g < efptNR; g++)
435 if (fepvals->n_lambda > 0)
437 snew(fepvals->all_lambda[g], fepvals->n_lambda);
441 gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
443 if (fepvals->init_lambda >= 0)
447 fepvals->separate_dvdl[efptFEP] = TRUE;
451 /* copy the contents of the efptFEP lambda component to all
452 the other components */
453 for (g = 0; g < efptNR; g++)
455 for (h = 0; h < fepvals->n_lambda; h++)
459 fepvals->all_lambda[g][h] =
460 fepvals->all_lambda[efptFEP][h];
469 fepvals->n_lambda = 0;
470 fepvals->all_lambda = NULL;
471 if (fepvals->init_lambda >= 0)
473 fepvals->separate_dvdl[efptFEP] = TRUE;
476 if (file_version >= 13)
478 gmx_fio_do_real(fio, fepvals->sc_alpha);
482 fepvals->sc_alpha = 0;
484 if (file_version >= 38)
486 gmx_fio_do_int(fio, fepvals->sc_power);
490 fepvals->sc_power = 2;
492 if (file_version >= 79)
494 gmx_fio_do_real(fio, fepvals->sc_r_power);
498 fepvals->sc_r_power = 6.0;
500 if (file_version >= 15)
502 gmx_fio_do_real(fio, fepvals->sc_sigma);
506 fepvals->sc_sigma = 0.3;
510 if (file_version >= 71)
512 fepvals->sc_sigma_min = fepvals->sc_sigma;
516 fepvals->sc_sigma_min = 0;
519 if (file_version >= 79)
521 gmx_fio_do_int(fio, fepvals->bScCoul);
525 fepvals->bScCoul = TRUE;
527 if (file_version >= 64)
529 gmx_fio_do_int(fio, fepvals->nstdhdl);
533 fepvals->nstdhdl = 1;
536 if (file_version >= 73)
538 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
539 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
543 fepvals->separate_dhdl_file = esepdhdlfileYES;
544 fepvals->dhdl_derivatives = edhdlderivativesYES;
546 if (file_version >= 71)
548 gmx_fio_do_int(fio, fepvals->dh_hist_size);
549 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
553 fepvals->dh_hist_size = 0;
554 fepvals->dh_hist_spacing = 0.1;
556 if (file_version >= 79)
558 gmx_fio_do_int(fio, fepvals->bPrintEnergy);
562 fepvals->bPrintEnergy = FALSE;
565 /* handle lambda_neighbors */
566 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
568 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
569 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
570 (fepvals->init_lambda < 0) )
572 fepvals->lambda_start_n = (fepvals->init_fep_state -
573 fepvals->lambda_neighbors);
574 fepvals->lambda_stop_n = (fepvals->init_fep_state +
575 fepvals->lambda_neighbors + 1);
576 if (fepvals->lambda_start_n < 0)
578 fepvals->lambda_start_n = 0;;
580 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
582 fepvals->lambda_stop_n = fepvals->n_lambda;
587 fepvals->lambda_start_n = 0;
588 fepvals->lambda_stop_n = fepvals->n_lambda;
593 fepvals->lambda_start_n = 0;
594 fepvals->lambda_stop_n = fepvals->n_lambda;
598 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
602 if (file_version >= 95)
604 gmx_fio_do_int(fio, pull->ngroup);
606 gmx_fio_do_int(fio, pull->ncoord);
607 if (file_version < 95)
609 pull->ngroup = pull->ncoord + 1;
611 gmx_fio_do_int(fio, pull->eGeom);
612 gmx_fio_do_ivec(fio, pull->dim);
613 gmx_fio_do_real(fio, pull->cyl_r1);
614 gmx_fio_do_real(fio, pull->cyl_r0);
615 gmx_fio_do_real(fio, pull->constr_tol);
616 if (file_version >= 95)
618 gmx_fio_do_int(fio, pull->bPrintRef);
620 gmx_fio_do_int(fio, pull->nstxout);
621 gmx_fio_do_int(fio, pull->nstfout);
624 snew(pull->group, pull->ngroup);
625 snew(pull->coord, pull->ncoord);
627 if (file_version < 95)
629 /* epullgPOS for position pulling, before epullgDIRPBC was removed */
630 if (pull->eGeom == epullgDIRPBC)
632 gmx_fatal(FARGS, "pull-geometry=position is no longer supported");
634 if (pull->eGeom > epullgDIRPBC)
639 for (g = 0; g < pull->ngroup; g++)
641 /* We read and ignore a pull coordinate for group 0 */
642 do_pullgrp_tpx_pre95(fio, &pull->group[g], &pull->coord[max(g-1,0)],
643 bRead, file_version);
646 pull->coord[g-1].group[0] = 0;
647 pull->coord[g-1].group[1] = g;
651 pull->bPrintRef = (pull->group[0].nat > 0);
655 for (g = 0; g < pull->ngroup; g++)
657 do_pull_group(fio, &pull->group[g], bRead);
659 for (g = 0; g < pull->ncoord; g++)
661 do_pull_coord(fio, &pull->coord[g]);
667 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead)
671 gmx_fio_do_int(fio, rotg->eType);
672 gmx_fio_do_int(fio, rotg->bMassW);
673 gmx_fio_do_int(fio, rotg->nat);
676 snew(rotg->ind, rotg->nat);
678 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
681 snew(rotg->x_ref, rotg->nat);
683 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
684 gmx_fio_do_rvec(fio, rotg->vec);
685 gmx_fio_do_rvec(fio, rotg->pivot);
686 gmx_fio_do_real(fio, rotg->rate);
687 gmx_fio_do_real(fio, rotg->k);
688 gmx_fio_do_real(fio, rotg->slab_dist);
689 gmx_fio_do_real(fio, rotg->min_gaussian);
690 gmx_fio_do_real(fio, rotg->eps);
691 gmx_fio_do_int(fio, rotg->eFittype);
692 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
693 gmx_fio_do_real(fio, rotg->PotAngle_step);
696 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead)
700 gmx_fio_do_int(fio, rot->ngrp);
701 gmx_fio_do_int(fio, rot->nstrout);
702 gmx_fio_do_int(fio, rot->nstsout);
705 snew(rot->grp, rot->ngrp);
707 for (g = 0; g < rot->ngrp; g++)
709 do_rotgrp(fio, &rot->grp[g], bRead);
714 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
715 int file_version, real *fudgeQQ)
717 int i, j, k, *tmp, idum = 0;
721 real zerotemptime, finish_t, init_temp, finish_temp;
723 if (file_version != tpx_version)
725 /* Give a warning about features that are not accessible */
726 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
727 file_version, tpx_version);
735 if (file_version == 0)
740 /* Basic inputrec stuff */
741 gmx_fio_do_int(fio, ir->eI);
742 if (file_version >= 62)
744 gmx_fio_do_int64(fio, ir->nsteps);
748 gmx_fio_do_int(fio, idum);
751 if (file_version > 25)
753 if (file_version >= 62)
755 gmx_fio_do_int64(fio, ir->init_step);
759 gmx_fio_do_int(fio, idum);
760 ir->init_step = idum;
768 if (file_version >= 58)
770 gmx_fio_do_int(fio, ir->simulation_part);
774 ir->simulation_part = 1;
777 if (file_version >= 67)
779 gmx_fio_do_int(fio, ir->nstcalcenergy);
783 ir->nstcalcenergy = 1;
785 if (file_version < 53)
787 /* The pbc info has been moved out of do_inputrec,
788 * since we always want it, also without reading the inputrec.
790 gmx_fio_do_int(fio, ir->ePBC);
791 if ((file_version <= 15) && (ir->ePBC == 2))
795 if (file_version >= 45)
797 gmx_fio_do_int(fio, ir->bPeriodicMols);
804 ir->bPeriodicMols = TRUE;
808 ir->bPeriodicMols = FALSE;
812 if (file_version >= 81)
814 gmx_fio_do_int(fio, ir->cutoff_scheme);
815 if (file_version < 94)
817 ir->cutoff_scheme = 1 - ir->cutoff_scheme;
822 ir->cutoff_scheme = ecutsGROUP;
824 gmx_fio_do_int(fio, ir->ns_type);
825 gmx_fio_do_int(fio, ir->nstlist);
826 gmx_fio_do_int(fio, ir->ndelta);
827 if (file_version < 41)
829 gmx_fio_do_int(fio, idum);
830 gmx_fio_do_int(fio, idum);
832 if (file_version >= 45)
834 gmx_fio_do_real(fio, ir->rtpi);
840 gmx_fio_do_int(fio, ir->nstcomm);
841 if (file_version > 34)
843 gmx_fio_do_int(fio, ir->comm_mode);
845 else if (ir->nstcomm < 0)
847 ir->comm_mode = ecmANGULAR;
851 ir->comm_mode = ecmLINEAR;
853 ir->nstcomm = abs(ir->nstcomm);
855 if (file_version > 25)
857 gmx_fio_do_int(fio, ir->nstcheckpoint);
861 ir->nstcheckpoint = 0;
864 gmx_fio_do_int(fio, ir->nstcgsteep);
866 if (file_version >= 30)
868 gmx_fio_do_int(fio, ir->nbfgscorr);
875 gmx_fio_do_int(fio, ir->nstlog);
876 gmx_fio_do_int(fio, ir->nstxout);
877 gmx_fio_do_int(fio, ir->nstvout);
878 gmx_fio_do_int(fio, ir->nstfout);
879 gmx_fio_do_int(fio, ir->nstenergy);
880 gmx_fio_do_int(fio, ir->nstxtcout);
881 if (file_version >= 59)
883 gmx_fio_do_double(fio, ir->init_t);
884 gmx_fio_do_double(fio, ir->delta_t);
888 gmx_fio_do_real(fio, rdum);
890 gmx_fio_do_real(fio, rdum);
893 gmx_fio_do_real(fio, ir->xtcprec);
894 if (file_version < 19)
896 gmx_fio_do_int(fio, idum);
897 gmx_fio_do_int(fio, idum);
899 if (file_version < 18)
901 gmx_fio_do_int(fio, idum);
903 if (file_version >= 81)
905 gmx_fio_do_real(fio, ir->verletbuf_tol);
909 ir->verletbuf_tol = 0;
911 gmx_fio_do_real(fio, ir->rlist);
912 if (file_version >= 67)
914 gmx_fio_do_real(fio, ir->rlistlong);
916 if (file_version >= 82 && file_version != 90)
918 gmx_fio_do_int(fio, ir->nstcalclr);
922 /* Calculate at NS steps */
923 ir->nstcalclr = ir->nstlist;
925 gmx_fio_do_int(fio, ir->coulombtype);
926 if (file_version < 32 && ir->coulombtype == eelRF)
928 ir->coulombtype = eelRF_NEC;
930 if (file_version >= 81)
932 gmx_fio_do_int(fio, ir->coulomb_modifier);
936 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
938 gmx_fio_do_real(fio, ir->rcoulomb_switch);
939 gmx_fio_do_real(fio, ir->rcoulomb);
940 gmx_fio_do_int(fio, ir->vdwtype);
941 if (file_version >= 81)
943 gmx_fio_do_int(fio, ir->vdw_modifier);
947 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
949 gmx_fio_do_real(fio, ir->rvdw_switch);
950 gmx_fio_do_real(fio, ir->rvdw);
951 if (file_version < 67)
953 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
955 gmx_fio_do_int(fio, ir->eDispCorr);
956 gmx_fio_do_real(fio, ir->epsilon_r);
957 if (file_version >= 37)
959 gmx_fio_do_real(fio, ir->epsilon_rf);
963 if (EEL_RF(ir->coulombtype))
965 ir->epsilon_rf = ir->epsilon_r;
970 ir->epsilon_rf = 1.0;
973 if (file_version >= 29)
975 gmx_fio_do_real(fio, ir->tabext);
982 if (file_version > 25)
984 gmx_fio_do_int(fio, ir->gb_algorithm);
985 gmx_fio_do_int(fio, ir->nstgbradii);
986 gmx_fio_do_real(fio, ir->rgbradii);
987 gmx_fio_do_real(fio, ir->gb_saltconc);
988 gmx_fio_do_int(fio, ir->implicit_solvent);
992 ir->gb_algorithm = egbSTILL;
996 ir->implicit_solvent = eisNO;
998 if (file_version >= 55)
1000 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
1001 gmx_fio_do_real(fio, ir->gb_obc_alpha);
1002 gmx_fio_do_real(fio, ir->gb_obc_beta);
1003 gmx_fio_do_real(fio, ir->gb_obc_gamma);
1004 if (file_version >= 60)
1006 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
1007 gmx_fio_do_int(fio, ir->sa_algorithm);
1011 ir->gb_dielectric_offset = 0.009;
1012 ir->sa_algorithm = esaAPPROX;
1014 gmx_fio_do_real(fio, ir->sa_surface_tension);
1016 /* Override sa_surface_tension if it is not changed in the mpd-file */
1017 if (ir->sa_surface_tension < 0)
1019 if (ir->gb_algorithm == egbSTILL)
1021 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
1023 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1025 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
1032 /* Better use sensible values than insane (0.0) ones... */
1033 ir->gb_epsilon_solvent = 80;
1034 ir->gb_obc_alpha = 1.0;
1035 ir->gb_obc_beta = 0.8;
1036 ir->gb_obc_gamma = 4.85;
1037 ir->sa_surface_tension = 2.092;
1041 if (file_version >= 81)
1043 gmx_fio_do_real(fio, ir->fourier_spacing);
1047 ir->fourier_spacing = 0.0;
1049 gmx_fio_do_int(fio, ir->nkx);
1050 gmx_fio_do_int(fio, ir->nky);
1051 gmx_fio_do_int(fio, ir->nkz);
1052 gmx_fio_do_int(fio, ir->pme_order);
1053 gmx_fio_do_real(fio, ir->ewald_rtol);
1055 if (file_version >= 93)
1057 gmx_fio_do_real(fio, ir->ewald_rtol_lj);
1061 ir->ewald_rtol_lj = ir->ewald_rtol;
1064 if (file_version >= 24)
1066 gmx_fio_do_int(fio, ir->ewald_geometry);
1070 ir->ewald_geometry = eewg3D;
1073 if (file_version <= 17)
1075 ir->epsilon_surface = 0;
1076 if (file_version == 17)
1078 gmx_fio_do_int(fio, idum);
1083 gmx_fio_do_real(fio, ir->epsilon_surface);
1086 gmx_fio_do_gmx_bool(fio, ir->bOptFFT);
1088 if (file_version >= 93)
1090 gmx_fio_do_int(fio, ir->ljpme_combination_rule);
1092 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1093 gmx_fio_do_int(fio, ir->etc);
1094 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1095 * but the values 0 and 1 still mean no and
1096 * berendsen temperature coupling, respectively.
1098 if (file_version >= 79)
1100 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1102 if (file_version >= 71)
1104 gmx_fio_do_int(fio, ir->nsttcouple);
1108 ir->nsttcouple = ir->nstcalcenergy;
1110 if (file_version <= 15)
1112 gmx_fio_do_int(fio, idum);
1114 if (file_version <= 17)
1116 gmx_fio_do_int(fio, ir->epct);
1117 if (file_version <= 15)
1121 ir->epct = epctSURFACETENSION;
1123 gmx_fio_do_int(fio, idum);
1126 /* we have removed the NO alternative at the beginning */
1130 ir->epct = epctISOTROPIC;
1134 ir->epc = epcBERENDSEN;
1139 gmx_fio_do_int(fio, ir->epc);
1140 gmx_fio_do_int(fio, ir->epct);
1142 if (file_version >= 71)
1144 gmx_fio_do_int(fio, ir->nstpcouple);
1148 ir->nstpcouple = ir->nstcalcenergy;
1150 gmx_fio_do_real(fio, ir->tau_p);
1151 if (file_version <= 15)
1153 gmx_fio_do_rvec(fio, vdum);
1154 clear_mat(ir->ref_p);
1155 for (i = 0; i < DIM; i++)
1157 ir->ref_p[i][i] = vdum[i];
1162 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1163 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1164 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1166 if (file_version <= 15)
1168 gmx_fio_do_rvec(fio, vdum);
1169 clear_mat(ir->compress);
1170 for (i = 0; i < DIM; i++)
1172 ir->compress[i][i] = vdum[i];
1177 gmx_fio_do_rvec(fio, ir->compress[XX]);
1178 gmx_fio_do_rvec(fio, ir->compress[YY]);
1179 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1181 if (file_version >= 47)
1183 gmx_fio_do_int(fio, ir->refcoord_scaling);
1184 gmx_fio_do_rvec(fio, ir->posres_com);
1185 gmx_fio_do_rvec(fio, ir->posres_comB);
1189 ir->refcoord_scaling = erscNO;
1190 clear_rvec(ir->posres_com);
1191 clear_rvec(ir->posres_comB);
1193 if ((file_version > 25) && (file_version < 79))
1195 gmx_fio_do_int(fio, ir->andersen_seed);
1199 ir->andersen_seed = 0;
1201 if (file_version < 26)
1203 gmx_fio_do_gmx_bool(fio, bSimAnn);
1204 gmx_fio_do_real(fio, zerotemptime);
1207 if (file_version < 37)
1209 gmx_fio_do_real(fio, rdum);
1212 gmx_fio_do_real(fio, ir->shake_tol);
1213 if (file_version < 54)
1215 gmx_fio_do_real(fio, *fudgeQQ);
1218 gmx_fio_do_int(fio, ir->efep);
1219 if (file_version <= 14 && ir->efep != efepNO)
1223 do_fepvals(fio, ir->fepvals, bRead, file_version);
1225 if (file_version >= 79)
1227 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1230 ir->bSimTemp = TRUE;
1235 ir->bSimTemp = FALSE;
1239 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1242 if (file_version >= 79)
1244 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1247 ir->bExpanded = TRUE;
1251 ir->bExpanded = FALSE;
1256 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1258 if (file_version >= 57)
1260 gmx_fio_do_int(fio, ir->eDisre);
1262 gmx_fio_do_int(fio, ir->eDisreWeighting);
1263 if (file_version < 22)
1265 if (ir->eDisreWeighting == 0)
1267 ir->eDisreWeighting = edrwEqual;
1271 ir->eDisreWeighting = edrwConservative;
1274 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1275 gmx_fio_do_real(fio, ir->dr_fc);
1276 gmx_fio_do_real(fio, ir->dr_tau);
1277 gmx_fio_do_int(fio, ir->nstdisreout);
1278 if (file_version >= 22)
1280 gmx_fio_do_real(fio, ir->orires_fc);
1281 gmx_fio_do_real(fio, ir->orires_tau);
1282 gmx_fio_do_int(fio, ir->nstorireout);
1288 ir->nstorireout = 0;
1290 if (file_version >= 26 && file_version < 79)
1292 gmx_fio_do_real(fio, ir->dihre_fc);
1293 if (file_version < 56)
1295 gmx_fio_do_real(fio, rdum);
1296 gmx_fio_do_int(fio, idum);
1304 gmx_fio_do_real(fio, ir->em_stepsize);
1305 gmx_fio_do_real(fio, ir->em_tol);
1306 if (file_version >= 22)
1308 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1312 ir->bShakeSOR = TRUE;
1314 if (file_version >= 11)
1316 gmx_fio_do_int(fio, ir->niter);
1321 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1324 if (file_version >= 21)
1326 gmx_fio_do_real(fio, ir->fc_stepsize);
1330 ir->fc_stepsize = 0;
1332 gmx_fio_do_int(fio, ir->eConstrAlg);
1333 gmx_fio_do_int(fio, ir->nProjOrder);
1334 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1335 if (file_version <= 14)
1337 gmx_fio_do_int(fio, idum);
1339 if (file_version >= 26)
1341 gmx_fio_do_int(fio, ir->nLincsIter);
1346 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1349 if (file_version < 33)
1351 gmx_fio_do_real(fio, bd_temp);
1353 gmx_fio_do_real(fio, ir->bd_fric);
1354 gmx_fio_do_int(fio, ir->ld_seed);
1355 if (file_version >= 33)
1357 for (i = 0; i < DIM; i++)
1359 gmx_fio_do_rvec(fio, ir->deform[i]);
1364 for (i = 0; i < DIM; i++)
1366 clear_rvec(ir->deform[i]);
1369 if (file_version >= 14)
1371 gmx_fio_do_real(fio, ir->cos_accel);
1377 gmx_fio_do_int(fio, ir->userint1);
1378 gmx_fio_do_int(fio, ir->userint2);
1379 gmx_fio_do_int(fio, ir->userint3);
1380 gmx_fio_do_int(fio, ir->userint4);
1381 gmx_fio_do_real(fio, ir->userreal1);
1382 gmx_fio_do_real(fio, ir->userreal2);
1383 gmx_fio_do_real(fio, ir->userreal3);
1384 gmx_fio_do_real(fio, ir->userreal4);
1387 if (file_version >= 77)
1389 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1394 snew(ir->adress, 1);
1396 gmx_fio_do_int(fio, ir->adress->type);
1397 gmx_fio_do_real(fio, ir->adress->const_wf);
1398 gmx_fio_do_real(fio, ir->adress->ex_width);
1399 gmx_fio_do_real(fio, ir->adress->hy_width);
1400 gmx_fio_do_int(fio, ir->adress->icor);
1401 gmx_fio_do_int(fio, ir->adress->site);
1402 gmx_fio_do_rvec(fio, ir->adress->refs);
1403 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1404 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1405 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1406 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1410 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1412 if (ir->adress->n_tf_grps > 0)
1414 gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1418 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1420 if (ir->adress->n_energy_grps > 0)
1422 gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1428 ir->bAdress = FALSE;
1432 if (file_version >= 48)
1434 gmx_fio_do_int(fio, ir->ePull);
1435 if (ir->ePull != epullNO)
1441 do_pull(fio, ir->pull, bRead, file_version);
1446 ir->ePull = epullNO;
1449 /* Enforced rotation */
1450 if (file_version >= 74)
1452 gmx_fio_do_int(fio, ir->bRot);
1453 if (ir->bRot == TRUE)
1459 do_rot(fio, ir->rot, bRead);
1468 gmx_fio_do_int(fio, ir->opts.ngtc);
1469 if (file_version >= 69)
1471 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1475 ir->opts.nhchainlength = 1;
1477 gmx_fio_do_int(fio, ir->opts.ngacc);
1478 gmx_fio_do_int(fio, ir->opts.ngfrz);
1479 gmx_fio_do_int(fio, ir->opts.ngener);
1483 snew(ir->opts.nrdf, ir->opts.ngtc);
1484 snew(ir->opts.ref_t, ir->opts.ngtc);
1485 snew(ir->opts.annealing, ir->opts.ngtc);
1486 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1487 snew(ir->opts.anneal_time, ir->opts.ngtc);
1488 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1489 snew(ir->opts.tau_t, ir->opts.ngtc);
1490 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1491 snew(ir->opts.acc, ir->opts.ngacc);
1492 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1494 if (ir->opts.ngtc > 0)
1496 if (bRead && file_version < 13)
1498 snew(tmp, ir->opts.ngtc);
1499 gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1500 for (i = 0; i < ir->opts.ngtc; i++)
1502 ir->opts.nrdf[i] = tmp[i];
1508 gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1510 gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1511 gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1512 if (file_version < 33 && ir->eI == eiBD)
1514 for (i = 0; i < ir->opts.ngtc; i++)
1516 ir->opts.tau_t[i] = bd_temp;
1520 if (ir->opts.ngfrz > 0)
1522 gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1524 if (ir->opts.ngacc > 0)
1526 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1528 if (file_version >= 12)
1530 gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1531 ir->opts.ngener*ir->opts.ngener);
1534 if (bRead && file_version < 26)
1536 for (i = 0; i < ir->opts.ngtc; i++)
1540 ir->opts.annealing[i] = eannSINGLE;
1541 ir->opts.anneal_npoints[i] = 2;
1542 snew(ir->opts.anneal_time[i], 2);
1543 snew(ir->opts.anneal_temp[i], 2);
1544 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1545 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1546 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1547 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1548 ir->opts.anneal_time[i][0] = ir->init_t;
1549 ir->opts.anneal_time[i][1] = finish_t;
1550 ir->opts.anneal_temp[i][0] = init_temp;
1551 ir->opts.anneal_temp[i][1] = finish_temp;
1555 ir->opts.annealing[i] = eannNO;
1556 ir->opts.anneal_npoints[i] = 0;
1562 /* file version 26 or later */
1563 /* First read the lists with annealing and npoints for each group */
1564 gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1565 gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1566 for (j = 0; j < (ir->opts.ngtc); j++)
1568 k = ir->opts.anneal_npoints[j];
1571 snew(ir->opts.anneal_time[j], k);
1572 snew(ir->opts.anneal_temp[j], k);
1574 gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1575 gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1579 if (file_version >= 45)
1581 gmx_fio_do_int(fio, ir->nwall);
1582 gmx_fio_do_int(fio, ir->wall_type);
1583 if (file_version >= 50)
1585 gmx_fio_do_real(fio, ir->wall_r_linpot);
1589 ir->wall_r_linpot = -1;
1591 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1592 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1593 gmx_fio_do_real(fio, ir->wall_density[0]);
1594 gmx_fio_do_real(fio, ir->wall_density[1]);
1595 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1601 ir->wall_atomtype[0] = -1;
1602 ir->wall_atomtype[1] = -1;
1603 ir->wall_density[0] = 0;
1604 ir->wall_density[1] = 0;
1605 ir->wall_ewald_zfac = 3;
1607 /* Cosine stuff for electric fields */
1608 for (j = 0; (j < DIM); j++)
1610 gmx_fio_do_int(fio, ir->ex[j].n);
1611 gmx_fio_do_int(fio, ir->et[j].n);
1614 snew(ir->ex[j].a, ir->ex[j].n);
1615 snew(ir->ex[j].phi, ir->ex[j].n);
1616 snew(ir->et[j].a, ir->et[j].n);
1617 snew(ir->et[j].phi, ir->et[j].n);
1619 gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1620 gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1621 gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1622 gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1626 if (file_version >= 39)
1628 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1629 gmx_fio_do_int(fio, ir->QMMMscheme);
1630 gmx_fio_do_real(fio, ir->scalefactor);
1631 gmx_fio_do_int(fio, ir->opts.ngQM);
1634 snew(ir->opts.QMmethod, ir->opts.ngQM);
1635 snew(ir->opts.QMbasis, ir->opts.ngQM);
1636 snew(ir->opts.QMcharge, ir->opts.ngQM);
1637 snew(ir->opts.QMmult, ir->opts.ngQM);
1638 snew(ir->opts.bSH, ir->opts.ngQM);
1639 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1640 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1641 snew(ir->opts.SAon, ir->opts.ngQM);
1642 snew(ir->opts.SAoff, ir->opts.ngQM);
1643 snew(ir->opts.SAsteps, ir->opts.ngQM);
1644 snew(ir->opts.bOPT, ir->opts.ngQM);
1645 snew(ir->opts.bTS, ir->opts.ngQM);
1647 if (ir->opts.ngQM > 0)
1649 gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1650 gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1651 gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1652 gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1653 gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1654 gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1655 gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1656 gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1657 gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1658 gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1659 gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1660 gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1662 /* end of QMMM stuff */
1667 static void do_harm(t_fileio *fio, t_iparams *iparams)
1669 gmx_fio_do_real(fio, iparams->harmonic.rA);
1670 gmx_fio_do_real(fio, iparams->harmonic.krA);
1671 gmx_fio_do_real(fio, iparams->harmonic.rB);
1672 gmx_fio_do_real(fio, iparams->harmonic.krB);
1675 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1676 gmx_bool bRead, int file_version)
1683 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1693 do_harm(fio, iparams);
1694 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1696 /* Correct incorrect storage of parameters */
1697 iparams->pdihs.phiB = iparams->pdihs.phiA;
1698 iparams->pdihs.cpB = iparams->pdihs.cpA;
1701 case F_LINEAR_ANGLES:
1702 gmx_fio_do_real(fio, iparams->linangle.klinA);
1703 gmx_fio_do_real(fio, iparams->linangle.aA);
1704 gmx_fio_do_real(fio, iparams->linangle.klinB);
1705 gmx_fio_do_real(fio, iparams->linangle.aB);
1708 gmx_fio_do_real(fio, iparams->fene.bm);
1709 gmx_fio_do_real(fio, iparams->fene.kb);
1712 gmx_fio_do_real(fio, iparams->restraint.lowA);
1713 gmx_fio_do_real(fio, iparams->restraint.up1A);
1714 gmx_fio_do_real(fio, iparams->restraint.up2A);
1715 gmx_fio_do_real(fio, iparams->restraint.kA);
1716 gmx_fio_do_real(fio, iparams->restraint.lowB);
1717 gmx_fio_do_real(fio, iparams->restraint.up1B);
1718 gmx_fio_do_real(fio, iparams->restraint.up2B);
1719 gmx_fio_do_real(fio, iparams->restraint.kB);
1725 gmx_fio_do_real(fio, iparams->tab.kA);
1726 gmx_fio_do_int(fio, iparams->tab.table);
1727 gmx_fio_do_real(fio, iparams->tab.kB);
1729 case F_CROSS_BOND_BONDS:
1730 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1731 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1732 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1734 case F_CROSS_BOND_ANGLES:
1735 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1736 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1737 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1738 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1740 case F_UREY_BRADLEY:
1741 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1742 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1743 gmx_fio_do_real(fio, iparams->u_b.r13A);
1744 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1745 if (file_version >= 79)
1747 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1748 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1749 gmx_fio_do_real(fio, iparams->u_b.r13B);
1750 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1754 iparams->u_b.thetaB = iparams->u_b.thetaA;
1755 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1756 iparams->u_b.r13B = iparams->u_b.r13A;
1757 iparams->u_b.kUBB = iparams->u_b.kUBA;
1760 case F_QUARTIC_ANGLES:
1761 gmx_fio_do_real(fio, iparams->qangle.theta);
1762 gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1765 gmx_fio_do_real(fio, iparams->bham.a);
1766 gmx_fio_do_real(fio, iparams->bham.b);
1767 gmx_fio_do_real(fio, iparams->bham.c);
1770 gmx_fio_do_real(fio, iparams->morse.b0A);
1771 gmx_fio_do_real(fio, iparams->morse.cbA);
1772 gmx_fio_do_real(fio, iparams->morse.betaA);
1773 if (file_version >= 79)
1775 gmx_fio_do_real(fio, iparams->morse.b0B);
1776 gmx_fio_do_real(fio, iparams->morse.cbB);
1777 gmx_fio_do_real(fio, iparams->morse.betaB);
1781 iparams->morse.b0B = iparams->morse.b0A;
1782 iparams->morse.cbB = iparams->morse.cbA;
1783 iparams->morse.betaB = iparams->morse.betaA;
1787 gmx_fio_do_real(fio, iparams->cubic.b0);
1788 gmx_fio_do_real(fio, iparams->cubic.kb);
1789 gmx_fio_do_real(fio, iparams->cubic.kcub);
1793 case F_POLARIZATION:
1794 gmx_fio_do_real(fio, iparams->polarize.alpha);
1797 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1798 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1799 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1802 if (file_version < 31)
1804 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1806 gmx_fio_do_real(fio, iparams->wpol.al_x);
1807 gmx_fio_do_real(fio, iparams->wpol.al_y);
1808 gmx_fio_do_real(fio, iparams->wpol.al_z);
1809 gmx_fio_do_real(fio, iparams->wpol.rOH);
1810 gmx_fio_do_real(fio, iparams->wpol.rHH);
1811 gmx_fio_do_real(fio, iparams->wpol.rOD);
1814 gmx_fio_do_real(fio, iparams->thole.a);
1815 gmx_fio_do_real(fio, iparams->thole.alpha1);
1816 gmx_fio_do_real(fio, iparams->thole.alpha2);
1817 gmx_fio_do_real(fio, iparams->thole.rfac);
1820 gmx_fio_do_real(fio, iparams->lj.c6);
1821 gmx_fio_do_real(fio, iparams->lj.c12);
1824 gmx_fio_do_real(fio, iparams->lj14.c6A);
1825 gmx_fio_do_real(fio, iparams->lj14.c12A);
1826 gmx_fio_do_real(fio, iparams->lj14.c6B);
1827 gmx_fio_do_real(fio, iparams->lj14.c12B);
1830 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1831 gmx_fio_do_real(fio, iparams->ljc14.qi);
1832 gmx_fio_do_real(fio, iparams->ljc14.qj);
1833 gmx_fio_do_real(fio, iparams->ljc14.c6);
1834 gmx_fio_do_real(fio, iparams->ljc14.c12);
1836 case F_LJC_PAIRS_NB:
1837 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1838 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1839 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1840 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1846 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1847 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1848 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1850 /* Read the incorrectly stored multiplicity */
1851 gmx_fio_do_real(fio, iparams->harmonic.rB);
1852 gmx_fio_do_real(fio, iparams->harmonic.krB);
1853 iparams->pdihs.phiB = iparams->pdihs.phiA;
1854 iparams->pdihs.cpB = iparams->pdihs.cpA;
1858 gmx_fio_do_real(fio, iparams->pdihs.phiB);
1859 gmx_fio_do_real(fio, iparams->pdihs.cpB);
1860 gmx_fio_do_int(fio, iparams->pdihs.mult);
1864 gmx_fio_do_int(fio, iparams->disres.label);
1865 gmx_fio_do_int(fio, iparams->disres.type);
1866 gmx_fio_do_real(fio, iparams->disres.low);
1867 gmx_fio_do_real(fio, iparams->disres.up1);
1868 gmx_fio_do_real(fio, iparams->disres.up2);
1869 gmx_fio_do_real(fio, iparams->disres.kfac);
1872 gmx_fio_do_int(fio, iparams->orires.ex);
1873 gmx_fio_do_int(fio, iparams->orires.label);
1874 gmx_fio_do_int(fio, iparams->orires.power);
1875 gmx_fio_do_real(fio, iparams->orires.c);
1876 gmx_fio_do_real(fio, iparams->orires.obs);
1877 gmx_fio_do_real(fio, iparams->orires.kfac);
1880 if (file_version < 82)
1882 gmx_fio_do_int(fio, idum);
1883 gmx_fio_do_int(fio, idum);
1885 gmx_fio_do_real(fio, iparams->dihres.phiA);
1886 gmx_fio_do_real(fio, iparams->dihres.dphiA);
1887 gmx_fio_do_real(fio, iparams->dihres.kfacA);
1888 if (file_version >= 82)
1890 gmx_fio_do_real(fio, iparams->dihres.phiB);
1891 gmx_fio_do_real(fio, iparams->dihres.dphiB);
1892 gmx_fio_do_real(fio, iparams->dihres.kfacB);
1896 iparams->dihres.phiB = iparams->dihres.phiA;
1897 iparams->dihres.dphiB = iparams->dihres.dphiA;
1898 iparams->dihres.kfacB = iparams->dihres.kfacA;
1902 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
1903 gmx_fio_do_rvec(fio, iparams->posres.fcA);
1904 if (bRead && file_version < 27)
1906 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
1907 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
1911 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
1912 gmx_fio_do_rvec(fio, iparams->posres.fcB);
1916 gmx_fio_do_int(fio, iparams->fbposres.geom);
1917 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
1918 gmx_fio_do_real(fio, iparams->fbposres.r);
1919 gmx_fio_do_real(fio, iparams->fbposres.k);
1922 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1923 if (file_version >= 25)
1925 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1929 /* Fourier dihedrals are internally represented
1930 * as Ryckaert-Bellemans since those are faster to compute.
1932 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1933 gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1937 gmx_fio_do_real(fio, iparams->constr.dA);
1938 gmx_fio_do_real(fio, iparams->constr.dB);
1941 gmx_fio_do_real(fio, iparams->settle.doh);
1942 gmx_fio_do_real(fio, iparams->settle.dhh);
1945 gmx_fio_do_real(fio, iparams->vsite.a);
1950 gmx_fio_do_real(fio, iparams->vsite.a);
1951 gmx_fio_do_real(fio, iparams->vsite.b);
1956 gmx_fio_do_real(fio, iparams->vsite.a);
1957 gmx_fio_do_real(fio, iparams->vsite.b);
1958 gmx_fio_do_real(fio, iparams->vsite.c);
1961 gmx_fio_do_int(fio, iparams->vsiten.n);
1962 gmx_fio_do_real(fio, iparams->vsiten.a);
1967 /* We got rid of some parameters in version 68 */
1968 if (bRead && file_version < 68)
1970 gmx_fio_do_real(fio, rdum);
1971 gmx_fio_do_real(fio, rdum);
1972 gmx_fio_do_real(fio, rdum);
1973 gmx_fio_do_real(fio, rdum);
1975 gmx_fio_do_real(fio, iparams->gb.sar);
1976 gmx_fio_do_real(fio, iparams->gb.st);
1977 gmx_fio_do_real(fio, iparams->gb.pi);
1978 gmx_fio_do_real(fio, iparams->gb.gbr);
1979 gmx_fio_do_real(fio, iparams->gb.bmlt);
1982 gmx_fio_do_int(fio, iparams->cmap.cmapA);
1983 gmx_fio_do_int(fio, iparams->cmap.cmapB);
1986 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1987 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1991 gmx_fio_unset_comment(fio);
1995 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
2002 gmx_fio_set_comment(fio, interaction_function[ftype].name);
2004 if (file_version < 44)
2006 for (i = 0; i < MAXNODES; i++)
2008 gmx_fio_do_int(fio, idum);
2011 gmx_fio_do_int(fio, ilist->nr);
2014 snew(ilist->iatoms, ilist->nr);
2016 gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
2019 gmx_fio_unset_comment(fio);
2023 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
2024 gmx_bool bRead, int file_version)
2029 gmx_fio_do_int(fio, ffparams->atnr);
2030 if (file_version < 57)
2032 gmx_fio_do_int(fio, idum);
2034 gmx_fio_do_int(fio, ffparams->ntypes);
2037 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
2038 ffparams->atnr, ffparams->ntypes);
2042 snew(ffparams->functype, ffparams->ntypes);
2043 snew(ffparams->iparams, ffparams->ntypes);
2045 /* Read/write all the function types */
2046 gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
2049 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
2052 if (file_version >= 66)
2054 gmx_fio_do_double(fio, ffparams->reppow);
2058 ffparams->reppow = 12.0;
2061 if (file_version >= 57)
2063 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2066 /* Check whether all these function types are supported by the code.
2067 * In practice the code is backwards compatible, which means that the
2068 * numbering may have to be altered from old numbering to new numbering
2070 for (i = 0; (i < ffparams->ntypes); i++)
2074 /* Loop over file versions */
2075 for (k = 0; (k < NFTUPD); k++)
2077 /* Compare the read file_version to the update table */
2078 if ((file_version < ftupd[k].fvnr) &&
2079 (ffparams->functype[i] >= ftupd[k].ftype))
2081 ffparams->functype[i] += 1;
2084 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
2085 i, ffparams->functype[i],
2086 interaction_function[ftupd[k].ftype].longname);
2093 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2097 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2102 static void add_settle_atoms(t_ilist *ilist)
2106 /* Settle used to only store the first atom: add the other two */
2107 srenew(ilist->iatoms, 2*ilist->nr);
2108 for (i = ilist->nr/2-1; i >= 0; i--)
2110 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2111 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2112 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2113 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2115 ilist->nr = 2*ilist->nr;
2118 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2121 int i, j, renum[F_NRE];
2125 for (j = 0; (j < F_NRE); j++)
2130 for (k = 0; k < NFTUPD; k++)
2132 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2141 ilist[j].iatoms = NULL;
2145 do_ilist(fio, &ilist[j], bRead, file_version, j);
2146 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2148 add_settle_atoms(&ilist[j]);
2152 if (bRead && gmx_debug_at)
2153 pr_ilist(debug,0,interaction_function[j].longname,
2154 functype,&ilist[j],TRUE);
2159 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2160 gmx_bool bRead, int file_version)
2162 do_ffparams(fio, ffparams, bRead, file_version);
2164 if (file_version >= 54)
2166 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2169 do_ilists(fio, molt->ilist, bRead, file_version);
2172 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2174 int i, idum, dum_nra, *dum_a;
2176 if (file_version < 44)
2178 for (i = 0; i < MAXNODES; i++)
2180 gmx_fio_do_int(fio, idum);
2183 gmx_fio_do_int(fio, block->nr);
2184 if (file_version < 51)
2186 gmx_fio_do_int(fio, dum_nra);
2190 if ((block->nalloc_index > 0) && (NULL != block->index))
2192 sfree(block->index);
2194 block->nalloc_index = block->nr+1;
2195 snew(block->index, block->nalloc_index);
2197 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2199 if (file_version < 51 && dum_nra > 0)
2201 snew(dum_a, dum_nra);
2202 gmx_fio_ndo_int(fio, dum_a, dum_nra);
2207 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2212 if (file_version < 44)
2214 for (i = 0; i < MAXNODES; i++)
2216 gmx_fio_do_int(fio, idum);
2219 gmx_fio_do_int(fio, block->nr);
2220 gmx_fio_do_int(fio, block->nra);
2223 block->nalloc_index = block->nr+1;
2224 snew(block->index, block->nalloc_index);
2225 block->nalloc_a = block->nra;
2226 snew(block->a, block->nalloc_a);
2228 gmx_fio_ndo_int(fio, block->index, block->nr+1);
2229 gmx_fio_ndo_int(fio, block->a, block->nra);
2232 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2233 int file_version, gmx_groups_t *groups, int atnr)
2237 gmx_fio_do_real(fio, atom->m);
2238 gmx_fio_do_real(fio, atom->q);
2239 gmx_fio_do_real(fio, atom->mB);
2240 gmx_fio_do_real(fio, atom->qB);
2241 gmx_fio_do_ushort(fio, atom->type);
2242 gmx_fio_do_ushort(fio, atom->typeB);
2243 gmx_fio_do_int(fio, atom->ptype);
2244 gmx_fio_do_int(fio, atom->resind);
2245 if (file_version >= 52)
2247 gmx_fio_do_int(fio, atom->atomnumber);
2251 atom->atomnumber = NOTSET;
2253 if (file_version < 23)
2257 else if (file_version < 39)
2266 if (file_version < 57)
2268 unsigned char uchar[egcNR];
2269 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2270 for (i = myngrp; (i < ngrp); i++)
2274 /* Copy the old data format to the groups struct */
2275 for (i = 0; i < ngrp; i++)
2277 groups->grpnr[i][atnr] = uchar[i];
2282 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2287 if (file_version < 23)
2291 else if (file_version < 39)
2300 for (j = 0; (j < ngrp); j++)
2304 gmx_fio_do_int(fio, grps[j].nr);
2307 snew(grps[j].nm_ind, grps[j].nr);
2309 gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2314 snew(grps[j].nm_ind, grps[j].nr);
2319 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2325 gmx_fio_do_int(fio, ls);
2326 *nm = get_symtab_handle(symtab, ls);
2330 ls = lookup_symtab(symtab, *nm);
2331 gmx_fio_do_int(fio, ls);
2335 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2340 for (j = 0; (j < nstr); j++)
2342 do_symstr(fio, &(nm[j]), bRead, symtab);
2346 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2347 t_symtab *symtab, int file_version)
2351 for (j = 0; (j < n); j++)
2353 do_symstr(fio, &(ri[j].name), bRead, symtab);
2354 if (file_version >= 63)
2356 gmx_fio_do_int(fio, ri[j].nr);
2357 gmx_fio_do_uchar(fio, ri[j].ic);
2367 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2369 gmx_groups_t *groups)
2373 gmx_fio_do_int(fio, atoms->nr);
2374 gmx_fio_do_int(fio, atoms->nres);
2375 if (file_version < 57)
2377 gmx_fio_do_int(fio, groups->ngrpname);
2378 for (i = 0; i < egcNR; i++)
2380 groups->ngrpnr[i] = atoms->nr;
2381 snew(groups->grpnr[i], groups->ngrpnr[i]);
2386 snew(atoms->atom, atoms->nr);
2387 snew(atoms->atomname, atoms->nr);
2388 snew(atoms->atomtype, atoms->nr);
2389 snew(atoms->atomtypeB, atoms->nr);
2390 snew(atoms->resinfo, atoms->nres);
2391 if (file_version < 57)
2393 snew(groups->grpname, groups->ngrpname);
2395 atoms->pdbinfo = NULL;
2397 for (i = 0; (i < atoms->nr); i++)
2399 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2401 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2402 if (bRead && (file_version <= 20))
2404 for (i = 0; i < atoms->nr; i++)
2406 atoms->atomtype[i] = put_symtab(symtab, "?");
2407 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2412 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2413 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2415 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2417 if (file_version < 57)
2419 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2421 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2425 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2426 gmx_bool bRead, t_symtab *symtab,
2431 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2432 gmx_fio_do_int(fio, groups->ngrpname);
2435 snew(groups->grpname, groups->ngrpname);
2437 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2438 for (g = 0; g < egcNR; g++)
2440 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2441 if (groups->ngrpnr[g] == 0)
2445 groups->grpnr[g] = NULL;
2452 snew(groups->grpnr[g], groups->ngrpnr[g]);
2454 gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2459 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2464 if (file_version > 25)
2466 gmx_fio_do_int(fio, atomtypes->nr);
2470 snew(atomtypes->radius, j);
2471 snew(atomtypes->vol, j);
2472 snew(atomtypes->surftens, j);
2473 snew(atomtypes->atomnumber, j);
2474 snew(atomtypes->gb_radius, j);
2475 snew(atomtypes->S_hct, j);
2477 gmx_fio_ndo_real(fio, atomtypes->radius, j);
2478 gmx_fio_ndo_real(fio, atomtypes->vol, j);
2479 gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2480 if (file_version >= 40)
2482 gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2484 if (file_version >= 60)
2486 gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2487 gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2492 /* File versions prior to 26 cannot do GBSA,
2493 * so they dont use this structure
2496 atomtypes->radius = NULL;
2497 atomtypes->vol = NULL;
2498 atomtypes->surftens = NULL;
2499 atomtypes->atomnumber = NULL;
2500 atomtypes->gb_radius = NULL;
2501 atomtypes->S_hct = NULL;
2505 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2511 gmx_fio_do_int(fio, symtab->nr);
2515 snew(symtab->symbuf, 1);
2516 symbuf = symtab->symbuf;
2517 symbuf->bufsize = nr;
2518 snew(symbuf->buf, nr);
2519 for (i = 0; (i < nr); i++)
2521 gmx_fio_do_string(fio, buf);
2522 symbuf->buf[i] = strdup(buf);
2527 symbuf = symtab->symbuf;
2528 while (symbuf != NULL)
2530 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2532 gmx_fio_do_string(fio, symbuf->buf[i]);
2535 symbuf = symbuf->next;
2539 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2544 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2546 int i, j, ngrid, gs, nelem;
2548 gmx_fio_do_int(fio, cmap_grid->ngrid);
2549 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2551 ngrid = cmap_grid->ngrid;
2552 gs = cmap_grid->grid_spacing;
2557 snew(cmap_grid->cmapdata, ngrid);
2559 for (i = 0; i < cmap_grid->ngrid; i++)
2561 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2565 for (i = 0; i < cmap_grid->ngrid; i++)
2567 for (j = 0; j < nelem; j++)
2569 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2570 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2571 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2572 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2578 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2580 int m, a, a0, a1, r;
2584 /* We always assign a new chain number, but save the chain id characters
2585 * for larger molecules.
2587 #define CHAIN_MIN_ATOMS 15
2591 for (m = 0; m < mols->nr; m++)
2593 a0 = mols->index[m];
2594 a1 = mols->index[m+1];
2595 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2604 for (a = a0; a < a1; a++)
2606 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2607 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2612 /* Blank out the chain id if there was only one chain */
2615 for (r = 0; r < atoms->nres; r++)
2617 atoms->resinfo[r].chainid = ' ';
2622 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2623 t_symtab *symtab, int file_version,
2624 gmx_groups_t *groups)
2628 if (file_version >= 57)
2630 do_symstr(fio, &(molt->name), bRead, symtab);
2633 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2635 if (bRead && gmx_debug_at)
2637 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2640 if (file_version >= 57)
2642 do_ilists(fio, molt->ilist, bRead, file_version);
2644 do_block(fio, &molt->cgs, bRead, file_version);
2645 if (bRead && gmx_debug_at)
2647 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2651 /* This used to be in the atoms struct */
2652 do_blocka(fio, &molt->excls, bRead, file_version);
2655 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead)
2659 gmx_fio_do_int(fio, molb->type);
2660 gmx_fio_do_int(fio, molb->nmol);
2661 gmx_fio_do_int(fio, molb->natoms_mol);
2662 /* Position restraint coordinates */
2663 gmx_fio_do_int(fio, molb->nposres_xA);
2664 if (molb->nposres_xA > 0)
2668 snew(molb->posres_xA, molb->nposres_xA);
2670 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2672 gmx_fio_do_int(fio, molb->nposres_xB);
2673 if (molb->nposres_xB > 0)
2677 snew(molb->posres_xB, molb->nposres_xB);
2679 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2684 static t_block mtop_mols(gmx_mtop_t *mtop)
2690 for (mb = 0; mb < mtop->nmolblock; mb++)
2692 mols.nr += mtop->molblock[mb].nmol;
2694 mols.nalloc_index = mols.nr + 1;
2695 snew(mols.index, mols.nalloc_index);
2700 for (mb = 0; mb < mtop->nmolblock; mb++)
2702 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2704 a += mtop->molblock[mb].natoms_mol;
2713 static void add_posres_molblock(gmx_mtop_t *mtop)
2718 gmx_molblock_t *molb;
2721 /* posres reference positions are stored in ip->posres (if present) and
2722 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2723 posres.pos0A are identical to fbposres.pos0. */
2724 il = &mtop->moltype[0].ilist[F_POSRES];
2725 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2726 if (il->nr == 0 && ilfb->nr == 0)
2732 for (i = 0; i < il->nr; i += 2)
2734 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2735 am = max(am, il->iatoms[i+1]);
2736 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2737 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2738 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2743 /* This loop is required if we have only flat-bottomed posres:
2745 - bFE == FALSE (no B-state for flat-bottomed posres) */
2748 for (i = 0; i < ilfb->nr; i += 2)
2750 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2751 am = max(am, ilfb->iatoms[i+1]);
2754 /* Make the posres coordinate block end at a molecule end */
2756 while (am >= mtop->mols.index[mol+1])
2760 molb = &mtop->molblock[0];
2761 molb->nposres_xA = mtop->mols.index[mol+1];
2762 snew(molb->posres_xA, molb->nposres_xA);
2765 molb->nposres_xB = molb->nposres_xA;
2766 snew(molb->posres_xB, molb->nposres_xB);
2770 molb->nposres_xB = 0;
2772 for (i = 0; i < il->nr; i += 2)
2774 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2775 a = il->iatoms[i+1];
2776 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2777 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2778 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2781 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2782 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2783 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2788 /* If only flat-bottomed posres are present, take reference pos from them.
2789 Here: bFE == FALSE */
2790 for (i = 0; i < ilfb->nr; i += 2)
2792 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2793 a = ilfb->iatoms[i+1];
2794 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2795 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2796 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2801 static void set_disres_npair(gmx_mtop_t *mtop)
2808 ip = mtop->ffparams.iparams;
2810 for (mt = 0; mt < mtop->nmoltype; mt++)
2812 il = &mtop->moltype[mt].ilist[F_DISRES];
2817 for (i = 0; i < il->nr; i += 3)
2820 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2822 ip[a[i]].disres.npair = npair;
2830 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2840 do_symtab(fio, &(mtop->symtab), bRead);
2843 pr_symtab(debug, 0, "symtab", &mtop->symtab);
2846 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2848 if (file_version >= 57)
2850 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
2852 gmx_fio_do_int(fio, mtop->nmoltype);
2860 snew(mtop->moltype, mtop->nmoltype);
2861 if (file_version < 57)
2863 mtop->moltype[0].name = mtop->name;
2866 for (mt = 0; mt < mtop->nmoltype; mt++)
2868 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
2872 if (file_version >= 57)
2874 gmx_fio_do_int(fio, mtop->nmolblock);
2878 mtop->nmolblock = 1;
2882 snew(mtop->molblock, mtop->nmolblock);
2884 if (file_version >= 57)
2886 for (mb = 0; mb < mtop->nmolblock; mb++)
2888 do_molblock(fio, &mtop->molblock[mb], bRead);
2890 gmx_fio_do_int(fio, mtop->natoms);
2894 mtop->molblock[0].type = 0;
2895 mtop->molblock[0].nmol = 1;
2896 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
2897 mtop->molblock[0].nposres_xA = 0;
2898 mtop->molblock[0].nposres_xB = 0;
2901 do_atomtypes (fio, &(mtop->atomtypes), bRead, file_version);
2904 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
2907 if (file_version < 57)
2909 /* Debug statements are inside do_idef */
2910 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
2911 mtop->natoms = mtop->moltype[0].atoms.nr;
2914 if (file_version >= 65)
2916 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
2920 mtop->ffparams.cmap_grid.ngrid = 0;
2921 mtop->ffparams.cmap_grid.grid_spacing = 0;
2922 mtop->ffparams.cmap_grid.cmapdata = NULL;
2925 if (file_version >= 57)
2927 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
2930 if (file_version < 57)
2932 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
2933 if (bRead && gmx_debug_at)
2935 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
2937 do_block(fio, &mtop->mols, bRead, file_version);
2938 /* Add the posres coordinates to the molblock */
2939 add_posres_molblock(mtop);
2943 if (file_version >= 57)
2945 done_block(&mtop->mols);
2946 mtop->mols = mtop_mols(mtop);
2950 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
2954 if (file_version < 51)
2956 /* Here used to be the shake blocks */
2957 do_blocka(fio, &dumb, bRead, file_version);
2970 close_symtab(&(mtop->symtab));
2974 /* If TopOnlyOK is TRUE then we can read even future versions
2975 * of tpx files, provided the file_generation hasn't changed.
2976 * If it is FALSE, we need the inputrecord too, and bail out
2977 * if the file is newer than the program.
2979 * The version and generation if the topology (see top of this file)
2980 * are returned in the two last arguments.
2982 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2984 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
2985 gmx_bool TopOnlyOK, int *file_version,
2986 int *file_generation)
2989 char file_tag[STRLEN];
2996 gmx_fio_checktype(fio);
2997 gmx_fio_setdebug(fio, bDebugMode());
2999 /* NEW! XDR tpb file */
3000 precision = sizeof(real);
3003 gmx_fio_do_string(fio, buf);
3004 if (strncmp(buf, "VERSION", 7))
3006 gmx_fatal(FARGS, "Can not read file %s,\n"
3007 " this file is from a Gromacs version which is older than 2.0\n"
3008 " Make a new one with grompp or use a gro or pdb file, if possible",
3009 gmx_fio_getname(fio));
3011 gmx_fio_do_int(fio, precision);
3012 bDouble = (precision == sizeof(double));
3013 if ((precision != sizeof(float)) && !bDouble)
3015 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
3016 "instead of %d or %d",
3017 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
3019 gmx_fio_setprecision(fio, bDouble);
3020 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
3021 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
3025 gmx_fio_write_string(fio, GromacsVersion());
3026 bDouble = (precision == sizeof(double));
3027 gmx_fio_setprecision(fio, bDouble);
3028 gmx_fio_do_int(fio, precision);
3030 sprintf(file_tag, "%s", tpx_tag);
3031 fgen = tpx_generation;
3034 /* Check versions! */
3035 gmx_fio_do_int(fio, fver);
3037 /* This is for backward compatibility with development versions 77-79
3038 * where the tag was, mistakenly, placed before the generation,
3039 * which would cause a segv instead of a proper error message
3040 * when reading the topology only from tpx with <77 code.
3042 if (fver >= 77 && fver <= 79)
3044 gmx_fio_do_string(fio, file_tag);
3049 gmx_fio_do_int(fio, fgen);
3058 gmx_fio_do_string(fio, file_tag);
3064 /* Versions before 77 don't have the tag, set it to release */
3065 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
3068 if (strcmp(file_tag, tpx_tag) != 0)
3070 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
3073 /* We only support reading tpx files with the same tag as the code
3074 * or tpx files with the release tag and with lower version number.
3076 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
3078 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
3079 gmx_fio_getname(fio), fver, file_tag,
3080 tpx_version, tpx_tag);
3085 if (file_version != NULL)
3087 *file_version = fver;
3089 if (file_generation != NULL)
3091 *file_generation = fgen;
3095 if ((fver <= tpx_incompatible_version) ||
3096 ((fver > tpx_version) && !TopOnlyOK) ||
3097 (fgen > tpx_generation) ||
3098 tpx_version == 80) /*80 was used by both 5.0-dev and 4.6-dev*/
3100 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3101 gmx_fio_getname(fio), fver, tpx_version);
3104 do_section(fio, eitemHEADER, bRead);
3105 gmx_fio_do_int(fio, tpx->natoms);
3108 gmx_fio_do_int(fio, tpx->ngtc);
3116 gmx_fio_do_int(fio, idum);
3117 gmx_fio_do_real(fio, rdum);
3119 /*a better decision will eventually (5.0 or later) need to be made
3120 on how to treat the alchemical state of the system, which can now
3121 vary through a simulation, and cannot be completely described
3122 though a single lambda variable, or even a single state
3123 index. Eventually, should probably be a vector. MRS*/
3126 gmx_fio_do_int(fio, tpx->fep_state);
3128 gmx_fio_do_real(fio, tpx->lambda);
3129 gmx_fio_do_int(fio, tpx->bIr);
3130 gmx_fio_do_int(fio, tpx->bTop);
3131 gmx_fio_do_int(fio, tpx->bX);
3132 gmx_fio_do_int(fio, tpx->bV);
3133 gmx_fio_do_int(fio, tpx->bF);
3134 gmx_fio_do_int(fio, tpx->bBox);
3136 if ((fgen > tpx_generation))
3138 /* This can only happen if TopOnlyOK=TRUE */
3143 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3144 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3145 gmx_bool bXVallocated)
3151 int file_version, file_generation;
3155 gmx_bool bPeriodicMols;
3159 tpx.natoms = state->natoms;
3160 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3161 tpx.fep_state = state->fep_state;
3162 tpx.lambda = state->lambda[efptFEP];
3163 tpx.bIr = (ir != NULL);
3164 tpx.bTop = (mtop != NULL);
3165 tpx.bX = (state->x != NULL);
3166 tpx.bV = (state->v != NULL);
3167 tpx.bF = (f != NULL);
3171 TopOnlyOK = (ir == NULL);
3173 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3178 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3179 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3184 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3185 state->natoms = tpx.natoms;
3186 state->nalloc = tpx.natoms;
3192 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3196 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3198 do_test(fio, tpx.bBox, state->box);
3199 do_section(fio, eitemBOX, bRead);
3202 gmx_fio_ndo_rvec(fio, state->box, DIM);
3203 if (file_version >= 51)
3205 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3209 /* We initialize box_rel after reading the inputrec */
3210 clear_mat(state->box_rel);
3212 if (file_version >= 28)
3214 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3215 if (file_version < 56)
3218 gmx_fio_ndo_rvec(fio, mdum, DIM);
3223 if (state->ngtc > 0 && file_version >= 28)
3226 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3227 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3228 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3229 snew(dumv, state->ngtc);
3230 if (file_version < 69)
3232 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3234 /* These used to be the Berendsen tcoupl_lambda's */
3235 gmx_fio_ndo_real(fio, dumv, state->ngtc);
3239 /* Prior to tpx version 26, the inputrec was here.
3240 * I moved it to enable partial forward-compatibility
3241 * for analysis/viewer programs.
3243 if (file_version < 26)
3245 do_test(fio, tpx.bIr, ir);
3246 do_section(fio, eitemIR, bRead);
3251 do_inputrec(fio, ir, bRead, file_version,
3252 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3255 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3260 do_inputrec(fio, &dum_ir, bRead, file_version,
3261 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3264 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3266 done_inputrec(&dum_ir);
3272 do_test(fio, tpx.bTop, mtop);
3273 do_section(fio, eitemTOP, bRead);
3278 do_mtop(fio, mtop, bRead, file_version);
3282 do_mtop(fio, &dum_top, bRead, file_version);
3283 done_mtop(&dum_top, TRUE);
3286 do_test(fio, tpx.bX, state->x);
3287 do_section(fio, eitemX, bRead);
3292 state->flags |= (1<<estX);
3294 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3297 do_test(fio, tpx.bV, state->v);
3298 do_section(fio, eitemV, bRead);
3303 state->flags |= (1<<estV);
3305 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3308 do_test(fio, tpx.bF, f);
3309 do_section(fio, eitemF, bRead);
3312 gmx_fio_ndo_rvec(fio, f, state->natoms);
3315 /* Starting with tpx version 26, we have the inputrec
3316 * at the end of the file, so we can ignore it
3317 * if the file is never than the software (but still the
3318 * same generation - see comments at the top of this file.
3323 bPeriodicMols = FALSE;
3324 if (file_version >= 26)
3326 do_test(fio, tpx.bIr, ir);
3327 do_section(fio, eitemIR, bRead);
3330 if (file_version >= 53)
3332 /* Removed the pbc info from do_inputrec, since we always want it */
3336 bPeriodicMols = ir->bPeriodicMols;
3338 gmx_fio_do_int(fio, ePBC);
3339 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3341 if (file_generation <= tpx_generation && ir)
3343 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3346 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3348 if (file_version < 51)
3350 set_box_rel(ir, state);
3352 if (file_version < 53)
3355 bPeriodicMols = ir->bPeriodicMols;
3358 if (bRead && ir && file_version >= 53)
3360 /* We need to do this after do_inputrec, since that initializes ir */
3362 ir->bPeriodicMols = bPeriodicMols;
3371 if (state->ngtc == 0)
3373 /* Reading old version without tcoupl state data: set it */
3374 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3376 if (tpx.bTop && mtop)
3378 if (file_version < 57)
3380 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3382 ir->eDisre = edrSimple;
3386 ir->eDisre = edrNone;
3389 set_disres_npair(mtop);
3393 if (tpx.bTop && mtop)
3395 gmx_mtop_finalize(mtop);
3398 if (file_version >= 57)
3402 env = getenv("GMX_NOCHARGEGROUPS");
3405 sscanf(env, "%d", &ienv);
3406 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3411 "Will make single atomic charge groups in non-solvent%s\n",
3412 ienv > 1 ? " and solvent" : "");
3413 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3415 fprintf(stderr, "\n");
3423 /************************************************************
3425 * The following routines are the exported ones
3427 ************************************************************/
3429 t_fileio *open_tpx(const char *fn, const char *mode)
3431 return gmx_fio_open(fn, mode);
3434 void close_tpx(t_fileio *fio)
3439 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3440 int *file_version, int *file_generation)
3444 fio = open_tpx(fn, "r");
3445 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3449 void write_tpx_state(const char *fn,
3450 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3454 fio = open_tpx(fn, "w");
3455 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3459 void read_tpx_state(const char *fn,
3460 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3464 fio = open_tpx(fn, "r");
3465 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3469 int read_tpx(const char *fn,
3470 t_inputrec *ir, matrix box, int *natoms,
3471 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3479 fio = open_tpx(fn, "r");
3480 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3482 *natoms = state.natoms;
3485 copy_mat(state.box, box);
3494 int read_tpx_top(const char *fn,
3495 t_inputrec *ir, matrix box, int *natoms,
3496 rvec *x, rvec *v, rvec *f, t_topology *top)
3502 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3504 *top = gmx_mtop_t_to_t_topology(&mtop);
3509 gmx_bool fn2bTPX(const char *file)
3511 switch (fn2ftp(file))
3522 static void done_gmx_groups_t(gmx_groups_t *g)
3526 for (i = 0; (i < egcNR); i++)
3528 if (NULL != g->grps[i].nm_ind)
3530 sfree(g->grps[i].nm_ind);
3531 g->grps[i].nm_ind = NULL;
3533 if (NULL != g->grpnr[i])
3539 /* The contents of this array is in symtab, don't free it here */
3543 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3544 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3547 int natoms, i, version, generation;
3548 gmx_bool bTop, bXNULL = FALSE;
3550 t_topology *topconv;
3553 bTop = fn2bTPX(infile);
3557 read_tpxheader(infile, &header, TRUE, &version, &generation);
3560 snew(*x, header.natoms);
3564 snew(*v, header.natoms);
3567 *ePBC = read_tpx(infile, NULL, box, &natoms,
3568 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3569 *top = gmx_mtop_t_to_t_topology(mtop);
3570 /* In this case we need to throw away the group data too */
3571 done_gmx_groups_t(&mtop->groups);
3573 strcpy(title, *top->name);
3574 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3578 get_stx_coordnum(infile, &natoms);
3579 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3590 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3598 aps = gmx_atomprop_init();
3599 for (i = 0; (i < natoms); i++)
3601 if (!gmx_atomprop_query(aps, epropMass,
3602 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3603 *top->atoms.atomname[i],
3604 &(top->atoms.atom[i].m)))
3608 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3609 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3610 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3611 *top->atoms.atomname[i]);
3615 gmx_atomprop_destroy(aps);
3617 top->idef.ntypes = -1;