1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROningen Mixture of Alchemy and Childrens' Stories
40 /* This file is completely threadsafe - keep it that way! */
42 #include <thread_mpi.h>
50 #include "gmx_fatal.h"
64 #include "mtop_util.h"
66 #define TPX_TAG_RELEASE "release"
68 /* This is the tag string which is stored in the tpx file.
69 * Change this if you want to change the tpx format in a feature branch.
70 * This ensures that there will not be different tpx formats around which
71 * can not be distinguished.
73 static const char *tpx_tag = TPX_TAG_RELEASE;
75 /* This number should be increased whenever the file format changes! */
76 static const int tpx_version = 92;
78 /* This number should only be increased when you edit the TOPOLOGY section
79 * or the HEADER of the tpx format.
80 * This way we can maintain forward compatibility too for all analysis tools
81 * and/or external programs that only need to know the atom/residue names,
82 * charges, and bond connectivity.
84 * It first appeared in tpx version 26, when I also moved the inputrecord
85 * to the end of the tpx file, so we can just skip it if we only
88 static const int tpx_generation = 25;
90 /* This number should be the most recent backwards incompatible version
91 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
93 static const int tpx_incompatible_version = 9;
97 /* Struct used to maintain tpx compatibility when function types are added */
99 int fvnr; /* file version number in which the function type first appeared */
100 int ftype; /* function type */
104 * The entries should be ordered in:
105 * 1. ascending file version number
106 * 2. ascending function type number
108 /*static const t_ftupd ftupd[] = {
109 { 20, F_CUBICBONDS },
113 { 22, F_DISRESVIOL },
119 { 26, F_DIHRESVIOL },
120 { 30, F_CROSS_BOND_BONDS },
121 { 30, F_CROSS_BOND_ANGLES },
122 { 30, F_UREY_BRADLEY },
123 { 30, F_POLARIZATION },
127 * The entries should be ordered in:
128 * 1. ascending function type number
129 * 2. ascending file version number
131 /* question; what is the purpose of the commented code above? */
132 static const t_ftupd ftupd[] = {
133 { 20, F_CUBICBONDS },
138 { 43, F_TABBONDSNC },
139 { 70, F_RESTRBONDS },
140 { 76, F_LINEAR_ANGLES },
141 { 30, F_CROSS_BOND_BONDS },
142 { 30, F_CROSS_BOND_ANGLES },
143 { 30, F_UREY_BRADLEY },
144 { 34, F_QUARTIC_ANGLES },
154 { 72, F_NPSOLVATION },
156 { 41, F_LJC_PAIRS_NB },
159 { 32, F_COUL_RECIP },
161 { 30, F_POLARIZATION },
164 { 22, F_DISRESVIOL },
168 { 26, F_DIHRESVIOL },
173 { 46, F_ECONSERVED },
174 { 69, F_VTEMP_NOLONGERUSED},
176 { 54, F_DVDL_CONSTR },
177 { 76, F_ANHARM_POL },
180 { 79, F_DVDL_BONDED, },
181 { 79, F_DVDL_RESTRAINT },
182 { 79, F_DVDL_TEMPERATURE },
184 #define NFTUPD asize(ftupd)
186 /* Needed for backward compatibility */
189 static void _do_section(t_fileio *fio, int key, gmx_bool bRead, const char *src,
195 if (gmx_fio_getftp(fio) == efTPA)
199 gmx_fio_write_string(fio, itemstr[key]);
200 bDbg = gmx_fio_getdebug(fio);
201 gmx_fio_setdebug(fio, FALSE);
202 gmx_fio_write_string(fio, comment_str[key]);
203 gmx_fio_setdebug(fio, bDbg);
207 if (gmx_fio_getdebug(fio))
209 fprintf(stderr, "Looking for section %s (%s, %d)",
210 itemstr[key], src, line);
215 gmx_fio_do_string(fio, buf);
217 while ((gmx_strcasecmp(buf, itemstr[key]) != 0));
219 if (gmx_strcasecmp(buf, itemstr[key]) != 0)
221 gmx_fatal(FARGS, "\nCould not find section heading %s", itemstr[key]);
223 else if (gmx_fio_getdebug(fio))
225 fprintf(stderr, " and found it\n");
231 #define do_section(fio, key, bRead) _do_section(fio, key, bRead, __FILE__, __LINE__)
233 /**************************************************************
235 * Now the higer level routines that do io of the structures and arrays
237 **************************************************************/
238 static void do_pullgrp(t_fileio *fio, t_pullgrp *pgrp, gmx_bool bRead,
241 gmx_bool bDum = TRUE;
244 gmx_fio_do_int(fio, pgrp->nat);
247 snew(pgrp->ind, pgrp->nat);
249 bDum = gmx_fio_ndo_int(fio, pgrp->ind, pgrp->nat);
250 gmx_fio_do_int(fio, pgrp->nweight);
253 snew(pgrp->weight, pgrp->nweight);
255 bDum = gmx_fio_ndo_real(fio, pgrp->weight, pgrp->nweight);
256 gmx_fio_do_int(fio, pgrp->pbcatom);
257 gmx_fio_do_rvec(fio, pgrp->vec);
258 gmx_fio_do_rvec(fio, pgrp->init);
259 gmx_fio_do_real(fio, pgrp->rate);
260 gmx_fio_do_real(fio, pgrp->k);
261 if (file_version >= 56)
263 gmx_fio_do_real(fio, pgrp->kB);
271 static void do_expandedvals(t_fileio *fio, t_expanded *expand, t_lambda *fepvals, gmx_bool bRead, int file_version)
273 /* i is used in the ndo_double macro*/
276 gmx_bool bDum = TRUE;
278 int n_lambda = fepvals->n_lambda;
280 /* reset the lambda calculation window */
281 fepvals->lambda_start_n = 0;
282 fepvals->lambda_stop_n = n_lambda;
283 if (file_version >= 79)
289 snew(expand->init_lambda_weights, n_lambda);
291 bDum = gmx_fio_ndo_real(fio, expand->init_lambda_weights, n_lambda);
292 gmx_fio_do_gmx_bool(fio, expand->bInit_weights);
295 gmx_fio_do_int(fio, expand->nstexpanded);
296 gmx_fio_do_int(fio, expand->elmcmove);
297 gmx_fio_do_int(fio, expand->elamstats);
298 gmx_fio_do_int(fio, expand->lmc_repeats);
299 gmx_fio_do_int(fio, expand->gibbsdeltalam);
300 gmx_fio_do_int(fio, expand->lmc_forced_nstart);
301 gmx_fio_do_int(fio, expand->lmc_seed);
302 gmx_fio_do_real(fio, expand->mc_temp);
303 gmx_fio_do_int(fio, expand->bSymmetrizedTMatrix);
304 gmx_fio_do_int(fio, expand->nstTij);
305 gmx_fio_do_int(fio, expand->minvarmin);
306 gmx_fio_do_int(fio, expand->c_range);
307 gmx_fio_do_real(fio, expand->wl_scale);
308 gmx_fio_do_real(fio, expand->wl_ratio);
309 gmx_fio_do_real(fio, expand->init_wl_delta);
310 gmx_fio_do_gmx_bool(fio, expand->bWLoneovert);
311 gmx_fio_do_int(fio, expand->elmceq);
312 gmx_fio_do_int(fio, expand->equil_steps);
313 gmx_fio_do_int(fio, expand->equil_samples);
314 gmx_fio_do_int(fio, expand->equil_n_at_lam);
315 gmx_fio_do_real(fio, expand->equil_wl_delta);
316 gmx_fio_do_real(fio, expand->equil_ratio);
320 static void do_simtempvals(t_fileio *fio, t_simtemp *simtemp, int n_lambda, gmx_bool bRead,
323 gmx_bool bDum = TRUE;
325 if (file_version >= 79)
327 gmx_fio_do_int(fio, simtemp->eSimTempScale);
328 gmx_fio_do_real(fio, simtemp->simtemp_high);
329 gmx_fio_do_real(fio, simtemp->simtemp_low);
334 snew(simtemp->temperatures, n_lambda);
336 bDum = gmx_fio_ndo_real(fio, simtemp->temperatures, n_lambda);
341 static void do_fepvals(t_fileio *fio, t_lambda *fepvals, gmx_bool bRead, int file_version)
343 /* i is defined in the ndo_double macro; use g to iterate. */
346 gmx_bool bDum = TRUE;
349 /* free energy values */
351 if (file_version >= 79)
353 gmx_fio_do_int(fio, fepvals->init_fep_state);
354 gmx_fio_do_double(fio, fepvals->init_lambda);
355 gmx_fio_do_double(fio, fepvals->delta_lambda);
357 else if (file_version >= 59)
359 gmx_fio_do_double(fio, fepvals->init_lambda);
360 gmx_fio_do_double(fio, fepvals->delta_lambda);
364 gmx_fio_do_real(fio, rdum);
365 fepvals->init_lambda = rdum;
366 gmx_fio_do_real(fio, rdum);
367 fepvals->delta_lambda = rdum;
369 if (file_version >= 79)
371 gmx_fio_do_int(fio, fepvals->n_lambda);
374 snew(fepvals->all_lambda, efptNR);
376 for (g = 0; g < efptNR; g++)
378 if (fepvals->n_lambda > 0)
382 snew(fepvals->all_lambda[g], fepvals->n_lambda);
384 bDum = gmx_fio_ndo_double(fio, fepvals->all_lambda[g], fepvals->n_lambda);
385 bDum = gmx_fio_ndo_int(fio, fepvals->separate_dvdl, efptNR);
387 else if (fepvals->init_lambda >= 0)
389 fepvals->separate_dvdl[efptFEP] = TRUE;
393 else if (file_version >= 64)
395 gmx_fio_do_int(fio, fepvals->n_lambda);
400 snew(fepvals->all_lambda, efptNR);
401 /* still allocate the all_lambda array's contents. */
402 for (g = 0; g < efptNR; g++)
404 if (fepvals->n_lambda > 0)
406 snew(fepvals->all_lambda[g], fepvals->n_lambda);
410 bDum = gmx_fio_ndo_double(fio, fepvals->all_lambda[efptFEP],
412 if (fepvals->init_lambda >= 0)
416 fepvals->separate_dvdl[efptFEP] = TRUE;
420 /* copy the contents of the efptFEP lambda component to all
421 the other components */
422 for (g = 0; g < efptNR; g++)
424 for (h = 0; h < fepvals->n_lambda; h++)
428 fepvals->all_lambda[g][h] =
429 fepvals->all_lambda[efptFEP][h];
438 fepvals->n_lambda = 0;
439 fepvals->all_lambda = NULL;
440 if (fepvals->init_lambda >= 0)
442 fepvals->separate_dvdl[efptFEP] = TRUE;
445 if (file_version >= 13)
447 gmx_fio_do_real(fio, fepvals->sc_alpha);
451 fepvals->sc_alpha = 0;
453 if (file_version >= 38)
455 gmx_fio_do_int(fio, fepvals->sc_power);
459 fepvals->sc_power = 2;
461 if (file_version >= 79)
463 gmx_fio_do_real(fio, fepvals->sc_r_power);
467 fepvals->sc_r_power = 6.0;
469 if (file_version >= 15)
471 gmx_fio_do_real(fio, fepvals->sc_sigma);
475 fepvals->sc_sigma = 0.3;
479 if (file_version >= 71)
481 fepvals->sc_sigma_min = fepvals->sc_sigma;
485 fepvals->sc_sigma_min = 0;
488 if (file_version >= 79)
490 gmx_fio_do_int(fio, fepvals->bScCoul);
494 fepvals->bScCoul = TRUE;
496 if (file_version >= 64)
498 gmx_fio_do_int(fio, fepvals->nstdhdl);
502 fepvals->nstdhdl = 1;
505 if (file_version >= 73)
507 gmx_fio_do_int(fio, fepvals->separate_dhdl_file);
508 gmx_fio_do_int(fio, fepvals->dhdl_derivatives);
512 fepvals->separate_dhdl_file = esepdhdlfileYES;
513 fepvals->dhdl_derivatives = edhdlderivativesYES;
515 if (file_version >= 71)
517 gmx_fio_do_int(fio, fepvals->dh_hist_size);
518 gmx_fio_do_double(fio, fepvals->dh_hist_spacing);
522 fepvals->dh_hist_size = 0;
523 fepvals->dh_hist_spacing = 0.1;
525 if (file_version >= 79)
527 gmx_fio_do_int(fio, fepvals->bPrintEnergy);
531 fepvals->bPrintEnergy = FALSE;
534 /* handle lambda_neighbors */
535 if ((file_version >= 83 && file_version < 90) || file_version >= 92)
537 gmx_fio_do_int(fio, fepvals->lambda_neighbors);
538 if ( (fepvals->lambda_neighbors >= 0) && (fepvals->init_fep_state >= 0) &&
539 (fepvals->init_lambda < 0) )
541 fepvals->lambda_start_n = (fepvals->init_fep_state -
542 fepvals->lambda_neighbors);
543 fepvals->lambda_stop_n = (fepvals->init_fep_state +
544 fepvals->lambda_neighbors + 1);
545 if (fepvals->lambda_start_n < 0)
547 fepvals->lambda_start_n = 0;;
549 if (fepvals->lambda_stop_n >= fepvals->n_lambda)
551 fepvals->lambda_stop_n = fepvals->n_lambda;
556 fepvals->lambda_start_n = 0;
557 fepvals->lambda_stop_n = fepvals->n_lambda;
562 fepvals->lambda_start_n = 0;
563 fepvals->lambda_stop_n = fepvals->n_lambda;
567 static void do_pull(t_fileio *fio, t_pull *pull, gmx_bool bRead, int file_version)
571 gmx_fio_do_int(fio, pull->ngrp);
572 gmx_fio_do_int(fio, pull->eGeom);
573 gmx_fio_do_ivec(fio, pull->dim);
574 gmx_fio_do_real(fio, pull->cyl_r1);
575 gmx_fio_do_real(fio, pull->cyl_r0);
576 gmx_fio_do_real(fio, pull->constr_tol);
577 gmx_fio_do_int(fio, pull->nstxout);
578 gmx_fio_do_int(fio, pull->nstfout);
581 snew(pull->grp, pull->ngrp+1);
583 for (g = 0; g < pull->ngrp+1; g++)
585 do_pullgrp(fio, &pull->grp[g], bRead, file_version);
590 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg, gmx_bool bRead, int file_version)
592 gmx_bool bDum = TRUE;
595 gmx_fio_do_int(fio, rotg->eType);
596 gmx_fio_do_int(fio, rotg->bMassW);
597 gmx_fio_do_int(fio, rotg->nat);
600 snew(rotg->ind, rotg->nat);
602 gmx_fio_ndo_int(fio, rotg->ind, rotg->nat);
605 snew(rotg->x_ref, rotg->nat);
607 gmx_fio_ndo_rvec(fio, rotg->x_ref, rotg->nat);
608 gmx_fio_do_rvec(fio, rotg->vec);
609 gmx_fio_do_rvec(fio, rotg->pivot);
610 gmx_fio_do_real(fio, rotg->rate);
611 gmx_fio_do_real(fio, rotg->k);
612 gmx_fio_do_real(fio, rotg->slab_dist);
613 gmx_fio_do_real(fio, rotg->min_gaussian);
614 gmx_fio_do_real(fio, rotg->eps);
615 gmx_fio_do_int(fio, rotg->eFittype);
616 gmx_fio_do_int(fio, rotg->PotAngle_nstep);
617 gmx_fio_do_real(fio, rotg->PotAngle_step);
620 static void do_rot(t_fileio *fio, t_rot *rot, gmx_bool bRead, int file_version)
624 gmx_fio_do_int(fio, rot->ngrp);
625 gmx_fio_do_int(fio, rot->nstrout);
626 gmx_fio_do_int(fio, rot->nstsout);
629 snew(rot->grp, rot->ngrp);
631 for (g = 0; g < rot->ngrp; g++)
633 do_rotgrp(fio, &rot->grp[g], bRead, file_version);
638 static void do_inputrec(t_fileio *fio, t_inputrec *ir, gmx_bool bRead,
639 int file_version, real *fudgeQQ)
641 int i, j, k, *tmp, idum = 0;
642 gmx_bool bDum = TRUE;
646 real zerotemptime, finish_t, init_temp, finish_temp;
648 if (file_version != tpx_version)
650 /* Give a warning about features that are not accessible */
651 fprintf(stderr, "Note: file tpx version %d, software tpx version %d\n",
652 file_version, tpx_version);
660 if (file_version == 0)
665 /* Basic inputrec stuff */
666 gmx_fio_do_int(fio, ir->eI);
667 if (file_version >= 62)
669 gmx_fio_do_gmx_large_int(fio, ir->nsteps);
673 gmx_fio_do_int(fio, idum);
676 if (file_version > 25)
678 if (file_version >= 62)
680 gmx_fio_do_gmx_large_int(fio, ir->init_step);
684 gmx_fio_do_int(fio, idum);
685 ir->init_step = idum;
693 if (file_version >= 58)
695 gmx_fio_do_int(fio, ir->simulation_part);
699 ir->simulation_part = 1;
702 if (file_version >= 67)
704 gmx_fio_do_int(fio, ir->nstcalcenergy);
708 ir->nstcalcenergy = 1;
710 if (file_version < 53)
712 /* The pbc info has been moved out of do_inputrec,
713 * since we always want it, also without reading the inputrec.
715 gmx_fio_do_int(fio, ir->ePBC);
716 if ((file_version <= 15) && (ir->ePBC == 2))
720 if (file_version >= 45)
722 gmx_fio_do_int(fio, ir->bPeriodicMols);
729 ir->bPeriodicMols = TRUE;
733 ir->bPeriodicMols = FALSE;
737 if (file_version >= 81)
739 gmx_fio_do_int(fio, ir->cutoff_scheme);
743 ir->cutoff_scheme = ecutsGROUP;
745 gmx_fio_do_int(fio, ir->ns_type);
746 gmx_fio_do_int(fio, ir->nstlist);
747 gmx_fio_do_int(fio, ir->ndelta);
748 if (file_version < 41)
750 gmx_fio_do_int(fio, idum);
751 gmx_fio_do_int(fio, idum);
753 if (file_version >= 45)
755 gmx_fio_do_real(fio, ir->rtpi);
761 gmx_fio_do_int(fio, ir->nstcomm);
762 if (file_version > 34)
764 gmx_fio_do_int(fio, ir->comm_mode);
766 else if (ir->nstcomm < 0)
768 ir->comm_mode = ecmANGULAR;
772 ir->comm_mode = ecmLINEAR;
774 ir->nstcomm = abs(ir->nstcomm);
776 if (file_version > 25)
778 gmx_fio_do_int(fio, ir->nstcheckpoint);
782 ir->nstcheckpoint = 0;
785 gmx_fio_do_int(fio, ir->nstcgsteep);
787 if (file_version >= 30)
789 gmx_fio_do_int(fio, ir->nbfgscorr);
796 gmx_fio_do_int(fio, ir->nstlog);
797 gmx_fio_do_int(fio, ir->nstxout);
798 gmx_fio_do_int(fio, ir->nstvout);
799 gmx_fio_do_int(fio, ir->nstfout);
800 gmx_fio_do_int(fio, ir->nstenergy);
801 gmx_fio_do_int(fio, ir->nstxtcout);
802 if (file_version >= 59)
804 gmx_fio_do_double(fio, ir->init_t);
805 gmx_fio_do_double(fio, ir->delta_t);
809 gmx_fio_do_real(fio, rdum);
811 gmx_fio_do_real(fio, rdum);
814 gmx_fio_do_real(fio, ir->xtcprec);
815 if (file_version < 19)
817 gmx_fio_do_int(fio, idum);
818 gmx_fio_do_int(fio, idum);
820 if (file_version < 18)
822 gmx_fio_do_int(fio, idum);
824 if (file_version >= 81)
826 gmx_fio_do_real(fio, ir->verletbuf_drift);
830 ir->verletbuf_drift = 0;
832 gmx_fio_do_real(fio, ir->rlist);
833 if (file_version >= 67)
835 gmx_fio_do_real(fio, ir->rlistlong);
837 if (file_version >= 82 && file_version != 90)
839 gmx_fio_do_int(fio, ir->nstcalclr);
843 /* Calculate at NS steps */
844 ir->nstcalclr = ir->nstlist;
846 gmx_fio_do_int(fio, ir->coulombtype);
847 if (file_version < 32 && ir->coulombtype == eelRF)
849 ir->coulombtype = eelRF_NEC;
851 if (file_version >= 81)
853 gmx_fio_do_int(fio, ir->coulomb_modifier);
857 ir->coulomb_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
859 gmx_fio_do_real(fio, ir->rcoulomb_switch);
860 gmx_fio_do_real(fio, ir->rcoulomb);
861 gmx_fio_do_int(fio, ir->vdwtype);
862 if (file_version >= 81)
864 gmx_fio_do_int(fio, ir->vdw_modifier);
868 ir->vdw_modifier = (ir->cutoff_scheme == ecutsVERLET ? eintmodPOTSHIFT : eintmodNONE);
870 gmx_fio_do_real(fio, ir->rvdw_switch);
871 gmx_fio_do_real(fio, ir->rvdw);
872 if (file_version < 67)
874 ir->rlistlong = max_cutoff(ir->rlist, max_cutoff(ir->rvdw, ir->rcoulomb));
876 gmx_fio_do_int(fio, ir->eDispCorr);
877 gmx_fio_do_real(fio, ir->epsilon_r);
878 if (file_version >= 37)
880 gmx_fio_do_real(fio, ir->epsilon_rf);
884 if (EEL_RF(ir->coulombtype))
886 ir->epsilon_rf = ir->epsilon_r;
891 ir->epsilon_rf = 1.0;
894 if (file_version >= 29)
896 gmx_fio_do_real(fio, ir->tabext);
903 if (file_version > 25)
905 gmx_fio_do_int(fio, ir->gb_algorithm);
906 gmx_fio_do_int(fio, ir->nstgbradii);
907 gmx_fio_do_real(fio, ir->rgbradii);
908 gmx_fio_do_real(fio, ir->gb_saltconc);
909 gmx_fio_do_int(fio, ir->implicit_solvent);
913 ir->gb_algorithm = egbSTILL;
917 ir->implicit_solvent = eisNO;
919 if (file_version >= 55)
921 gmx_fio_do_real(fio, ir->gb_epsilon_solvent);
922 gmx_fio_do_real(fio, ir->gb_obc_alpha);
923 gmx_fio_do_real(fio, ir->gb_obc_beta);
924 gmx_fio_do_real(fio, ir->gb_obc_gamma);
925 if (file_version >= 60)
927 gmx_fio_do_real(fio, ir->gb_dielectric_offset);
928 gmx_fio_do_int(fio, ir->sa_algorithm);
932 ir->gb_dielectric_offset = 0.009;
933 ir->sa_algorithm = esaAPPROX;
935 gmx_fio_do_real(fio, ir->sa_surface_tension);
937 /* Override sa_surface_tension if it is not changed in the mpd-file */
938 if (ir->sa_surface_tension < 0)
940 if (ir->gb_algorithm == egbSTILL)
942 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
944 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
946 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
953 /* Better use sensible values than insane (0.0) ones... */
954 ir->gb_epsilon_solvent = 80;
955 ir->gb_obc_alpha = 1.0;
956 ir->gb_obc_beta = 0.8;
957 ir->gb_obc_gamma = 4.85;
958 ir->sa_surface_tension = 2.092;
962 if (file_version >= 81)
964 gmx_fio_do_real(fio, ir->fourier_spacing);
968 ir->fourier_spacing = 0.0;
970 gmx_fio_do_int(fio, ir->nkx);
971 gmx_fio_do_int(fio, ir->nky);
972 gmx_fio_do_int(fio, ir->nkz);
973 gmx_fio_do_int(fio, ir->pme_order);
974 gmx_fio_do_real(fio, ir->ewald_rtol);
976 if (file_version >= 24)
978 gmx_fio_do_int(fio, ir->ewald_geometry);
982 ir->ewald_geometry = eewg3D;
985 if (file_version <= 17)
987 ir->epsilon_surface = 0;
988 if (file_version == 17)
990 gmx_fio_do_int(fio, idum);
995 gmx_fio_do_real(fio, ir->epsilon_surface);
998 gmx_fio_do_gmx_bool(fio, ir->bOptFFT);
1000 gmx_fio_do_gmx_bool(fio, ir->bContinuation);
1001 gmx_fio_do_int(fio, ir->etc);
1002 /* before version 18, ir->etc was a gmx_bool (ir->btc),
1003 * but the values 0 and 1 still mean no and
1004 * berendsen temperature coupling, respectively.
1006 if (file_version >= 79)
1008 gmx_fio_do_gmx_bool(fio, ir->bPrintNHChains);
1010 if (file_version >= 71)
1012 gmx_fio_do_int(fio, ir->nsttcouple);
1016 ir->nsttcouple = ir->nstcalcenergy;
1018 if (file_version <= 15)
1020 gmx_fio_do_int(fio, idum);
1022 if (file_version <= 17)
1024 gmx_fio_do_int(fio, ir->epct);
1025 if (file_version <= 15)
1029 ir->epct = epctSURFACETENSION;
1031 gmx_fio_do_int(fio, idum);
1034 /* we have removed the NO alternative at the beginning */
1038 ir->epct = epctISOTROPIC;
1042 ir->epc = epcBERENDSEN;
1047 gmx_fio_do_int(fio, ir->epc);
1048 gmx_fio_do_int(fio, ir->epct);
1050 if (file_version >= 71)
1052 gmx_fio_do_int(fio, ir->nstpcouple);
1056 ir->nstpcouple = ir->nstcalcenergy;
1058 gmx_fio_do_real(fio, ir->tau_p);
1059 if (file_version <= 15)
1061 gmx_fio_do_rvec(fio, vdum);
1062 clear_mat(ir->ref_p);
1063 for (i = 0; i < DIM; i++)
1065 ir->ref_p[i][i] = vdum[i];
1070 gmx_fio_do_rvec(fio, ir->ref_p[XX]);
1071 gmx_fio_do_rvec(fio, ir->ref_p[YY]);
1072 gmx_fio_do_rvec(fio, ir->ref_p[ZZ]);
1074 if (file_version <= 15)
1076 gmx_fio_do_rvec(fio, vdum);
1077 clear_mat(ir->compress);
1078 for (i = 0; i < DIM; i++)
1080 ir->compress[i][i] = vdum[i];
1085 gmx_fio_do_rvec(fio, ir->compress[XX]);
1086 gmx_fio_do_rvec(fio, ir->compress[YY]);
1087 gmx_fio_do_rvec(fio, ir->compress[ZZ]);
1089 if (file_version >= 47)
1091 gmx_fio_do_int(fio, ir->refcoord_scaling);
1092 gmx_fio_do_rvec(fio, ir->posres_com);
1093 gmx_fio_do_rvec(fio, ir->posres_comB);
1097 ir->refcoord_scaling = erscNO;
1098 clear_rvec(ir->posres_com);
1099 clear_rvec(ir->posres_comB);
1101 if ((file_version > 25) && (file_version < 79))
1103 gmx_fio_do_int(fio, ir->andersen_seed);
1107 ir->andersen_seed = 0;
1109 if (file_version < 26)
1111 gmx_fio_do_gmx_bool(fio, bSimAnn);
1112 gmx_fio_do_real(fio, zerotemptime);
1115 if (file_version < 37)
1117 gmx_fio_do_real(fio, rdum);
1120 gmx_fio_do_real(fio, ir->shake_tol);
1121 if (file_version < 54)
1123 gmx_fio_do_real(fio, *fudgeQQ);
1126 gmx_fio_do_int(fio, ir->efep);
1127 if (file_version <= 14 && ir->efep != efepNO)
1131 do_fepvals(fio, ir->fepvals, bRead, file_version);
1133 if (file_version >= 79)
1135 gmx_fio_do_gmx_bool(fio, ir->bSimTemp);
1138 ir->bSimTemp = TRUE;
1143 ir->bSimTemp = FALSE;
1147 do_simtempvals(fio, ir->simtempvals, ir->fepvals->n_lambda, bRead, file_version);
1150 if (file_version >= 79)
1152 gmx_fio_do_gmx_bool(fio, ir->bExpanded);
1155 ir->bExpanded = TRUE;
1159 ir->bExpanded = FALSE;
1164 do_expandedvals(fio, ir->expandedvals, ir->fepvals, bRead, file_version);
1166 if (file_version >= 57)
1168 gmx_fio_do_int(fio, ir->eDisre);
1170 gmx_fio_do_int(fio, ir->eDisreWeighting);
1171 if (file_version < 22)
1173 if (ir->eDisreWeighting == 0)
1175 ir->eDisreWeighting = edrwEqual;
1179 ir->eDisreWeighting = edrwConservative;
1182 gmx_fio_do_gmx_bool(fio, ir->bDisreMixed);
1183 gmx_fio_do_real(fio, ir->dr_fc);
1184 gmx_fio_do_real(fio, ir->dr_tau);
1185 gmx_fio_do_int(fio, ir->nstdisreout);
1186 if (file_version >= 22)
1188 gmx_fio_do_real(fio, ir->orires_fc);
1189 gmx_fio_do_real(fio, ir->orires_tau);
1190 gmx_fio_do_int(fio, ir->nstorireout);
1196 ir->nstorireout = 0;
1198 if (file_version >= 26 && file_version < 79)
1200 gmx_fio_do_real(fio, ir->dihre_fc);
1201 if (file_version < 56)
1203 gmx_fio_do_real(fio, rdum);
1204 gmx_fio_do_int(fio, idum);
1212 gmx_fio_do_real(fio, ir->em_stepsize);
1213 gmx_fio_do_real(fio, ir->em_tol);
1214 if (file_version >= 22)
1216 gmx_fio_do_gmx_bool(fio, ir->bShakeSOR);
1220 ir->bShakeSOR = TRUE;
1222 if (file_version >= 11)
1224 gmx_fio_do_int(fio, ir->niter);
1229 fprintf(stderr, "Note: niter not in run input file, setting it to %d\n",
1232 if (file_version >= 21)
1234 gmx_fio_do_real(fio, ir->fc_stepsize);
1238 ir->fc_stepsize = 0;
1240 gmx_fio_do_int(fio, ir->eConstrAlg);
1241 gmx_fio_do_int(fio, ir->nProjOrder);
1242 gmx_fio_do_real(fio, ir->LincsWarnAngle);
1243 if (file_version <= 14)
1245 gmx_fio_do_int(fio, idum);
1247 if (file_version >= 26)
1249 gmx_fio_do_int(fio, ir->nLincsIter);
1254 fprintf(stderr, "Note: nLincsIter not in run input file, setting it to %d\n",
1257 if (file_version < 33)
1259 gmx_fio_do_real(fio, bd_temp);
1261 gmx_fio_do_real(fio, ir->bd_fric);
1262 gmx_fio_do_int(fio, ir->ld_seed);
1263 if (file_version >= 33)
1265 for (i = 0; i < DIM; i++)
1267 gmx_fio_do_rvec(fio, ir->deform[i]);
1272 for (i = 0; i < DIM; i++)
1274 clear_rvec(ir->deform[i]);
1277 if (file_version >= 14)
1279 gmx_fio_do_real(fio, ir->cos_accel);
1285 gmx_fio_do_int(fio, ir->userint1);
1286 gmx_fio_do_int(fio, ir->userint2);
1287 gmx_fio_do_int(fio, ir->userint3);
1288 gmx_fio_do_int(fio, ir->userint4);
1289 gmx_fio_do_real(fio, ir->userreal1);
1290 gmx_fio_do_real(fio, ir->userreal2);
1291 gmx_fio_do_real(fio, ir->userreal3);
1292 gmx_fio_do_real(fio, ir->userreal4);
1295 if (file_version >= 77)
1297 gmx_fio_do_gmx_bool(fio, ir->bAdress);
1302 snew(ir->adress, 1);
1304 gmx_fio_do_int(fio, ir->adress->type);
1305 gmx_fio_do_real(fio, ir->adress->const_wf);
1306 gmx_fio_do_real(fio, ir->adress->ex_width);
1307 gmx_fio_do_real(fio, ir->adress->hy_width);
1308 gmx_fio_do_int(fio, ir->adress->icor);
1309 gmx_fio_do_int(fio, ir->adress->site);
1310 gmx_fio_do_rvec(fio, ir->adress->refs);
1311 gmx_fio_do_int(fio, ir->adress->n_tf_grps);
1312 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
1313 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
1314 gmx_fio_do_int(fio, ir->adress->do_hybridpairs);
1318 snew(ir->adress->tf_table_index, ir->adress->n_tf_grps);
1320 if (ir->adress->n_tf_grps > 0)
1322 bDum = gmx_fio_ndo_int(fio, ir->adress->tf_table_index, ir->adress->n_tf_grps);
1326 snew(ir->adress->group_explicit, ir->adress->n_energy_grps);
1328 if (ir->adress->n_energy_grps > 0)
1330 bDum = gmx_fio_ndo_int(fio, ir->adress->group_explicit, ir->adress->n_energy_grps);
1336 ir->bAdress = FALSE;
1340 if (file_version >= 48)
1342 gmx_fio_do_int(fio, ir->ePull);
1343 if (ir->ePull != epullNO)
1349 do_pull(fio, ir->pull, bRead, file_version);
1354 ir->ePull = epullNO;
1357 /* Enforced rotation */
1358 if (file_version >= 74)
1360 gmx_fio_do_int(fio, ir->bRot);
1361 if (ir->bRot == TRUE)
1367 do_rot(fio, ir->rot, bRead, file_version);
1376 gmx_fio_do_int(fio, ir->opts.ngtc);
1377 if (file_version >= 69)
1379 gmx_fio_do_int(fio, ir->opts.nhchainlength);
1383 ir->opts.nhchainlength = 1;
1385 gmx_fio_do_int(fio, ir->opts.ngacc);
1386 gmx_fio_do_int(fio, ir->opts.ngfrz);
1387 gmx_fio_do_int(fio, ir->opts.ngener);
1391 snew(ir->opts.nrdf, ir->opts.ngtc);
1392 snew(ir->opts.ref_t, ir->opts.ngtc);
1393 snew(ir->opts.annealing, ir->opts.ngtc);
1394 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
1395 snew(ir->opts.anneal_time, ir->opts.ngtc);
1396 snew(ir->opts.anneal_temp, ir->opts.ngtc);
1397 snew(ir->opts.tau_t, ir->opts.ngtc);
1398 snew(ir->opts.nFreeze, ir->opts.ngfrz);
1399 snew(ir->opts.acc, ir->opts.ngacc);
1400 snew(ir->opts.egp_flags, ir->opts.ngener*ir->opts.ngener);
1402 if (ir->opts.ngtc > 0)
1404 if (bRead && file_version < 13)
1406 snew(tmp, ir->opts.ngtc);
1407 bDum = gmx_fio_ndo_int(fio, tmp, ir->opts.ngtc);
1408 for (i = 0; i < ir->opts.ngtc; i++)
1410 ir->opts.nrdf[i] = tmp[i];
1416 bDum = gmx_fio_ndo_real(fio, ir->opts.nrdf, ir->opts.ngtc);
1418 bDum = gmx_fio_ndo_real(fio, ir->opts.ref_t, ir->opts.ngtc);
1419 bDum = gmx_fio_ndo_real(fio, ir->opts.tau_t, ir->opts.ngtc);
1420 if (file_version < 33 && ir->eI == eiBD)
1422 for (i = 0; i < ir->opts.ngtc; i++)
1424 ir->opts.tau_t[i] = bd_temp;
1428 if (ir->opts.ngfrz > 0)
1430 bDum = gmx_fio_ndo_ivec(fio, ir->opts.nFreeze, ir->opts.ngfrz);
1432 if (ir->opts.ngacc > 0)
1434 gmx_fio_ndo_rvec(fio, ir->opts.acc, ir->opts.ngacc);
1436 if (file_version >= 12)
1438 bDum = gmx_fio_ndo_int(fio, ir->opts.egp_flags,
1439 ir->opts.ngener*ir->opts.ngener);
1442 if (bRead && file_version < 26)
1444 for (i = 0; i < ir->opts.ngtc; i++)
1448 ir->opts.annealing[i] = eannSINGLE;
1449 ir->opts.anneal_npoints[i] = 2;
1450 snew(ir->opts.anneal_time[i], 2);
1451 snew(ir->opts.anneal_temp[i], 2);
1452 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
1453 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
1454 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
1455 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
1456 ir->opts.anneal_time[i][0] = ir->init_t;
1457 ir->opts.anneal_time[i][1] = finish_t;
1458 ir->opts.anneal_temp[i][0] = init_temp;
1459 ir->opts.anneal_temp[i][1] = finish_temp;
1463 ir->opts.annealing[i] = eannNO;
1464 ir->opts.anneal_npoints[i] = 0;
1470 /* file version 26 or later */
1471 /* First read the lists with annealing and npoints for each group */
1472 bDum = gmx_fio_ndo_int(fio, ir->opts.annealing, ir->opts.ngtc);
1473 bDum = gmx_fio_ndo_int(fio, ir->opts.anneal_npoints, ir->opts.ngtc);
1474 for (j = 0; j < (ir->opts.ngtc); j++)
1476 k = ir->opts.anneal_npoints[j];
1479 snew(ir->opts.anneal_time[j], k);
1480 snew(ir->opts.anneal_temp[j], k);
1482 bDum = gmx_fio_ndo_real(fio, ir->opts.anneal_time[j], k);
1483 bDum = gmx_fio_ndo_real(fio, ir->opts.anneal_temp[j], k);
1487 if (file_version >= 45)
1489 gmx_fio_do_int(fio, ir->nwall);
1490 gmx_fio_do_int(fio, ir->wall_type);
1491 if (file_version >= 50)
1493 gmx_fio_do_real(fio, ir->wall_r_linpot);
1497 ir->wall_r_linpot = -1;
1499 gmx_fio_do_int(fio, ir->wall_atomtype[0]);
1500 gmx_fio_do_int(fio, ir->wall_atomtype[1]);
1501 gmx_fio_do_real(fio, ir->wall_density[0]);
1502 gmx_fio_do_real(fio, ir->wall_density[1]);
1503 gmx_fio_do_real(fio, ir->wall_ewald_zfac);
1509 ir->wall_atomtype[0] = -1;
1510 ir->wall_atomtype[1] = -1;
1511 ir->wall_density[0] = 0;
1512 ir->wall_density[1] = 0;
1513 ir->wall_ewald_zfac = 3;
1515 /* Cosine stuff for electric fields */
1516 for (j = 0; (j < DIM); j++)
1518 gmx_fio_do_int(fio, ir->ex[j].n);
1519 gmx_fio_do_int(fio, ir->et[j].n);
1522 snew(ir->ex[j].a, ir->ex[j].n);
1523 snew(ir->ex[j].phi, ir->ex[j].n);
1524 snew(ir->et[j].a, ir->et[j].n);
1525 snew(ir->et[j].phi, ir->et[j].n);
1527 bDum = gmx_fio_ndo_real(fio, ir->ex[j].a, ir->ex[j].n);
1528 bDum = gmx_fio_ndo_real(fio, ir->ex[j].phi, ir->ex[j].n);
1529 bDum = gmx_fio_ndo_real(fio, ir->et[j].a, ir->et[j].n);
1530 bDum = gmx_fio_ndo_real(fio, ir->et[j].phi, ir->et[j].n);
1534 if (file_version >= 39)
1536 gmx_fio_do_gmx_bool(fio, ir->bQMMM);
1537 gmx_fio_do_int(fio, ir->QMMMscheme);
1538 gmx_fio_do_real(fio, ir->scalefactor);
1539 gmx_fio_do_int(fio, ir->opts.ngQM);
1542 snew(ir->opts.QMmethod, ir->opts.ngQM);
1543 snew(ir->opts.QMbasis, ir->opts.ngQM);
1544 snew(ir->opts.QMcharge, ir->opts.ngQM);
1545 snew(ir->opts.QMmult, ir->opts.ngQM);
1546 snew(ir->opts.bSH, ir->opts.ngQM);
1547 snew(ir->opts.CASorbitals, ir->opts.ngQM);
1548 snew(ir->opts.CASelectrons, ir->opts.ngQM);
1549 snew(ir->opts.SAon, ir->opts.ngQM);
1550 snew(ir->opts.SAoff, ir->opts.ngQM);
1551 snew(ir->opts.SAsteps, ir->opts.ngQM);
1552 snew(ir->opts.bOPT, ir->opts.ngQM);
1553 snew(ir->opts.bTS, ir->opts.ngQM);
1555 if (ir->opts.ngQM > 0)
1557 bDum = gmx_fio_ndo_int(fio, ir->opts.QMmethod, ir->opts.ngQM);
1558 bDum = gmx_fio_ndo_int(fio, ir->opts.QMbasis, ir->opts.ngQM);
1559 bDum = gmx_fio_ndo_int(fio, ir->opts.QMcharge, ir->opts.ngQM);
1560 bDum = gmx_fio_ndo_int(fio, ir->opts.QMmult, ir->opts.ngQM);
1561 bDum = gmx_fio_ndo_gmx_bool(fio, ir->opts.bSH, ir->opts.ngQM);
1562 bDum = gmx_fio_ndo_int(fio, ir->opts.CASorbitals, ir->opts.ngQM);
1563 bDum = gmx_fio_ndo_int(fio, ir->opts.CASelectrons, ir->opts.ngQM);
1564 bDum = gmx_fio_ndo_real(fio, ir->opts.SAon, ir->opts.ngQM);
1565 bDum = gmx_fio_ndo_real(fio, ir->opts.SAoff, ir->opts.ngQM);
1566 bDum = gmx_fio_ndo_int(fio, ir->opts.SAsteps, ir->opts.ngQM);
1567 bDum = gmx_fio_ndo_gmx_bool(fio, ir->opts.bOPT, ir->opts.ngQM);
1568 bDum = gmx_fio_ndo_gmx_bool(fio, ir->opts.bTS, ir->opts.ngQM);
1570 /* end of QMMM stuff */
1575 static void do_harm(t_fileio *fio, t_iparams *iparams, gmx_bool bRead)
1577 gmx_fio_do_real(fio, iparams->harmonic.rA);
1578 gmx_fio_do_real(fio, iparams->harmonic.krA);
1579 gmx_fio_do_real(fio, iparams->harmonic.rB);
1580 gmx_fio_do_real(fio, iparams->harmonic.krB);
1583 void do_iparams(t_fileio *fio, t_functype ftype, t_iparams *iparams,
1584 gmx_bool bRead, int file_version)
1592 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1602 do_harm(fio, iparams, bRead);
1603 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead)
1605 /* Correct incorrect storage of parameters */
1606 iparams->pdihs.phiB = iparams->pdihs.phiA;
1607 iparams->pdihs.cpB = iparams->pdihs.cpA;
1610 case F_LINEAR_ANGLES:
1611 gmx_fio_do_real(fio, iparams->linangle.klinA);
1612 gmx_fio_do_real(fio, iparams->linangle.aA);
1613 gmx_fio_do_real(fio, iparams->linangle.klinB);
1614 gmx_fio_do_real(fio, iparams->linangle.aB);
1617 gmx_fio_do_real(fio, iparams->fene.bm);
1618 gmx_fio_do_real(fio, iparams->fene.kb);
1621 gmx_fio_do_real(fio, iparams->restraint.lowA);
1622 gmx_fio_do_real(fio, iparams->restraint.up1A);
1623 gmx_fio_do_real(fio, iparams->restraint.up2A);
1624 gmx_fio_do_real(fio, iparams->restraint.kA);
1625 gmx_fio_do_real(fio, iparams->restraint.lowB);
1626 gmx_fio_do_real(fio, iparams->restraint.up1B);
1627 gmx_fio_do_real(fio, iparams->restraint.up2B);
1628 gmx_fio_do_real(fio, iparams->restraint.kB);
1634 gmx_fio_do_real(fio, iparams->tab.kA);
1635 gmx_fio_do_int(fio, iparams->tab.table);
1636 gmx_fio_do_real(fio, iparams->tab.kB);
1638 case F_CROSS_BOND_BONDS:
1639 gmx_fio_do_real(fio, iparams->cross_bb.r1e);
1640 gmx_fio_do_real(fio, iparams->cross_bb.r2e);
1641 gmx_fio_do_real(fio, iparams->cross_bb.krr);
1643 case F_CROSS_BOND_ANGLES:
1644 gmx_fio_do_real(fio, iparams->cross_ba.r1e);
1645 gmx_fio_do_real(fio, iparams->cross_ba.r2e);
1646 gmx_fio_do_real(fio, iparams->cross_ba.r3e);
1647 gmx_fio_do_real(fio, iparams->cross_ba.krt);
1649 case F_UREY_BRADLEY:
1650 gmx_fio_do_real(fio, iparams->u_b.thetaA);
1651 gmx_fio_do_real(fio, iparams->u_b.kthetaA);
1652 gmx_fio_do_real(fio, iparams->u_b.r13A);
1653 gmx_fio_do_real(fio, iparams->u_b.kUBA);
1654 if (file_version >= 79)
1656 gmx_fio_do_real(fio, iparams->u_b.thetaB);
1657 gmx_fio_do_real(fio, iparams->u_b.kthetaB);
1658 gmx_fio_do_real(fio, iparams->u_b.r13B);
1659 gmx_fio_do_real(fio, iparams->u_b.kUBB);
1663 iparams->u_b.thetaB = iparams->u_b.thetaA;
1664 iparams->u_b.kthetaB = iparams->u_b.kthetaA;
1665 iparams->u_b.r13B = iparams->u_b.r13A;
1666 iparams->u_b.kUBB = iparams->u_b.kUBA;
1669 case F_QUARTIC_ANGLES:
1670 gmx_fio_do_real(fio, iparams->qangle.theta);
1671 bDum = gmx_fio_ndo_real(fio, iparams->qangle.c, 5);
1674 gmx_fio_do_real(fio, iparams->bham.a);
1675 gmx_fio_do_real(fio, iparams->bham.b);
1676 gmx_fio_do_real(fio, iparams->bham.c);
1679 gmx_fio_do_real(fio, iparams->morse.b0A);
1680 gmx_fio_do_real(fio, iparams->morse.cbA);
1681 gmx_fio_do_real(fio, iparams->morse.betaA);
1682 if (file_version >= 79)
1684 gmx_fio_do_real(fio, iparams->morse.b0B);
1685 gmx_fio_do_real(fio, iparams->morse.cbB);
1686 gmx_fio_do_real(fio, iparams->morse.betaB);
1690 iparams->morse.b0B = iparams->morse.b0A;
1691 iparams->morse.cbB = iparams->morse.cbA;
1692 iparams->morse.betaB = iparams->morse.betaA;
1696 gmx_fio_do_real(fio, iparams->cubic.b0);
1697 gmx_fio_do_real(fio, iparams->cubic.kb);
1698 gmx_fio_do_real(fio, iparams->cubic.kcub);
1702 case F_POLARIZATION:
1703 gmx_fio_do_real(fio, iparams->polarize.alpha);
1706 gmx_fio_do_real(fio, iparams->anharm_polarize.alpha);
1707 gmx_fio_do_real(fio, iparams->anharm_polarize.drcut);
1708 gmx_fio_do_real(fio, iparams->anharm_polarize.khyp);
1711 if (file_version < 31)
1713 gmx_fatal(FARGS, "Old tpr files with water_polarization not supported. Make a new.");
1715 gmx_fio_do_real(fio, iparams->wpol.al_x);
1716 gmx_fio_do_real(fio, iparams->wpol.al_y);
1717 gmx_fio_do_real(fio, iparams->wpol.al_z);
1718 gmx_fio_do_real(fio, iparams->wpol.rOH);
1719 gmx_fio_do_real(fio, iparams->wpol.rHH);
1720 gmx_fio_do_real(fio, iparams->wpol.rOD);
1723 gmx_fio_do_real(fio, iparams->thole.a);
1724 gmx_fio_do_real(fio, iparams->thole.alpha1);
1725 gmx_fio_do_real(fio, iparams->thole.alpha2);
1726 gmx_fio_do_real(fio, iparams->thole.rfac);
1729 gmx_fio_do_real(fio, iparams->lj.c6);
1730 gmx_fio_do_real(fio, iparams->lj.c12);
1733 gmx_fio_do_real(fio, iparams->lj14.c6A);
1734 gmx_fio_do_real(fio, iparams->lj14.c12A);
1735 gmx_fio_do_real(fio, iparams->lj14.c6B);
1736 gmx_fio_do_real(fio, iparams->lj14.c12B);
1739 gmx_fio_do_real(fio, iparams->ljc14.fqq);
1740 gmx_fio_do_real(fio, iparams->ljc14.qi);
1741 gmx_fio_do_real(fio, iparams->ljc14.qj);
1742 gmx_fio_do_real(fio, iparams->ljc14.c6);
1743 gmx_fio_do_real(fio, iparams->ljc14.c12);
1745 case F_LJC_PAIRS_NB:
1746 gmx_fio_do_real(fio, iparams->ljcnb.qi);
1747 gmx_fio_do_real(fio, iparams->ljcnb.qj);
1748 gmx_fio_do_real(fio, iparams->ljcnb.c6);
1749 gmx_fio_do_real(fio, iparams->ljcnb.c12);
1755 gmx_fio_do_real(fio, iparams->pdihs.phiA);
1756 gmx_fio_do_real(fio, iparams->pdihs.cpA);
1757 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42)
1759 /* Read the incorrectly stored multiplicity */
1760 gmx_fio_do_real(fio, iparams->harmonic.rB);
1761 gmx_fio_do_real(fio, iparams->harmonic.krB);
1762 iparams->pdihs.phiB = iparams->pdihs.phiA;
1763 iparams->pdihs.cpB = iparams->pdihs.cpA;
1767 gmx_fio_do_real(fio, iparams->pdihs.phiB);
1768 gmx_fio_do_real(fio, iparams->pdihs.cpB);
1769 gmx_fio_do_int(fio, iparams->pdihs.mult);
1773 gmx_fio_do_int(fio, iparams->disres.label);
1774 gmx_fio_do_int(fio, iparams->disres.type);
1775 gmx_fio_do_real(fio, iparams->disres.low);
1776 gmx_fio_do_real(fio, iparams->disres.up1);
1777 gmx_fio_do_real(fio, iparams->disres.up2);
1778 gmx_fio_do_real(fio, iparams->disres.kfac);
1781 gmx_fio_do_int(fio, iparams->orires.ex);
1782 gmx_fio_do_int(fio, iparams->orires.label);
1783 gmx_fio_do_int(fio, iparams->orires.power);
1784 gmx_fio_do_real(fio, iparams->orires.c);
1785 gmx_fio_do_real(fio, iparams->orires.obs);
1786 gmx_fio_do_real(fio, iparams->orires.kfac);
1789 if (file_version < 72)
1791 gmx_fio_do_int(fio, idum);
1792 gmx_fio_do_int(fio, idum);
1794 gmx_fio_do_real(fio, iparams->dihres.phiA);
1795 gmx_fio_do_real(fio, iparams->dihres.dphiA);
1796 gmx_fio_do_real(fio, iparams->dihres.kfacA);
1797 if (file_version >= 72)
1799 gmx_fio_do_real(fio, iparams->dihres.phiB);
1800 gmx_fio_do_real(fio, iparams->dihres.dphiB);
1801 gmx_fio_do_real(fio, iparams->dihres.kfacB);
1805 iparams->dihres.phiB = iparams->dihres.phiA;
1806 iparams->dihres.dphiB = iparams->dihres.dphiA;
1807 iparams->dihres.kfacB = iparams->dihres.kfacA;
1811 gmx_fio_do_rvec(fio, iparams->posres.pos0A);
1812 gmx_fio_do_rvec(fio, iparams->posres.fcA);
1813 if (bRead && file_version < 27)
1815 copy_rvec(iparams->posres.pos0A, iparams->posres.pos0B);
1816 copy_rvec(iparams->posres.fcA, iparams->posres.fcB);
1820 gmx_fio_do_rvec(fio, iparams->posres.pos0B);
1821 gmx_fio_do_rvec(fio, iparams->posres.fcB);
1825 gmx_fio_do_int(fio, iparams->fbposres.geom);
1826 gmx_fio_do_rvec(fio, iparams->fbposres.pos0);
1827 gmx_fio_do_real(fio, iparams->fbposres.r);
1828 gmx_fio_do_real(fio, iparams->fbposres.k);
1831 bDum = gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1832 if (file_version >= 25)
1834 bDum = gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1838 /* Fourier dihedrals are internally represented
1839 * as Ryckaert-Bellemans since those are faster to compute.
1841 bDum = gmx_fio_ndo_real(fio, iparams->rbdihs.rbcA, NR_RBDIHS);
1842 bDum = gmx_fio_ndo_real(fio, iparams->rbdihs.rbcB, NR_RBDIHS);
1846 gmx_fio_do_real(fio, iparams->constr.dA);
1847 gmx_fio_do_real(fio, iparams->constr.dB);
1850 gmx_fio_do_real(fio, iparams->settle.doh);
1851 gmx_fio_do_real(fio, iparams->settle.dhh);
1854 gmx_fio_do_real(fio, iparams->vsite.a);
1859 gmx_fio_do_real(fio, iparams->vsite.a);
1860 gmx_fio_do_real(fio, iparams->vsite.b);
1865 gmx_fio_do_real(fio, iparams->vsite.a);
1866 gmx_fio_do_real(fio, iparams->vsite.b);
1867 gmx_fio_do_real(fio, iparams->vsite.c);
1870 gmx_fio_do_int(fio, iparams->vsiten.n);
1871 gmx_fio_do_real(fio, iparams->vsiten.a);
1876 /* We got rid of some parameters in version 68 */
1877 if (bRead && file_version < 68)
1879 gmx_fio_do_real(fio, rdum);
1880 gmx_fio_do_real(fio, rdum);
1881 gmx_fio_do_real(fio, rdum);
1882 gmx_fio_do_real(fio, rdum);
1884 gmx_fio_do_real(fio, iparams->gb.sar);
1885 gmx_fio_do_real(fio, iparams->gb.st);
1886 gmx_fio_do_real(fio, iparams->gb.pi);
1887 gmx_fio_do_real(fio, iparams->gb.gbr);
1888 gmx_fio_do_real(fio, iparams->gb.bmlt);
1891 gmx_fio_do_int(fio, iparams->cmap.cmapA);
1892 gmx_fio_do_int(fio, iparams->cmap.cmapB);
1895 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1896 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1900 gmx_fio_unset_comment(fio);
1904 static void do_ilist(t_fileio *fio, t_ilist *ilist, gmx_bool bRead, int file_version,
1908 gmx_bool bDum = TRUE;
1912 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1914 if (file_version < 44)
1916 for (i = 0; i < MAXNODES; i++)
1918 gmx_fio_do_int(fio, idum);
1921 gmx_fio_do_int(fio, ilist->nr);
1924 snew(ilist->iatoms, ilist->nr);
1926 bDum = gmx_fio_ndo_int(fio, ilist->iatoms, ilist->nr);
1929 gmx_fio_unset_comment(fio);
1933 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1934 gmx_bool bRead, int file_version)
1937 gmx_bool bDum = TRUE;
1940 gmx_fio_do_int(fio, ffparams->atnr);
1941 if (file_version < 57)
1943 gmx_fio_do_int(fio, idum);
1945 gmx_fio_do_int(fio, ffparams->ntypes);
1948 fprintf(debug, "ffparams->atnr = %d, ntypes = %d\n",
1949 ffparams->atnr, ffparams->ntypes);
1953 snew(ffparams->functype, ffparams->ntypes);
1954 snew(ffparams->iparams, ffparams->ntypes);
1956 /* Read/write all the function types */
1957 bDum = gmx_fio_ndo_int(fio, ffparams->functype, ffparams->ntypes);
1960 pr_ivec(debug, 0, "functype", ffparams->functype, ffparams->ntypes, TRUE);
1963 if (file_version >= 66)
1965 gmx_fio_do_double(fio, ffparams->reppow);
1969 ffparams->reppow = 12.0;
1972 if (file_version >= 57)
1974 gmx_fio_do_real(fio, ffparams->fudgeQQ);
1977 /* Check whether all these function types are supported by the code.
1978 * In practice the code is backwards compatible, which means that the
1979 * numbering may have to be altered from old numbering to new numbering
1981 for (i = 0; (i < ffparams->ntypes); i++)
1985 /* Loop over file versions */
1986 for (k = 0; (k < NFTUPD); k++)
1988 /* Compare the read file_version to the update table */
1989 if ((file_version < ftupd[k].fvnr) &&
1990 (ffparams->functype[i] >= ftupd[k].ftype))
1992 ffparams->functype[i] += 1;
1995 fprintf(debug, "Incrementing function type %d to %d (due to %s)\n",
1996 i, ffparams->functype[i],
1997 interaction_function[ftupd[k].ftype].longname);
2004 do_iparams(fio, ffparams->functype[i], &ffparams->iparams[i], bRead,
2008 pr_iparams(debug, ffparams->functype[i], &ffparams->iparams[i]);
2013 static void add_settle_atoms(t_ilist *ilist)
2017 /* Settle used to only store the first atom: add the other two */
2018 srenew(ilist->iatoms, 2*ilist->nr);
2019 for (i = ilist->nr/2-1; i >= 0; i--)
2021 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
2022 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
2023 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
2024 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
2026 ilist->nr = 2*ilist->nr;
2029 static void do_ilists(t_fileio *fio, t_ilist *ilist, gmx_bool bRead,
2032 int i, j, renum[F_NRE];
2033 gmx_bool bDum = TRUE, bClear;
2036 for (j = 0; (j < F_NRE); j++)
2041 for (k = 0; k < NFTUPD; k++)
2043 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
2052 ilist[j].iatoms = NULL;
2056 do_ilist(fio, &ilist[j], bRead, file_version, j);
2057 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
2059 add_settle_atoms(&ilist[j]);
2063 if (bRead && gmx_debug_at)
2064 pr_ilist(debug,0,interaction_function[j].longname,
2065 functype,&ilist[j],TRUE);
2070 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams, gmx_moltype_t *molt,
2071 gmx_bool bRead, int file_version)
2073 do_ffparams(fio, ffparams, bRead, file_version);
2075 if (file_version >= 54)
2077 gmx_fio_do_real(fio, ffparams->fudgeQQ);
2080 do_ilists(fio, molt->ilist, bRead, file_version);
2083 static void do_block(t_fileio *fio, t_block *block, gmx_bool bRead, int file_version)
2085 int i, idum, dum_nra, *dum_a;
2086 gmx_bool bDum = TRUE;
2088 if (file_version < 44)
2090 for (i = 0; i < MAXNODES; i++)
2092 gmx_fio_do_int(fio, idum);
2095 gmx_fio_do_int(fio, block->nr);
2096 if (file_version < 51)
2098 gmx_fio_do_int(fio, dum_nra);
2102 block->nalloc_index = block->nr+1;
2103 snew(block->index, block->nalloc_index);
2105 bDum = gmx_fio_ndo_int(fio, block->index, block->nr+1);
2107 if (file_version < 51 && dum_nra > 0)
2109 snew(dum_a, dum_nra);
2110 bDum = gmx_fio_ndo_int(fio, dum_a, dum_nra);
2115 static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
2119 gmx_bool bDum = TRUE;
2121 if (file_version < 44)
2123 for (i = 0; i < MAXNODES; i++)
2125 gmx_fio_do_int(fio, idum);
2128 gmx_fio_do_int(fio, block->nr);
2129 gmx_fio_do_int(fio, block->nra);
2132 block->nalloc_index = block->nr+1;
2133 snew(block->index, block->nalloc_index);
2134 block->nalloc_a = block->nra;
2135 snew(block->a, block->nalloc_a);
2137 bDum = gmx_fio_ndo_int(fio, block->index, block->nr+1);
2138 bDum = gmx_fio_ndo_int(fio, block->a, block->nra);
2141 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
2142 int file_version, gmx_groups_t *groups, int atnr)
2146 gmx_fio_do_real(fio, atom->m);
2147 gmx_fio_do_real(fio, atom->q);
2148 gmx_fio_do_real(fio, atom->mB);
2149 gmx_fio_do_real(fio, atom->qB);
2150 gmx_fio_do_ushort(fio, atom->type);
2151 gmx_fio_do_ushort(fio, atom->typeB);
2152 gmx_fio_do_int(fio, atom->ptype);
2153 gmx_fio_do_int(fio, atom->resind);
2154 if (file_version >= 52)
2156 gmx_fio_do_int(fio, atom->atomnumber);
2160 atom->atomnumber = NOTSET;
2162 if (file_version < 23)
2166 else if (file_version < 39)
2175 if (file_version < 57)
2177 unsigned char uchar[egcNR];
2178 gmx_fio_ndo_uchar(fio, uchar, myngrp);
2179 for (i = myngrp; (i < ngrp); i++)
2183 /* Copy the old data format to the groups struct */
2184 for (i = 0; i < ngrp; i++)
2186 groups->grpnr[i][atnr] = uchar[i];
2191 static void do_grps(t_fileio *fio, int ngrp, t_grps grps[], gmx_bool bRead,
2195 gmx_bool bDum = TRUE;
2197 if (file_version < 23)
2201 else if (file_version < 39)
2210 for (j = 0; (j < ngrp); j++)
2214 gmx_fio_do_int(fio, grps[j].nr);
2217 snew(grps[j].nm_ind, grps[j].nr);
2219 bDum = gmx_fio_ndo_int(fio, grps[j].nm_ind, grps[j].nr);
2224 snew(grps[j].nm_ind, grps[j].nr);
2229 static void do_symstr(t_fileio *fio, char ***nm, gmx_bool bRead, t_symtab *symtab)
2235 gmx_fio_do_int(fio, ls);
2236 *nm = get_symtab_handle(symtab, ls);
2240 ls = lookup_symtab(symtab, *nm);
2241 gmx_fio_do_int(fio, ls);
2245 static void do_strstr(t_fileio *fio, int nstr, char ***nm, gmx_bool bRead,
2250 for (j = 0; (j < nstr); j++)
2252 do_symstr(fio, &(nm[j]), bRead, symtab);
2256 static void do_resinfo(t_fileio *fio, int n, t_resinfo *ri, gmx_bool bRead,
2257 t_symtab *symtab, int file_version)
2261 for (j = 0; (j < n); j++)
2263 do_symstr(fio, &(ri[j].name), bRead, symtab);
2264 if (file_version >= 63)
2266 gmx_fio_do_int(fio, ri[j].nr);
2267 gmx_fio_do_uchar(fio, ri[j].ic);
2277 static void do_atoms(t_fileio *fio, t_atoms *atoms, gmx_bool bRead, t_symtab *symtab,
2279 gmx_groups_t *groups)
2283 gmx_fio_do_int(fio, atoms->nr);
2284 gmx_fio_do_int(fio, atoms->nres);
2285 if (file_version < 57)
2287 gmx_fio_do_int(fio, groups->ngrpname);
2288 for (i = 0; i < egcNR; i++)
2290 groups->ngrpnr[i] = atoms->nr;
2291 snew(groups->grpnr[i], groups->ngrpnr[i]);
2296 snew(atoms->atom, atoms->nr);
2297 snew(atoms->atomname, atoms->nr);
2298 snew(atoms->atomtype, atoms->nr);
2299 snew(atoms->atomtypeB, atoms->nr);
2300 snew(atoms->resinfo, atoms->nres);
2301 if (file_version < 57)
2303 snew(groups->grpname, groups->ngrpname);
2305 atoms->pdbinfo = NULL;
2307 for (i = 0; (i < atoms->nr); i++)
2309 do_atom(fio, &atoms->atom[i], egcNR, bRead, file_version, groups, i);
2311 do_strstr(fio, atoms->nr, atoms->atomname, bRead, symtab);
2312 if (bRead && (file_version <= 20))
2314 for (i = 0; i < atoms->nr; i++)
2316 atoms->atomtype[i] = put_symtab(symtab, "?");
2317 atoms->atomtypeB[i] = put_symtab(symtab, "?");
2322 do_strstr(fio, atoms->nr, atoms->atomtype, bRead, symtab);
2323 do_strstr(fio, atoms->nr, atoms->atomtypeB, bRead, symtab);
2325 do_resinfo(fio, atoms->nres, atoms->resinfo, bRead, symtab, file_version);
2327 if (file_version < 57)
2329 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2331 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2335 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
2336 gmx_bool bRead, t_symtab *symtab,
2340 gmx_bool bDum = TRUE;
2342 do_grps(fio, egcNR, groups->grps, bRead, file_version);
2343 gmx_fio_do_int(fio, groups->ngrpname);
2346 snew(groups->grpname, groups->ngrpname);
2348 do_strstr(fio, groups->ngrpname, groups->grpname, bRead, symtab);
2349 for (g = 0; g < egcNR; g++)
2351 gmx_fio_do_int(fio, groups->ngrpnr[g]);
2352 if (groups->ngrpnr[g] == 0)
2356 groups->grpnr[g] = NULL;
2363 snew(groups->grpnr[g], groups->ngrpnr[g]);
2365 bDum = gmx_fio_ndo_uchar(fio, groups->grpnr[g], groups->ngrpnr[g]);
2370 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes, gmx_bool bRead,
2371 t_symtab *symtab, int file_version)
2374 gmx_bool bDum = TRUE;
2376 if (file_version > 25)
2378 gmx_fio_do_int(fio, atomtypes->nr);
2382 snew(atomtypes->radius, j);
2383 snew(atomtypes->vol, j);
2384 snew(atomtypes->surftens, j);
2385 snew(atomtypes->atomnumber, j);
2386 snew(atomtypes->gb_radius, j);
2387 snew(atomtypes->S_hct, j);
2389 bDum = gmx_fio_ndo_real(fio, atomtypes->radius, j);
2390 bDum = gmx_fio_ndo_real(fio, atomtypes->vol, j);
2391 bDum = gmx_fio_ndo_real(fio, atomtypes->surftens, j);
2392 if (file_version >= 40)
2394 bDum = gmx_fio_ndo_int(fio, atomtypes->atomnumber, j);
2396 if (file_version >= 60)
2398 bDum = gmx_fio_ndo_real(fio, atomtypes->gb_radius, j);
2399 bDum = gmx_fio_ndo_real(fio, atomtypes->S_hct, j);
2404 /* File versions prior to 26 cannot do GBSA,
2405 * so they dont use this structure
2408 atomtypes->radius = NULL;
2409 atomtypes->vol = NULL;
2410 atomtypes->surftens = NULL;
2411 atomtypes->atomnumber = NULL;
2412 atomtypes->gb_radius = NULL;
2413 atomtypes->S_hct = NULL;
2417 static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
2423 gmx_fio_do_int(fio, symtab->nr);
2427 snew(symtab->symbuf, 1);
2428 symbuf = symtab->symbuf;
2429 symbuf->bufsize = nr;
2430 snew(symbuf->buf, nr);
2431 for (i = 0; (i < nr); i++)
2433 gmx_fio_do_string(fio, buf);
2434 symbuf->buf[i] = strdup(buf);
2439 symbuf = symtab->symbuf;
2440 while (symbuf != NULL)
2442 for (i = 0; (i < symbuf->bufsize) && (i < nr); i++)
2444 gmx_fio_do_string(fio, symbuf->buf[i]);
2447 symbuf = symbuf->next;
2451 gmx_fatal(FARGS, "nr of symtab strings left: %d", nr);
2456 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
2458 int i, j, ngrid, gs, nelem;
2460 gmx_fio_do_int(fio, cmap_grid->ngrid);
2461 gmx_fio_do_int(fio, cmap_grid->grid_spacing);
2463 ngrid = cmap_grid->ngrid;
2464 gs = cmap_grid->grid_spacing;
2469 snew(cmap_grid->cmapdata, ngrid);
2471 for (i = 0; i < cmap_grid->ngrid; i++)
2473 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
2477 for (i = 0; i < cmap_grid->ngrid; i++)
2479 for (j = 0; j < nelem; j++)
2481 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
2482 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
2483 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+2]);
2484 gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+3]);
2490 void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
2492 int m, a, a0, a1, r;
2496 /* We always assign a new chain number, but save the chain id characters
2497 * for larger molecules.
2499 #define CHAIN_MIN_ATOMS 15
2503 for (m = 0; m < mols->nr; m++)
2505 a0 = mols->index[m];
2506 a1 = mols->index[m+1];
2507 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
2516 for (a = a0; a < a1; a++)
2518 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
2519 atoms->resinfo[atoms->atom[a].resind].chainid = c;
2524 /* Blank out the chain id if there was only one chain */
2527 for (r = 0; r < atoms->nres; r++)
2529 atoms->resinfo[r].chainid = ' ';
2534 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt, gmx_bool bRead,
2535 t_symtab *symtab, int file_version,
2536 gmx_groups_t *groups)
2540 if (file_version >= 57)
2542 do_symstr(fio, &(molt->name), bRead, symtab);
2545 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
2547 if (bRead && gmx_debug_at)
2549 pr_atoms(debug, 0, "atoms", &molt->atoms, TRUE);
2552 if (file_version >= 57)
2554 do_ilists(fio, molt->ilist, bRead, file_version);
2556 do_block(fio, &molt->cgs, bRead, file_version);
2557 if (bRead && gmx_debug_at)
2559 pr_block(debug, 0, "cgs", &molt->cgs, TRUE);
2563 /* This used to be in the atoms struct */
2564 do_blocka(fio, &molt->excls, bRead, file_version);
2567 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb, gmx_bool bRead,
2572 gmx_fio_do_int(fio, molb->type);
2573 gmx_fio_do_int(fio, molb->nmol);
2574 gmx_fio_do_int(fio, molb->natoms_mol);
2575 /* Position restraint coordinates */
2576 gmx_fio_do_int(fio, molb->nposres_xA);
2577 if (molb->nposres_xA > 0)
2581 snew(molb->posres_xA, molb->nposres_xA);
2583 gmx_fio_ndo_rvec(fio, molb->posres_xA, molb->nposres_xA);
2585 gmx_fio_do_int(fio, molb->nposres_xB);
2586 if (molb->nposres_xB > 0)
2590 snew(molb->posres_xB, molb->nposres_xB);
2592 gmx_fio_ndo_rvec(fio, molb->posres_xB, molb->nposres_xB);
2597 static t_block mtop_mols(gmx_mtop_t *mtop)
2603 for (mb = 0; mb < mtop->nmolblock; mb++)
2605 mols.nr += mtop->molblock[mb].nmol;
2607 mols.nalloc_index = mols.nr + 1;
2608 snew(mols.index, mols.nalloc_index);
2613 for (mb = 0; mb < mtop->nmolblock; mb++)
2615 for (mol = 0; mol < mtop->molblock[mb].nmol; mol++)
2617 a += mtop->molblock[mb].natoms_mol;
2626 static void add_posres_molblock(gmx_mtop_t *mtop)
2631 gmx_molblock_t *molb;
2634 /* posres reference positions are stored in ip->posres (if present) and
2635 in ip->fbposres (if present). If normal and flat-bottomed posres are present,
2636 posres.pos0A are identical to fbposres.pos0. */
2637 il = &mtop->moltype[0].ilist[F_POSRES];
2638 ilfb = &mtop->moltype[0].ilist[F_FBPOSRES];
2639 if (il->nr == 0 && ilfb->nr == 0)
2645 for (i = 0; i < il->nr; i += 2)
2647 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2648 am = max(am, il->iatoms[i+1]);
2649 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
2650 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
2651 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ])
2656 /* This loop is required if we have only flat-bottomed posres:
2658 - bFE == FALSE (no B-state for flat-bottomed posres) */
2661 for (i = 0; i < ilfb->nr; i += 2)
2663 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2664 am = max(am, ilfb->iatoms[i+1]);
2667 /* Make the posres coordinate block end at a molecule end */
2669 while (am >= mtop->mols.index[mol+1])
2673 molb = &mtop->molblock[0];
2674 molb->nposres_xA = mtop->mols.index[mol+1];
2675 snew(molb->posres_xA, molb->nposres_xA);
2678 molb->nposres_xB = molb->nposres_xA;
2679 snew(molb->posres_xB, molb->nposres_xB);
2683 molb->nposres_xB = 0;
2685 for (i = 0; i < il->nr; i += 2)
2687 ip = &mtop->ffparams.iparams[il->iatoms[i]];
2688 a = il->iatoms[i+1];
2689 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
2690 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
2691 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
2694 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
2695 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
2696 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
2701 /* If only flat-bottomed posres are present, take reference pos from them.
2702 Here: bFE == FALSE */
2703 for (i = 0; i < ilfb->nr; i += 2)
2705 ip = &mtop->ffparams.iparams[ilfb->iatoms[i]];
2706 a = ilfb->iatoms[i+1];
2707 molb->posres_xA[a][XX] = ip->fbposres.pos0[XX];
2708 molb->posres_xA[a][YY] = ip->fbposres.pos0[YY];
2709 molb->posres_xA[a][ZZ] = ip->fbposres.pos0[ZZ];
2714 static void set_disres_npair(gmx_mtop_t *mtop)
2721 ip = mtop->ffparams.iparams;
2723 for (mt = 0; mt < mtop->nmoltype; mt++)
2725 il = &mtop->moltype[mt].ilist[F_DISRES];
2730 for (i = 0; i < il->nr; i += 3)
2733 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label)
2735 ip[a[i]].disres.npair = npair;
2743 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
2753 do_symtab(fio, &(mtop->symtab), bRead);
2756 pr_symtab(debug, 0, "symtab", &mtop->symtab);
2759 do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
2761 if (file_version >= 57)
2763 do_ffparams(fio, &mtop->ffparams, bRead, file_version);
2765 gmx_fio_do_int(fio, mtop->nmoltype);
2773 snew(mtop->moltype, mtop->nmoltype);
2774 if (file_version < 57)
2776 mtop->moltype[0].name = mtop->name;
2779 for (mt = 0; mt < mtop->nmoltype; mt++)
2781 do_moltype(fio, &mtop->moltype[mt], bRead, &mtop->symtab, file_version,
2785 if (file_version >= 57)
2787 gmx_fio_do_int(fio, mtop->nmolblock);
2791 mtop->nmolblock = 1;
2795 snew(mtop->molblock, mtop->nmolblock);
2797 if (file_version >= 57)
2799 for (mb = 0; mb < mtop->nmolblock; mb++)
2801 do_molblock(fio, &mtop->molblock[mb], bRead, file_version);
2803 gmx_fio_do_int(fio, mtop->natoms);
2807 mtop->molblock[0].type = 0;
2808 mtop->molblock[0].nmol = 1;
2809 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
2810 mtop->molblock[0].nposres_xA = 0;
2811 mtop->molblock[0].nposres_xB = 0;
2814 do_atomtypes (fio, &(mtop->atomtypes), bRead, &(mtop->symtab), file_version);
2817 pr_atomtypes(debug, 0, "atomtypes", &mtop->atomtypes, TRUE);
2820 if (file_version < 57)
2822 /* Debug statements are inside do_idef */
2823 do_idef (fio, &mtop->ffparams, &mtop->moltype[0], bRead, file_version);
2824 mtop->natoms = mtop->moltype[0].atoms.nr;
2827 if (file_version >= 65)
2829 do_cmap(fio, &mtop->ffparams.cmap_grid, bRead);
2833 mtop->ffparams.cmap_grid.ngrid = 0;
2834 mtop->ffparams.cmap_grid.grid_spacing = 0;
2835 mtop->ffparams.cmap_grid.cmapdata = NULL;
2838 if (file_version >= 57)
2840 do_groups(fio, &mtop->groups, bRead, &(mtop->symtab), file_version);
2843 if (file_version < 57)
2845 do_block(fio, &mtop->moltype[0].cgs, bRead, file_version);
2846 if (bRead && gmx_debug_at)
2848 pr_block(debug, 0, "cgs", &mtop->moltype[0].cgs, TRUE);
2850 do_block(fio, &mtop->mols, bRead, file_version);
2851 /* Add the posres coordinates to the molblock */
2852 add_posres_molblock(mtop);
2856 if (file_version >= 57)
2858 mtop->mols = mtop_mols(mtop);
2862 pr_block(debug, 0, "mols", &mtop->mols, TRUE);
2866 if (file_version < 51)
2868 /* Here used to be the shake blocks */
2869 do_blocka(fio, &dumb, bRead, file_version);
2882 close_symtab(&(mtop->symtab));
2886 /* If TopOnlyOK is TRUE then we can read even future versions
2887 * of tpx files, provided the file_generation hasn't changed.
2888 * If it is FALSE, we need the inputrecord too, and bail out
2889 * if the file is newer than the program.
2891 * The version and generation if the topology (see top of this file)
2892 * are returned in the two last arguments.
2894 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2896 static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx,
2897 gmx_bool TopOnlyOK, int *file_version,
2898 int *file_generation)
2901 char file_tag[STRLEN];
2908 gmx_fio_checktype(fio);
2909 gmx_fio_setdebug(fio, bDebugMode());
2911 /* NEW! XDR tpb file */
2912 precision = sizeof(real);
2915 gmx_fio_do_string(fio, buf);
2916 if (strncmp(buf, "VERSION", 7))
2918 gmx_fatal(FARGS, "Can not read file %s,\n"
2919 " this file is from a Gromacs version which is older than 2.0\n"
2920 " Make a new one with grompp or use a gro or pdb file, if possible",
2921 gmx_fio_getname(fio));
2923 gmx_fio_do_int(fio, precision);
2924 bDouble = (precision == sizeof(double));
2925 if ((precision != sizeof(float)) && !bDouble)
2927 gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes "
2928 "instead of %d or %d",
2929 gmx_fio_getname(fio), precision, sizeof(float), sizeof(double));
2931 gmx_fio_setprecision(fio, bDouble);
2932 fprintf(stderr, "Reading file %s, %s (%s precision)\n",
2933 gmx_fio_getname(fio), buf, bDouble ? "double" : "single");
2937 gmx_fio_write_string(fio, GromacsVersion());
2938 bDouble = (precision == sizeof(double));
2939 gmx_fio_setprecision(fio, bDouble);
2940 gmx_fio_do_int(fio, precision);
2942 sprintf(file_tag, "%s", tpx_tag);
2943 fgen = tpx_generation;
2946 /* Check versions! */
2947 gmx_fio_do_int(fio, fver);
2949 /* This is for backward compatibility with development versions 77-79
2950 * where the tag was, mistakenly, placed before the generation,
2951 * which would cause a segv instead of a proper error message
2952 * when reading the topology only from tpx with <77 code.
2954 if (fver >= 77 && fver <= 79)
2956 gmx_fio_do_string(fio, file_tag);
2961 gmx_fio_do_int(fio, fgen);
2970 gmx_fio_do_string(fio, file_tag);
2976 /* Versions before 77 don't have the tag, set it to release */
2977 sprintf(file_tag, "%s", TPX_TAG_RELEASE);
2980 if (strcmp(file_tag, tpx_tag) != 0)
2982 fprintf(stderr, "Note: file tpx tag '%s', software tpx tag '%s'\n",
2985 /* We only support reading tpx files with the same tag as the code
2986 * or tpx files with the release tag and with lower version number.
2988 if (!strcmp(file_tag, TPX_TAG_RELEASE) == 0 && fver < tpx_version)
2990 gmx_fatal(FARGS, "tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2991 gmx_fio_getname(fio), fver, file_tag,
2992 tpx_version, tpx_tag);
2997 if (file_version != NULL)
2999 *file_version = fver;
3001 if (file_generation != NULL)
3003 *file_generation = fgen;
3007 if ((fver <= tpx_incompatible_version) ||
3008 ((fver > tpx_version) && !TopOnlyOK) ||
3009 (fgen > tpx_generation))
3011 gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program",
3012 gmx_fio_getname(fio), fver, tpx_version);
3015 do_section(fio, eitemHEADER, bRead);
3016 gmx_fio_do_int(fio, tpx->natoms);
3019 gmx_fio_do_int(fio, tpx->ngtc);
3027 gmx_fio_do_int(fio, idum);
3028 gmx_fio_do_real(fio, rdum);
3030 /*a better decision will eventually (5.0 or later) need to be made
3031 on how to treat the alchemical state of the system, which can now
3032 vary through a simulation, and cannot be completely described
3033 though a single lambda variable, or even a single state
3034 index. Eventually, should probably be a vector. MRS*/
3037 gmx_fio_do_int(fio, tpx->fep_state);
3039 gmx_fio_do_real(fio, tpx->lambda);
3040 gmx_fio_do_int(fio, tpx->bIr);
3041 gmx_fio_do_int(fio, tpx->bTop);
3042 gmx_fio_do_int(fio, tpx->bX);
3043 gmx_fio_do_int(fio, tpx->bV);
3044 gmx_fio_do_int(fio, tpx->bF);
3045 gmx_fio_do_int(fio, tpx->bBox);
3047 if ((fgen > tpx_generation))
3049 /* This can only happen if TopOnlyOK=TRUE */
3054 static int do_tpx(t_fileio *fio, gmx_bool bRead,
3055 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop,
3056 gmx_bool bXVallocated)
3061 gmx_bool TopOnlyOK, bDum = TRUE;
3062 int file_version, file_generation;
3066 gmx_bool bPeriodicMols;
3070 tpx.natoms = state->natoms;
3071 tpx.ngtc = state->ngtc; /* need to add nnhpres here? */
3072 tpx.fep_state = state->fep_state;
3073 tpx.lambda = state->lambda[efptFEP];
3074 tpx.bIr = (ir != NULL);
3075 tpx.bTop = (mtop != NULL);
3076 tpx.bX = (state->x != NULL);
3077 tpx.bV = (state->v != NULL);
3078 tpx.bF = (f != NULL);
3082 TopOnlyOK = (ir == NULL);
3084 do_tpxheader(fio, bRead, &tpx, TopOnlyOK, &file_version, &file_generation);
3089 /* state->lambda = tpx.lambda;*/ /*remove this eventually? */
3090 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
3095 init_state(state, 0, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
3096 state->natoms = tpx.natoms;
3097 state->nalloc = tpx.natoms;
3103 init_state(state, tpx.natoms, tpx.ngtc, 0, 0, 0); /* nose-hoover chains */
3107 #define do_test(fio, b, p) if (bRead && (p != NULL) && !b) gmx_fatal(FARGS, "No %s in %s",#p, gmx_fio_getname(fio))
3109 do_test(fio, tpx.bBox, state->box);
3110 do_section(fio, eitemBOX, bRead);
3113 gmx_fio_ndo_rvec(fio, state->box, DIM);
3114 if (file_version >= 51)
3116 gmx_fio_ndo_rvec(fio, state->box_rel, DIM);
3120 /* We initialize box_rel after reading the inputrec */
3121 clear_mat(state->box_rel);
3123 if (file_version >= 28)
3125 gmx_fio_ndo_rvec(fio, state->boxv, DIM);
3126 if (file_version < 56)
3129 gmx_fio_ndo_rvec(fio, mdum, DIM);
3134 if (state->ngtc > 0 && file_version >= 28)
3137 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
3138 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
3139 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
3140 snew(dumv, state->ngtc);
3141 if (file_version < 69)
3143 bDum = gmx_fio_ndo_real(fio, dumv, state->ngtc);
3145 /* These used to be the Berendsen tcoupl_lambda's */
3146 bDum = gmx_fio_ndo_real(fio, dumv, state->ngtc);
3150 /* Prior to tpx version 26, the inputrec was here.
3151 * I moved it to enable partial forward-compatibility
3152 * for analysis/viewer programs.
3154 if (file_version < 26)
3156 do_test(fio, tpx.bIr, ir);
3157 do_section(fio, eitemIR, bRead);
3162 do_inputrec(fio, ir, bRead, file_version,
3163 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3166 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3171 do_inputrec(fio, &dum_ir, bRead, file_version,
3172 mtop ? &mtop->ffparams.fudgeQQ : NULL);
3175 pr_inputrec(debug, 0, "inputrec", &dum_ir, FALSE);
3177 done_inputrec(&dum_ir);
3183 do_test(fio, tpx.bTop, mtop);
3184 do_section(fio, eitemTOP, bRead);
3187 int mtop_file_version = file_version;
3188 /*allow reading of Gromacs 4.6 files*/
3189 if (mtop_file_version > 80 && mtop_file_version < 90)
3191 mtop_file_version = 79;
3195 do_mtop(fio, mtop, bRead, mtop_file_version);
3199 do_mtop(fio, &dum_top, bRead, mtop_file_version);
3200 done_mtop(&dum_top, TRUE);
3203 do_test(fio, tpx.bX, state->x);
3204 do_section(fio, eitemX, bRead);
3209 state->flags |= (1<<estX);
3211 gmx_fio_ndo_rvec(fio, state->x, state->natoms);
3214 do_test(fio, tpx.bV, state->v);
3215 do_section(fio, eitemV, bRead);
3220 state->flags |= (1<<estV);
3222 gmx_fio_ndo_rvec(fio, state->v, state->natoms);
3225 do_test(fio, tpx.bF, f);
3226 do_section(fio, eitemF, bRead);
3229 gmx_fio_ndo_rvec(fio, f, state->natoms);
3232 /* Starting with tpx version 26, we have the inputrec
3233 * at the end of the file, so we can ignore it
3234 * if the file is never than the software (but still the
3235 * same generation - see comments at the top of this file.
3240 bPeriodicMols = FALSE;
3241 if (file_version >= 26)
3243 do_test(fio, tpx.bIr, ir);
3244 do_section(fio, eitemIR, bRead);
3247 if (file_version >= 53)
3249 /* Removed the pbc info from do_inputrec, since we always want it */
3253 bPeriodicMols = ir->bPeriodicMols;
3255 gmx_fio_do_int(fio, ePBC);
3256 gmx_fio_do_gmx_bool(fio, bPeriodicMols);
3258 if (file_generation <= tpx_generation && ir)
3260 do_inputrec(fio, ir, bRead, file_version, mtop ? &mtop->ffparams.fudgeQQ : NULL);
3263 pr_inputrec(debug, 0, "inputrec", ir, FALSE);
3265 if (file_version < 51)
3267 set_box_rel(ir, state);
3269 if (file_version < 53)
3272 bPeriodicMols = ir->bPeriodicMols;
3275 if (bRead && ir && file_version >= 53)
3277 /* We need to do this after do_inputrec, since that initializes ir */
3279 ir->bPeriodicMols = bPeriodicMols;
3288 if (state->ngtc == 0)
3290 /* Reading old version without tcoupl state data: set it */
3291 init_gtc_state(state, ir->opts.ngtc, 0, ir->opts.nhchainlength);
3293 if (tpx.bTop && mtop)
3295 if (file_version < 57)
3297 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
3299 ir->eDisre = edrSimple;
3303 ir->eDisre = edrNone;
3306 set_disres_npair(mtop);
3310 if (tpx.bTop && mtop)
3312 gmx_mtop_finalize(mtop);
3315 if (file_version >= 57)
3319 env = getenv("GMX_NOCHARGEGROUPS");
3322 sscanf(env, "%d", &ienv);
3323 fprintf(stderr, "\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
3328 "Will make single atomic charge groups in non-solvent%s\n",
3329 ienv > 1 ? " and solvent" : "");
3330 gmx_mtop_make_atomic_charge_groups(mtop, ienv == 1);
3332 fprintf(stderr, "\n");
3340 /************************************************************
3342 * The following routines are the exported ones
3344 ************************************************************/
3346 t_fileio *open_tpx(const char *fn, const char *mode)
3348 return gmx_fio_open(fn, mode);
3351 void close_tpx(t_fileio *fio)
3356 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
3357 int *file_version, int *file_generation)
3361 fio = open_tpx(fn, "r");
3362 do_tpxheader(fio, TRUE, tpx, TopOnlyOK, file_version, file_generation);
3366 void write_tpx_state(const char *fn,
3367 t_inputrec *ir, t_state *state, gmx_mtop_t *mtop)
3371 fio = open_tpx(fn, "w");
3372 do_tpx(fio, FALSE, ir, state, NULL, mtop, FALSE);
3376 void read_tpx_state(const char *fn,
3377 t_inputrec *ir, t_state *state, rvec *f, gmx_mtop_t *mtop)
3381 fio = open_tpx(fn, "r");
3382 do_tpx(fio, TRUE, ir, state, f, mtop, FALSE);
3386 int read_tpx(const char *fn,
3387 t_inputrec *ir, matrix box, int *natoms,
3388 rvec *x, rvec *v, rvec *f, gmx_mtop_t *mtop)
3396 fio = open_tpx(fn, "r");
3397 ePBC = do_tpx(fio, TRUE, ir, &state, f, mtop, TRUE);
3399 *natoms = state.natoms;
3402 copy_mat(state.box, box);
3411 int read_tpx_top(const char *fn,
3412 t_inputrec *ir, matrix box, int *natoms,
3413 rvec *x, rvec *v, rvec *f, t_topology *top)
3419 ePBC = read_tpx(fn, ir, box, natoms, x, v, f, &mtop);
3421 *top = gmx_mtop_t_to_t_topology(&mtop);
3426 gmx_bool fn2bTPX(const char *file)
3428 switch (fn2ftp(file))
3439 gmx_bool read_tps_conf(const char *infile, char *title, t_topology *top, int *ePBC,
3440 rvec **x, rvec **v, matrix box, gmx_bool bMass)
3443 int natoms, i, version, generation;
3444 gmx_bool bTop, bXNULL = FALSE;
3446 t_topology *topconv;
3449 bTop = fn2bTPX(infile);
3453 read_tpxheader(infile, &header, TRUE, &version, &generation);
3456 snew(*x, header.natoms);
3460 snew(*v, header.natoms);
3463 *ePBC = read_tpx(infile, NULL, box, &natoms,
3464 (x == NULL) ? NULL : *x, (v == NULL) ? NULL : *v, NULL, mtop);
3465 *top = gmx_mtop_t_to_t_topology(mtop);
3467 strcpy(title, *top->name);
3468 tpx_make_chain_identifiers(&top->atoms, &top->mols);
3472 get_stx_coordnum(infile, &natoms);
3473 init_t_atoms(&top->atoms, natoms, (fn2ftp(infile) == efPDB));
3484 read_stx_conf(infile, title, &top->atoms, *x, (v == NULL) ? NULL : *v, ePBC, box);
3492 aps = gmx_atomprop_init();
3493 for (i = 0; (i < natoms); i++)
3495 if (!gmx_atomprop_query(aps, epropMass,
3496 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3497 *top->atoms.atomname[i],
3498 &(top->atoms.atom[i].m)))
3502 fprintf(debug, "Can not find mass for atom %s %d %s, setting to 1\n",
3503 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
3504 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
3505 *top->atoms.atomname[i]);
3509 gmx_atomprop_destroy(aps);
3511 top->idef.ntypes = -1;