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 /* This number should be increased whenever the file format changes! */
67 static const int tpx_version = 74;
69 /* This number should only be increased when you edit the TOPOLOGY section
70 * of the tpx format. This way we can maintain forward compatibility too
71 * for all analysis tools and/or external programs that only need to
72 * know the atom/residue names, charges, and bond connectivity.
74 * It first appeared in tpx version 26, when I also moved the inputrecord
75 * to the end of the tpx file, so we can just skip it if we only
78 static const int tpx_generation = 23;
80 /* This number should be the most recent backwards incompatible version
81 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
83 static const int tpx_incompatible_version = 9;
87 /* Struct used to maintain tpx compatibility when function types are added */
89 int fvnr; /* file version number in which the function type first appeared */
90 int ftype; /* function type */
94 *The entries should be ordered in:
95 * 1. ascending file version number
96 * 2. ascending function type number
98 /*static const t_ftupd ftupd[] = {
103 { 22, F_DISRESVIOL },
109 { 26, F_DIHRESVIOL },
110 { 30, F_CROSS_BOND_BONDS },
111 { 30, F_CROSS_BOND_ANGLES },
112 { 30, F_UREY_BRADLEY },
113 { 30, F_POLARIZATION }
117 *The entries should be ordered in:
118 * 1. ascending function type number
119 * 2. ascending file version number
121 static const t_ftupd ftupd[] = {
122 { 20, F_CUBICBONDS },
127 { 43, F_TABBONDSNC },
128 { 70, F_RESTRBONDS },
129 { 30, F_CROSS_BOND_BONDS },
130 { 30, F_CROSS_BOND_ANGLES },
131 { 30, F_UREY_BRADLEY },
132 { 34, F_QUARTIC_ANGLES },
142 { 72, F_NPSOLVATION },
144 { 41, F_LJC_PAIRS_NB },
147 { 32, F_COUL_RECIP },
149 { 30, F_POLARIZATION },
151 { 22, F_DISRESVIOL },
155 { 26, F_DIHRESVIOL },
160 { 46, F_ECONSERVED },
165 #define NFTUPD asize(ftupd)
167 /* Needed for backward compatibility */
170 static void _do_section(t_fileio *fio,int key,gmx_bool bRead,const char *src,
176 if (gmx_fio_getftp(fio) == efTPA) {
178 gmx_fio_write_string(fio,itemstr[key]);
179 bDbg = gmx_fio_getdebug(fio);
180 gmx_fio_setdebug(fio,FALSE);
181 gmx_fio_write_string(fio,comment_str[key]);
182 gmx_fio_setdebug(fio,bDbg);
185 if (gmx_fio_getdebug(fio))
186 fprintf(stderr,"Looking for section %s (%s, %d)",
187 itemstr[key],src,line);
190 gmx_fio_do_string(fio,buf);
191 } while ((gmx_strcasecmp(buf,itemstr[key]) != 0));
193 if (gmx_strcasecmp(buf,itemstr[key]) != 0)
194 gmx_fatal(FARGS,"\nCould not find section heading %s",itemstr[key]);
195 else if (gmx_fio_getdebug(fio))
196 fprintf(stderr," and found it\n");
201 #define do_section(fio,key,bRead) _do_section(fio,key,bRead,__FILE__,__LINE__)
203 /**************************************************************
205 * Now the higer level routines that do io of the structures and arrays
207 **************************************************************/
208 static void do_pullgrp(t_fileio *fio, t_pullgrp *pgrp, gmx_bool bRead,
214 gmx_fio_do_int(fio,pgrp->nat);
216 snew(pgrp->ind,pgrp->nat);
217 bDum=gmx_fio_ndo_int(fio,pgrp->ind,pgrp->nat);
218 gmx_fio_do_int(fio,pgrp->nweight);
220 snew(pgrp->weight,pgrp->nweight);
221 bDum=gmx_fio_ndo_real(fio,pgrp->weight,pgrp->nweight);
222 gmx_fio_do_int(fio,pgrp->pbcatom);
223 gmx_fio_do_rvec(fio,pgrp->vec);
224 gmx_fio_do_rvec(fio,pgrp->init);
225 gmx_fio_do_real(fio,pgrp->rate);
226 gmx_fio_do_real(fio,pgrp->k);
227 if (file_version >= 56) {
228 gmx_fio_do_real(fio,pgrp->kB);
234 static void do_pull(t_fileio *fio, t_pull *pull,gmx_bool bRead, int file_version)
238 gmx_fio_do_int(fio,pull->ngrp);
239 gmx_fio_do_int(fio,pull->eGeom);
240 gmx_fio_do_ivec(fio,pull->dim);
241 gmx_fio_do_real(fio,pull->cyl_r1);
242 gmx_fio_do_real(fio,pull->cyl_r0);
243 gmx_fio_do_real(fio,pull->constr_tol);
244 gmx_fio_do_int(fio,pull->nstxout);
245 gmx_fio_do_int(fio,pull->nstfout);
247 snew(pull->grp,pull->ngrp+1);
248 for(g=0; g<pull->ngrp+1; g++)
249 do_pullgrp(fio,&pull->grp[g],bRead,file_version);
253 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg,gmx_bool bRead, int file_version)
258 gmx_fio_do_int(fio,rotg->eType);
259 gmx_fio_do_int(fio,rotg->bMassW);
260 gmx_fio_do_int(fio,rotg->nat);
262 snew(rotg->ind,rotg->nat);
263 gmx_fio_ndo_int(fio,rotg->ind,rotg->nat);
265 snew(rotg->x_ref,rotg->nat);
266 gmx_fio_ndo_rvec(fio,rotg->x_ref,rotg->nat);
267 gmx_fio_do_rvec(fio,rotg->vec);
268 gmx_fio_do_rvec(fio,rotg->pivot);
269 gmx_fio_do_real(fio,rotg->rate);
270 gmx_fio_do_real(fio,rotg->k);
271 gmx_fio_do_real(fio,rotg->slab_dist);
272 gmx_fio_do_real(fio,rotg->min_gaussian);
273 gmx_fio_do_real(fio,rotg->eps);
274 gmx_fio_do_int(fio,rotg->eFittype);
277 static void do_rot(t_fileio *fio, t_rot *rot,gmx_bool bRead, int file_version)
281 gmx_fio_do_int(fio,rot->ngrp);
282 gmx_fio_do_int(fio,rot->nstrout);
283 gmx_fio_do_int(fio,rot->nstsout);
285 snew(rot->grp,rot->ngrp);
286 for(g=0; g<rot->ngrp; g++)
287 do_rotgrp(fio, &rot->grp[g],bRead,file_version);
291 static void do_inputrec(t_fileio *fio, t_inputrec *ir,gmx_bool bRead,
292 int file_version, real *fudgeQQ)
294 int i,j,k,*tmp,idum=0;
299 real zerotemptime,finish_t,init_temp,finish_temp;
301 if (file_version != tpx_version)
303 /* Give a warning about features that are not accessible */
304 fprintf(stderr,"Note: tpx file_version %d, software version %d\n",
305 file_version,tpx_version);
313 if (file_version == 0)
318 /* Basic inputrec stuff */
319 gmx_fio_do_int(fio,ir->eI);
320 if (file_version >= 62) {
321 gmx_fio_do_gmx_large_int(fio, ir->nsteps);
323 gmx_fio_do_int(fio,idum);
326 if(file_version > 25) {
327 if (file_version >= 62) {
328 gmx_fio_do_gmx_large_int(fio, ir->init_step);
330 gmx_fio_do_int(fio,idum);
331 ir->init_step = idum;
337 if(file_version >= 58)
338 gmx_fio_do_int(fio,ir->simulation_part);
340 ir->simulation_part=1;
342 if (file_version >= 67) {
343 gmx_fio_do_int(fio,ir->nstcalcenergy);
345 ir->nstcalcenergy = 1;
347 if (file_version < 53) {
348 /* The pbc info has been moved out of do_inputrec,
349 * since we always want it, also without reading the inputrec.
351 gmx_fio_do_int(fio,ir->ePBC);
352 if ((file_version <= 15) && (ir->ePBC == 2))
354 if (file_version >= 45) {
355 gmx_fio_do_int(fio,ir->bPeriodicMols);
359 ir->bPeriodicMols = TRUE;
361 ir->bPeriodicMols = FALSE;
365 gmx_fio_do_int(fio,ir->ns_type);
366 gmx_fio_do_int(fio,ir->nstlist);
367 gmx_fio_do_int(fio,ir->ndelta);
368 if (file_version < 41) {
369 gmx_fio_do_int(fio,idum);
370 gmx_fio_do_int(fio,idum);
372 if (file_version >= 45)
373 gmx_fio_do_real(fio,ir->rtpi);
376 gmx_fio_do_int(fio,ir->nstcomm);
377 if (file_version > 34)
378 gmx_fio_do_int(fio,ir->comm_mode);
379 else if (ir->nstcomm < 0)
380 ir->comm_mode = ecmANGULAR;
382 ir->comm_mode = ecmLINEAR;
383 ir->nstcomm = abs(ir->nstcomm);
385 if(file_version > 25)
386 gmx_fio_do_int(fio,ir->nstcheckpoint);
390 gmx_fio_do_int(fio,ir->nstcgsteep);
393 gmx_fio_do_int(fio,ir->nbfgscorr);
397 gmx_fio_do_int(fio,ir->nstlog);
398 gmx_fio_do_int(fio,ir->nstxout);
399 gmx_fio_do_int(fio,ir->nstvout);
400 gmx_fio_do_int(fio,ir->nstfout);
401 gmx_fio_do_int(fio,ir->nstenergy);
402 gmx_fio_do_int(fio,ir->nstxtcout);
403 if (file_version >= 59) {
404 gmx_fio_do_double(fio,ir->init_t);
405 gmx_fio_do_double(fio,ir->delta_t);
407 gmx_fio_do_real(fio,rdum);
409 gmx_fio_do_real(fio,rdum);
412 gmx_fio_do_real(fio,ir->xtcprec);
413 if (file_version < 19) {
414 gmx_fio_do_int(fio,idum);
415 gmx_fio_do_int(fio,idum);
417 if(file_version < 18)
418 gmx_fio_do_int(fio,idum);
419 gmx_fio_do_real(fio,ir->rlist);
420 if (file_version >= 67) {
421 gmx_fio_do_real(fio,ir->rlistlong);
423 gmx_fio_do_int(fio,ir->coulombtype);
424 if (file_version < 32 && ir->coulombtype == eelRF)
425 ir->coulombtype = eelRF_NEC;
426 gmx_fio_do_real(fio,ir->rcoulomb_switch);
427 gmx_fio_do_real(fio,ir->rcoulomb);
428 gmx_fio_do_int(fio,ir->vdwtype);
429 gmx_fio_do_real(fio,ir->rvdw_switch);
430 gmx_fio_do_real(fio,ir->rvdw);
431 if (file_version < 67) {
432 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
434 gmx_fio_do_int(fio,ir->eDispCorr);
435 gmx_fio_do_real(fio,ir->epsilon_r);
436 if (file_version >= 37) {
437 gmx_fio_do_real(fio,ir->epsilon_rf);
439 if (EEL_RF(ir->coulombtype)) {
440 ir->epsilon_rf = ir->epsilon_r;
443 ir->epsilon_rf = 1.0;
446 if (file_version >= 29)
447 gmx_fio_do_real(fio,ir->tabext);
451 if(file_version > 25) {
452 gmx_fio_do_int(fio,ir->gb_algorithm);
453 gmx_fio_do_int(fio,ir->nstgbradii);
454 gmx_fio_do_real(fio,ir->rgbradii);
455 gmx_fio_do_real(fio,ir->gb_saltconc);
456 gmx_fio_do_int(fio,ir->implicit_solvent);
458 ir->gb_algorithm=egbSTILL;
462 ir->implicit_solvent=eisNO;
466 gmx_fio_do_real(fio,ir->gb_epsilon_solvent);
467 gmx_fio_do_real(fio,ir->gb_obc_alpha);
468 gmx_fio_do_real(fio,ir->gb_obc_beta);
469 gmx_fio_do_real(fio,ir->gb_obc_gamma);
472 gmx_fio_do_real(fio,ir->gb_dielectric_offset);
473 gmx_fio_do_int(fio,ir->sa_algorithm);
477 ir->gb_dielectric_offset = 0.009;
478 ir->sa_algorithm = esaAPPROX;
480 gmx_fio_do_real(fio,ir->sa_surface_tension);
482 /* Override sa_surface_tension if it is not changed in the mpd-file */
483 if(ir->sa_surface_tension<0)
485 if(ir->gb_algorithm==egbSTILL)
487 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
489 else if(ir->gb_algorithm==egbHCT || ir->gb_algorithm==egbOBC)
491 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
498 /* Better use sensible values than insane (0.0) ones... */
499 ir->gb_epsilon_solvent = 80;
500 ir->gb_obc_alpha = 1.0;
501 ir->gb_obc_beta = 0.8;
502 ir->gb_obc_gamma = 4.85;
503 ir->sa_surface_tension = 2.092;
507 gmx_fio_do_int(fio,ir->nkx);
508 gmx_fio_do_int(fio,ir->nky);
509 gmx_fio_do_int(fio,ir->nkz);
510 gmx_fio_do_int(fio,ir->pme_order);
511 gmx_fio_do_real(fio,ir->ewald_rtol);
513 if (file_version >=24)
514 gmx_fio_do_int(fio,ir->ewald_geometry);
516 ir->ewald_geometry=eewg3D;
518 if (file_version <=17) {
519 ir->epsilon_surface=0;
520 if (file_version==17)
521 gmx_fio_do_int(fio,idum);
524 gmx_fio_do_real(fio,ir->epsilon_surface);
526 gmx_fio_do_gmx_bool(fio,ir->bOptFFT);
528 gmx_fio_do_gmx_bool(fio,ir->bContinuation);
529 gmx_fio_do_int(fio,ir->etc);
530 /* before version 18, ir->etc was a gmx_bool (ir->btc),
531 * but the values 0 and 1 still mean no and
532 * berendsen temperature coupling, respectively.
534 if (file_version >= 71)
536 gmx_fio_do_int(fio,ir->nsttcouple);
540 ir->nsttcouple = ir->nstcalcenergy;
542 if (file_version <= 15)
544 gmx_fio_do_int(fio,idum);
546 if (file_version <=17)
548 gmx_fio_do_int(fio,ir->epct);
549 if (file_version <= 15)
553 ir->epct = epctSURFACETENSION;
555 gmx_fio_do_int(fio,idum);
558 /* we have removed the NO alternative at the beginning */
562 ir->epct=epctISOTROPIC;
566 ir->epc=epcBERENDSEN;
571 gmx_fio_do_int(fio,ir->epc);
572 gmx_fio_do_int(fio,ir->epct);
574 if (file_version >= 71)
576 gmx_fio_do_int(fio,ir->nstpcouple);
580 ir->nstpcouple = ir->nstcalcenergy;
582 gmx_fio_do_real(fio,ir->tau_p);
583 if (file_version <= 15) {
584 gmx_fio_do_rvec(fio,vdum);
585 clear_mat(ir->ref_p);
587 ir->ref_p[i][i] = vdum[i];
589 gmx_fio_do_rvec(fio,ir->ref_p[XX]);
590 gmx_fio_do_rvec(fio,ir->ref_p[YY]);
591 gmx_fio_do_rvec(fio,ir->ref_p[ZZ]);
593 if (file_version <= 15) {
594 gmx_fio_do_rvec(fio,vdum);
595 clear_mat(ir->compress);
597 ir->compress[i][i] = vdum[i];
600 gmx_fio_do_rvec(fio,ir->compress[XX]);
601 gmx_fio_do_rvec(fio,ir->compress[YY]);
602 gmx_fio_do_rvec(fio,ir->compress[ZZ]);
604 if (file_version >= 47) {
605 gmx_fio_do_int(fio,ir->refcoord_scaling);
606 gmx_fio_do_rvec(fio,ir->posres_com);
607 gmx_fio_do_rvec(fio,ir->posres_comB);
609 ir->refcoord_scaling = erscNO;
610 clear_rvec(ir->posres_com);
611 clear_rvec(ir->posres_comB);
613 if(file_version > 25)
614 gmx_fio_do_int(fio,ir->andersen_seed);
618 if(file_version < 26) {
619 gmx_fio_do_gmx_bool(fio,bSimAnn);
620 gmx_fio_do_real(fio,zerotemptime);
623 if (file_version < 37)
624 gmx_fio_do_real(fio,rdum);
626 gmx_fio_do_real(fio,ir->shake_tol);
627 if (file_version < 54)
628 gmx_fio_do_real(fio,*fudgeQQ);
629 gmx_fio_do_int(fio,ir->efep);
630 if (file_version <= 14 && ir->efep > efepNO)
632 if (file_version >= 59) {
633 gmx_fio_do_double(fio, ir->init_lambda);
634 gmx_fio_do_double(fio, ir->delta_lambda);
636 gmx_fio_do_real(fio,rdum);
637 ir->init_lambda = rdum;
638 gmx_fio_do_real(fio,rdum);
639 ir->delta_lambda = rdum;
641 if (file_version >= 64) {
642 gmx_fio_do_int(fio,ir->n_flambda);
644 snew(ir->flambda,ir->n_flambda);
646 bDum=gmx_fio_ndo_double(fio,ir->flambda,ir->n_flambda);
651 if (file_version >= 13)
652 gmx_fio_do_real(fio,ir->sc_alpha);
655 if (file_version >= 38)
656 gmx_fio_do_int(fio,ir->sc_power);
659 if (file_version >= 15)
660 gmx_fio_do_real(fio,ir->sc_sigma);
665 if (file_version >= 71)
667 ir->sc_sigma_min = ir->sc_sigma;
671 ir->sc_sigma_min = 0;
674 if (file_version >= 64) {
675 gmx_fio_do_int(fio,ir->nstdhdl);
680 if (file_version >= 73)
682 gmx_fio_do_int(fio, ir->separate_dhdl_file);
683 gmx_fio_do_int(fio, ir->dhdl_derivatives);
687 ir->separate_dhdl_file = sepdhdlfileYES;
688 ir->dhdl_derivatives = dhdlderivativesYES;
691 if (file_version >= 71)
693 gmx_fio_do_int(fio,ir->dh_hist_size);
694 gmx_fio_do_double(fio,ir->dh_hist_spacing);
698 ir->dh_hist_size = 0;
699 ir->dh_hist_spacing = 0.1;
701 if (file_version >= 57) {
702 gmx_fio_do_int(fio,ir->eDisre);
704 gmx_fio_do_int(fio,ir->eDisreWeighting);
705 if (file_version < 22) {
706 if (ir->eDisreWeighting == 0)
707 ir->eDisreWeighting = edrwEqual;
709 ir->eDisreWeighting = edrwConservative;
711 gmx_fio_do_gmx_bool(fio,ir->bDisreMixed);
712 gmx_fio_do_real(fio,ir->dr_fc);
713 gmx_fio_do_real(fio,ir->dr_tau);
714 gmx_fio_do_int(fio,ir->nstdisreout);
715 if (file_version >= 22) {
716 gmx_fio_do_real(fio,ir->orires_fc);
717 gmx_fio_do_real(fio,ir->orires_tau);
718 gmx_fio_do_int(fio,ir->nstorireout);
724 if(file_version >= 26) {
725 gmx_fio_do_real(fio,ir->dihre_fc);
726 if (file_version < 56) {
727 gmx_fio_do_real(fio,rdum);
728 gmx_fio_do_int(fio,idum);
734 gmx_fio_do_real(fio,ir->em_stepsize);
735 gmx_fio_do_real(fio,ir->em_tol);
736 if (file_version >= 22)
737 gmx_fio_do_gmx_bool(fio,ir->bShakeSOR);
739 ir->bShakeSOR = TRUE;
740 if (file_version >= 11)
741 gmx_fio_do_int(fio,ir->niter);
744 fprintf(stderr,"Note: niter not in run input file, setting it to %d\n",
747 if (file_version >= 21)
748 gmx_fio_do_real(fio,ir->fc_stepsize);
751 gmx_fio_do_int(fio,ir->eConstrAlg);
752 gmx_fio_do_int(fio,ir->nProjOrder);
753 gmx_fio_do_real(fio,ir->LincsWarnAngle);
754 if (file_version <= 14)
755 gmx_fio_do_int(fio,idum);
756 if (file_version >=26)
757 gmx_fio_do_int(fio,ir->nLincsIter);
760 fprintf(stderr,"Note: nLincsIter not in run input file, setting it to %d\n",
763 if (file_version < 33)
764 gmx_fio_do_real(fio,bd_temp);
765 gmx_fio_do_real(fio,ir->bd_fric);
766 gmx_fio_do_int(fio,ir->ld_seed);
767 if (file_version >= 33) {
769 gmx_fio_do_rvec(fio,ir->deform[i]);
772 clear_rvec(ir->deform[i]);
774 if (file_version >= 14)
775 gmx_fio_do_real(fio,ir->cos_accel);
778 gmx_fio_do_int(fio,ir->userint1);
779 gmx_fio_do_int(fio,ir->userint2);
780 gmx_fio_do_int(fio,ir->userint3);
781 gmx_fio_do_int(fio,ir->userint4);
782 gmx_fio_do_real(fio,ir->userreal1);
783 gmx_fio_do_real(fio,ir->userreal2);
784 gmx_fio_do_real(fio,ir->userreal3);
785 gmx_fio_do_real(fio,ir->userreal4);
788 if (file_version >= 48) {
789 gmx_fio_do_int(fio,ir->ePull);
790 if (ir->ePull != epullNO) {
793 do_pull(fio, ir->pull,bRead,file_version);
799 /* Enforced rotation */
800 if (file_version >= 74) {
801 gmx_fio_do_int(fio,ir->bRot);
802 if (ir->bRot == TRUE) {
805 do_rot(fio, ir->rot,bRead,file_version);
812 gmx_fio_do_int(fio,ir->opts.ngtc);
813 if (file_version >= 69) {
814 gmx_fio_do_int(fio,ir->opts.nhchainlength);
816 ir->opts.nhchainlength = 1;
818 gmx_fio_do_int(fio,ir->opts.ngacc);
819 gmx_fio_do_int(fio,ir->opts.ngfrz);
820 gmx_fio_do_int(fio,ir->opts.ngener);
823 snew(ir->opts.nrdf, ir->opts.ngtc);
824 snew(ir->opts.ref_t, ir->opts.ngtc);
825 snew(ir->opts.annealing, ir->opts.ngtc);
826 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
827 snew(ir->opts.anneal_time, ir->opts.ngtc);
828 snew(ir->opts.anneal_temp, ir->opts.ngtc);
829 snew(ir->opts.tau_t, ir->opts.ngtc);
830 snew(ir->opts.nFreeze,ir->opts.ngfrz);
831 snew(ir->opts.acc, ir->opts.ngacc);
832 snew(ir->opts.egp_flags,ir->opts.ngener*ir->opts.ngener);
834 if (ir->opts.ngtc > 0) {
835 if (bRead && file_version<13) {
836 snew(tmp,ir->opts.ngtc);
837 bDum=gmx_fio_ndo_int(fio,tmp, ir->opts.ngtc);
838 for(i=0; i<ir->opts.ngtc; i++)
839 ir->opts.nrdf[i] = tmp[i];
842 bDum=gmx_fio_ndo_real(fio,ir->opts.nrdf, ir->opts.ngtc);
844 bDum=gmx_fio_ndo_real(fio,ir->opts.ref_t,ir->opts.ngtc);
845 bDum=gmx_fio_ndo_real(fio,ir->opts.tau_t,ir->opts.ngtc);
846 if (file_version<33 && ir->eI==eiBD) {
847 for(i=0; i<ir->opts.ngtc; i++)
848 ir->opts.tau_t[i] = bd_temp;
851 if (ir->opts.ngfrz > 0)
852 bDum=gmx_fio_ndo_ivec(fio,ir->opts.nFreeze,ir->opts.ngfrz);
853 if (ir->opts.ngacc > 0)
854 gmx_fio_ndo_rvec(fio,ir->opts.acc,ir->opts.ngacc);
855 if (file_version >= 12)
856 bDum=gmx_fio_ndo_int(fio,ir->opts.egp_flags,
857 ir->opts.ngener*ir->opts.ngener);
859 if(bRead && file_version < 26) {
860 for(i=0;i<ir->opts.ngtc;i++) {
862 ir->opts.annealing[i] = eannSINGLE;
863 ir->opts.anneal_npoints[i] = 2;
864 snew(ir->opts.anneal_time[i],2);
865 snew(ir->opts.anneal_temp[i],2);
866 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
867 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
868 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
869 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
870 ir->opts.anneal_time[i][0] = ir->init_t;
871 ir->opts.anneal_time[i][1] = finish_t;
872 ir->opts.anneal_temp[i][0] = init_temp;
873 ir->opts.anneal_temp[i][1] = finish_temp;
875 ir->opts.annealing[i] = eannNO;
876 ir->opts.anneal_npoints[i] = 0;
880 /* file version 26 or later */
881 /* First read the lists with annealing and npoints for each group */
882 bDum=gmx_fio_ndo_int(fio,ir->opts.annealing,ir->opts.ngtc);
883 bDum=gmx_fio_ndo_int(fio,ir->opts.anneal_npoints,ir->opts.ngtc);
884 for(j=0;j<(ir->opts.ngtc);j++) {
885 k=ir->opts.anneal_npoints[j];
887 snew(ir->opts.anneal_time[j],k);
888 snew(ir->opts.anneal_temp[j],k);
890 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_time[j],k);
891 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_temp[j],k);
895 if (file_version >= 45) {
896 gmx_fio_do_int(fio,ir->nwall);
897 gmx_fio_do_int(fio,ir->wall_type);
898 if (file_version >= 50)
899 gmx_fio_do_real(fio,ir->wall_r_linpot);
901 ir->wall_r_linpot = -1;
902 gmx_fio_do_int(fio,ir->wall_atomtype[0]);
903 gmx_fio_do_int(fio,ir->wall_atomtype[1]);
904 gmx_fio_do_real(fio,ir->wall_density[0]);
905 gmx_fio_do_real(fio,ir->wall_density[1]);
906 gmx_fio_do_real(fio,ir->wall_ewald_zfac);
910 ir->wall_atomtype[0] = -1;
911 ir->wall_atomtype[1] = -1;
912 ir->wall_density[0] = 0;
913 ir->wall_density[1] = 0;
914 ir->wall_ewald_zfac = 3;
916 /* Cosine stuff for electric fields */
917 for(j=0; (j<DIM); j++) {
918 gmx_fio_do_int(fio,ir->ex[j].n);
919 gmx_fio_do_int(fio,ir->et[j].n);
921 snew(ir->ex[j].a, ir->ex[j].n);
922 snew(ir->ex[j].phi,ir->ex[j].n);
923 snew(ir->et[j].a, ir->et[j].n);
924 snew(ir->et[j].phi,ir->et[j].n);
926 bDum=gmx_fio_ndo_real(fio,ir->ex[j].a, ir->ex[j].n);
927 bDum=gmx_fio_ndo_real(fio,ir->ex[j].phi,ir->ex[j].n);
928 bDum=gmx_fio_ndo_real(fio,ir->et[j].a, ir->et[j].n);
929 bDum=gmx_fio_ndo_real(fio,ir->et[j].phi,ir->et[j].n);
933 if(file_version>=39){
934 gmx_fio_do_gmx_bool(fio,ir->bQMMM);
935 gmx_fio_do_int(fio,ir->QMMMscheme);
936 gmx_fio_do_real(fio,ir->scalefactor);
937 gmx_fio_do_int(fio,ir->opts.ngQM);
939 snew(ir->opts.QMmethod, ir->opts.ngQM);
940 snew(ir->opts.QMbasis, ir->opts.ngQM);
941 snew(ir->opts.QMcharge, ir->opts.ngQM);
942 snew(ir->opts.QMmult, ir->opts.ngQM);
943 snew(ir->opts.bSH, ir->opts.ngQM);
944 snew(ir->opts.CASorbitals, ir->opts.ngQM);
945 snew(ir->opts.CASelectrons,ir->opts.ngQM);
946 snew(ir->opts.SAon, ir->opts.ngQM);
947 snew(ir->opts.SAoff, ir->opts.ngQM);
948 snew(ir->opts.SAsteps, ir->opts.ngQM);
949 snew(ir->opts.bOPT, ir->opts.ngQM);
950 snew(ir->opts.bTS, ir->opts.ngQM);
952 if (ir->opts.ngQM > 0) {
953 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmethod,ir->opts.ngQM);
954 bDum=gmx_fio_ndo_int(fio,ir->opts.QMbasis,ir->opts.ngQM);
955 bDum=gmx_fio_ndo_int(fio,ir->opts.QMcharge,ir->opts.ngQM);
956 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmult,ir->opts.ngQM);
957 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bSH,ir->opts.ngQM);
958 bDum=gmx_fio_ndo_int(fio,ir->opts.CASorbitals,ir->opts.ngQM);
959 bDum=gmx_fio_ndo_int(fio,ir->opts.CASelectrons,ir->opts.ngQM);
960 bDum=gmx_fio_ndo_real(fio,ir->opts.SAon,ir->opts.ngQM);
961 bDum=gmx_fio_ndo_real(fio,ir->opts.SAoff,ir->opts.ngQM);
962 bDum=gmx_fio_ndo_int(fio,ir->opts.SAsteps,ir->opts.ngQM);
963 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bOPT,ir->opts.ngQM);
964 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bTS,ir->opts.ngQM);
966 /* end of QMMM stuff */
971 static void do_harm(t_fileio *fio, t_iparams *iparams,gmx_bool bRead)
973 gmx_fio_do_real(fio,iparams->harmonic.rA);
974 gmx_fio_do_real(fio,iparams->harmonic.krA);
975 gmx_fio_do_real(fio,iparams->harmonic.rB);
976 gmx_fio_do_real(fio,iparams->harmonic.krB);
979 void do_iparams(t_fileio *fio, t_functype ftype,t_iparams *iparams,
980 gmx_bool bRead, int file_version)
987 gmx_fio_set_comment(fio, interaction_function[ftype].name);
995 do_harm(fio, iparams,bRead);
996 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead) {
997 /* Correct incorrect storage of parameters */
998 iparams->pdihs.phiB = iparams->pdihs.phiA;
999 iparams->pdihs.cpB = iparams->pdihs.cpA;
1003 gmx_fio_do_real(fio,iparams->fene.bm);
1004 gmx_fio_do_real(fio,iparams->fene.kb);
1007 gmx_fio_do_real(fio,iparams->restraint.lowA);
1008 gmx_fio_do_real(fio,iparams->restraint.up1A);
1009 gmx_fio_do_real(fio,iparams->restraint.up2A);
1010 gmx_fio_do_real(fio,iparams->restraint.kA);
1011 gmx_fio_do_real(fio,iparams->restraint.lowB);
1012 gmx_fio_do_real(fio,iparams->restraint.up1B);
1013 gmx_fio_do_real(fio,iparams->restraint.up2B);
1014 gmx_fio_do_real(fio,iparams->restraint.kB);
1020 gmx_fio_do_real(fio,iparams->tab.kA);
1021 gmx_fio_do_int(fio,iparams->tab.table);
1022 gmx_fio_do_real(fio,iparams->tab.kB);
1024 case F_CROSS_BOND_BONDS:
1025 gmx_fio_do_real(fio,iparams->cross_bb.r1e);
1026 gmx_fio_do_real(fio,iparams->cross_bb.r2e);
1027 gmx_fio_do_real(fio,iparams->cross_bb.krr);
1029 case F_CROSS_BOND_ANGLES:
1030 gmx_fio_do_real(fio,iparams->cross_ba.r1e);
1031 gmx_fio_do_real(fio,iparams->cross_ba.r2e);
1032 gmx_fio_do_real(fio,iparams->cross_ba.r3e);
1033 gmx_fio_do_real(fio,iparams->cross_ba.krt);
1035 case F_UREY_BRADLEY:
1036 gmx_fio_do_real(fio,iparams->u_b.theta);
1037 gmx_fio_do_real(fio,iparams->u_b.ktheta);
1038 gmx_fio_do_real(fio,iparams->u_b.r13);
1039 gmx_fio_do_real(fio,iparams->u_b.kUB);
1041 case F_QUARTIC_ANGLES:
1042 gmx_fio_do_real(fio,iparams->qangle.theta);
1043 bDum=gmx_fio_ndo_real(fio,iparams->qangle.c,5);
1046 gmx_fio_do_real(fio,iparams->bham.a);
1047 gmx_fio_do_real(fio,iparams->bham.b);
1048 gmx_fio_do_real(fio,iparams->bham.c);
1051 gmx_fio_do_real(fio,iparams->morse.b0);
1052 gmx_fio_do_real(fio,iparams->morse.cb);
1053 gmx_fio_do_real(fio,iparams->morse.beta);
1056 gmx_fio_do_real(fio,iparams->cubic.b0);
1057 gmx_fio_do_real(fio,iparams->cubic.kb);
1058 gmx_fio_do_real(fio,iparams->cubic.kcub);
1062 case F_POLARIZATION:
1063 gmx_fio_do_real(fio,iparams->polarize.alpha);
1066 if (file_version < 31)
1067 gmx_fatal(FARGS,"Old tpr files with water_polarization not supported. Make a new.");
1068 gmx_fio_do_real(fio,iparams->wpol.al_x);
1069 gmx_fio_do_real(fio,iparams->wpol.al_y);
1070 gmx_fio_do_real(fio,iparams->wpol.al_z);
1071 gmx_fio_do_real(fio,iparams->wpol.rOH);
1072 gmx_fio_do_real(fio,iparams->wpol.rHH);
1073 gmx_fio_do_real(fio,iparams->wpol.rOD);
1076 gmx_fio_do_real(fio,iparams->thole.a);
1077 gmx_fio_do_real(fio,iparams->thole.alpha1);
1078 gmx_fio_do_real(fio,iparams->thole.alpha2);
1079 gmx_fio_do_real(fio,iparams->thole.rfac);
1082 gmx_fio_do_real(fio,iparams->lj.c6);
1083 gmx_fio_do_real(fio,iparams->lj.c12);
1086 gmx_fio_do_real(fio,iparams->lj14.c6A);
1087 gmx_fio_do_real(fio,iparams->lj14.c12A);
1088 gmx_fio_do_real(fio,iparams->lj14.c6B);
1089 gmx_fio_do_real(fio,iparams->lj14.c12B);
1092 gmx_fio_do_real(fio,iparams->ljc14.fqq);
1093 gmx_fio_do_real(fio,iparams->ljc14.qi);
1094 gmx_fio_do_real(fio,iparams->ljc14.qj);
1095 gmx_fio_do_real(fio,iparams->ljc14.c6);
1096 gmx_fio_do_real(fio,iparams->ljc14.c12);
1098 case F_LJC_PAIRS_NB:
1099 gmx_fio_do_real(fio,iparams->ljcnb.qi);
1100 gmx_fio_do_real(fio,iparams->ljcnb.qj);
1101 gmx_fio_do_real(fio,iparams->ljcnb.c6);
1102 gmx_fio_do_real(fio,iparams->ljcnb.c12);
1108 gmx_fio_do_real(fio,iparams->pdihs.phiA);
1109 gmx_fio_do_real(fio,iparams->pdihs.cpA);
1110 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42) {
1111 /* Read the incorrectly stored multiplicity */
1112 gmx_fio_do_real(fio,iparams->harmonic.rB);
1113 gmx_fio_do_real(fio,iparams->harmonic.krB);
1114 iparams->pdihs.phiB = iparams->pdihs.phiA;
1115 iparams->pdihs.cpB = iparams->pdihs.cpA;
1117 gmx_fio_do_real(fio,iparams->pdihs.phiB);
1118 gmx_fio_do_real(fio,iparams->pdihs.cpB);
1119 gmx_fio_do_int(fio,iparams->pdihs.mult);
1123 gmx_fio_do_int(fio,iparams->disres.label);
1124 gmx_fio_do_int(fio,iparams->disres.type);
1125 gmx_fio_do_real(fio,iparams->disres.low);
1126 gmx_fio_do_real(fio,iparams->disres.up1);
1127 gmx_fio_do_real(fio,iparams->disres.up2);
1128 gmx_fio_do_real(fio,iparams->disres.kfac);
1131 gmx_fio_do_int(fio,iparams->orires.ex);
1132 gmx_fio_do_int(fio,iparams->orires.label);
1133 gmx_fio_do_int(fio,iparams->orires.power);
1134 gmx_fio_do_real(fio,iparams->orires.c);
1135 gmx_fio_do_real(fio,iparams->orires.obs);
1136 gmx_fio_do_real(fio,iparams->orires.kfac);
1139 gmx_fio_do_int(fio,iparams->dihres.power);
1140 gmx_fio_do_int(fio,iparams->dihres.label);
1141 gmx_fio_do_real(fio,iparams->dihres.phi);
1142 gmx_fio_do_real(fio,iparams->dihres.dphi);
1143 gmx_fio_do_real(fio,iparams->dihres.kfac);
1146 gmx_fio_do_rvec(fio,iparams->posres.pos0A);
1147 gmx_fio_do_rvec(fio,iparams->posres.fcA);
1148 if (bRead && file_version < 27) {
1149 copy_rvec(iparams->posres.pos0A,iparams->posres.pos0B);
1150 copy_rvec(iparams->posres.fcA,iparams->posres.fcB);
1152 gmx_fio_do_rvec(fio,iparams->posres.pos0B);
1153 gmx_fio_do_rvec(fio,iparams->posres.fcB);
1157 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA,NR_RBDIHS);
1158 if(file_version>=25)
1159 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB,NR_RBDIHS);
1162 /* Fourier dihedrals are internally represented
1163 * as Ryckaert-Bellemans since those are faster to compute.
1165 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA, NR_RBDIHS);
1166 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB, NR_RBDIHS);
1170 gmx_fio_do_real(fio,iparams->constr.dA);
1171 gmx_fio_do_real(fio,iparams->constr.dB);
1174 gmx_fio_do_real(fio,iparams->settle.doh);
1175 gmx_fio_do_real(fio,iparams->settle.dhh);
1178 gmx_fio_do_real(fio,iparams->vsite.a);
1183 gmx_fio_do_real(fio,iparams->vsite.a);
1184 gmx_fio_do_real(fio,iparams->vsite.b);
1189 gmx_fio_do_real(fio,iparams->vsite.a);
1190 gmx_fio_do_real(fio,iparams->vsite.b);
1191 gmx_fio_do_real(fio,iparams->vsite.c);
1194 gmx_fio_do_int(fio,iparams->vsiten.n);
1195 gmx_fio_do_real(fio,iparams->vsiten.a);
1200 /* We got rid of some parameters in version 68 */
1201 if(bRead && file_version<68)
1203 gmx_fio_do_real(fio,rdum);
1204 gmx_fio_do_real(fio,rdum);
1205 gmx_fio_do_real(fio,rdum);
1206 gmx_fio_do_real(fio,rdum);
1208 gmx_fio_do_real(fio,iparams->gb.sar);
1209 gmx_fio_do_real(fio,iparams->gb.st);
1210 gmx_fio_do_real(fio,iparams->gb.pi);
1211 gmx_fio_do_real(fio,iparams->gb.gbr);
1212 gmx_fio_do_real(fio,iparams->gb.bmlt);
1215 gmx_fio_do_int(fio,iparams->cmap.cmapA);
1216 gmx_fio_do_int(fio,iparams->cmap.cmapB);
1219 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
1221 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
1224 gmx_fio_unset_comment(fio);
1227 static void do_ilist(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,int file_version,
1234 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1236 if (file_version < 44) {
1237 for(i=0; i<MAXNODES; i++)
1238 gmx_fio_do_int(fio,idum);
1240 gmx_fio_do_int(fio,ilist->nr);
1242 snew(ilist->iatoms,ilist->nr);
1243 bDum=gmx_fio_ndo_int(fio,ilist->iatoms,ilist->nr);
1245 gmx_fio_unset_comment(fio);
1248 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1249 gmx_bool bRead, int file_version)
1255 gmx_fio_do_int(fio,ffparams->atnr);
1256 if (file_version < 57) {
1257 gmx_fio_do_int(fio,idum);
1259 gmx_fio_do_int(fio,ffparams->ntypes);
1261 fprintf(debug,"ffparams->atnr = %d, ntypes = %d\n",
1262 ffparams->atnr,ffparams->ntypes);
1264 snew(ffparams->functype,ffparams->ntypes);
1265 snew(ffparams->iparams,ffparams->ntypes);
1267 /* Read/write all the function types */
1268 bDum=gmx_fio_ndo_int(fio,ffparams->functype,ffparams->ntypes);
1270 pr_ivec(debug,0,"functype",ffparams->functype,ffparams->ntypes,TRUE);
1272 if (file_version >= 66) {
1273 gmx_fio_do_double(fio,ffparams->reppow);
1275 ffparams->reppow = 12.0;
1278 if (file_version >= 57) {
1279 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1282 /* Check whether all these function types are supported by the code.
1283 * In practice the code is backwards compatible, which means that the
1284 * numbering may have to be altered from old numbering to new numbering
1286 for (i=0; (i<ffparams->ntypes); i++) {
1288 /* Loop over file versions */
1289 for (k=0; (k<NFTUPD); k++)
1290 /* Compare the read file_version to the update table */
1291 if ((file_version < ftupd[k].fvnr) &&
1292 (ffparams->functype[i] >= ftupd[k].ftype)) {
1293 ffparams->functype[i] += 1;
1295 fprintf(debug,"Incrementing function type %d to %d (due to %s)\n",
1296 i,ffparams->functype[i],
1297 interaction_function[ftupd[k].ftype].longname);
1302 do_iparams(fio, ffparams->functype[i],&ffparams->iparams[i],bRead,
1305 pr_iparams(debug,ffparams->functype[i],&ffparams->iparams[i]);
1309 static void do_ilists(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,
1312 int i,j,renum[F_NRE];
1313 gmx_bool bDum=TRUE,bClear;
1316 for(j=0; (j<F_NRE); j++) {
1319 for (k=0; k<NFTUPD; k++)
1320 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
1324 ilist[j].iatoms = NULL;
1326 do_ilist(fio, &ilist[j],bRead,file_version,j);
1329 if (bRead && gmx_debug_at)
1330 pr_ilist(debug,0,interaction_function[j].longname,
1331 functype,&ilist[j],TRUE);
1336 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams,gmx_moltype_t *molt,
1337 gmx_bool bRead, int file_version)
1339 do_ffparams(fio, ffparams,bRead,file_version);
1341 if (file_version >= 54) {
1342 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1345 do_ilists(fio, molt->ilist,bRead,file_version);
1348 static void do_block(t_fileio *fio, t_block *block,gmx_bool bRead,int file_version)
1350 int i,idum,dum_nra,*dum_a;
1353 if (file_version < 44)
1354 for(i=0; i<MAXNODES; i++)
1355 gmx_fio_do_int(fio,idum);
1356 gmx_fio_do_int(fio,block->nr);
1357 if (file_version < 51)
1358 gmx_fio_do_int(fio,dum_nra);
1360 block->nalloc_index = block->nr+1;
1361 snew(block->index,block->nalloc_index);
1363 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1365 if (file_version < 51 && dum_nra > 0) {
1366 snew(dum_a,dum_nra);
1367 bDum=gmx_fio_ndo_int(fio,dum_a,dum_nra);
1372 static void do_blocka(t_fileio *fio, t_blocka *block,gmx_bool bRead,
1378 if (file_version < 44)
1379 for(i=0; i<MAXNODES; i++)
1380 gmx_fio_do_int(fio,idum);
1381 gmx_fio_do_int(fio,block->nr);
1382 gmx_fio_do_int(fio,block->nra);
1384 block->nalloc_index = block->nr+1;
1385 snew(block->index,block->nalloc_index);
1386 block->nalloc_a = block->nra;
1387 snew(block->a,block->nalloc_a);
1389 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1390 bDum=gmx_fio_ndo_int(fio,block->a,block->nra);
1393 static void do_atom(t_fileio *fio, t_atom *atom,int ngrp,gmx_bool bRead,
1394 int file_version, gmx_groups_t *groups,int atnr)
1398 gmx_fio_do_real(fio,atom->m);
1399 gmx_fio_do_real(fio,atom->q);
1400 gmx_fio_do_real(fio,atom->mB);
1401 gmx_fio_do_real(fio,atom->qB);
1402 gmx_fio_do_ushort(fio, atom->type);
1403 gmx_fio_do_ushort(fio, atom->typeB);
1404 gmx_fio_do_int(fio,atom->ptype);
1405 gmx_fio_do_int(fio,atom->resind);
1406 if (file_version >= 52)
1407 gmx_fio_do_int(fio,atom->atomnumber);
1409 atom->atomnumber = NOTSET;
1410 if (file_version < 23)
1412 else if (file_version < 39)
1417 if (file_version < 57) {
1418 unsigned char uchar[egcNR];
1419 gmx_fio_ndo_uchar(fio,uchar,myngrp);
1420 for(i=myngrp; (i<ngrp); i++) {
1423 /* Copy the old data format to the groups struct */
1424 for(i=0; i<ngrp; i++) {
1425 groups->grpnr[i][atnr] = uchar[i];
1430 static void do_grps(t_fileio *fio, int ngrp,t_grps grps[],gmx_bool bRead,
1436 if (file_version < 23)
1438 else if (file_version < 39)
1443 for(j=0; (j<ngrp); j++) {
1445 gmx_fio_do_int(fio,grps[j].nr);
1447 snew(grps[j].nm_ind,grps[j].nr);
1448 bDum=gmx_fio_ndo_int(fio,grps[j].nm_ind,grps[j].nr);
1452 snew(grps[j].nm_ind,grps[j].nr);
1457 static void do_symstr(t_fileio *fio, char ***nm,gmx_bool bRead,t_symtab *symtab)
1462 gmx_fio_do_int(fio,ls);
1463 *nm = get_symtab_handle(symtab,ls);
1466 ls = lookup_symtab(symtab,*nm);
1467 gmx_fio_do_int(fio,ls);
1471 static void do_strstr(t_fileio *fio, int nstr,char ***nm,gmx_bool bRead,
1476 for (j=0; (j<nstr); j++)
1477 do_symstr(fio, &(nm[j]),bRead,symtab);
1480 static void do_resinfo(t_fileio *fio, int n,t_resinfo *ri,gmx_bool bRead,
1481 t_symtab *symtab, int file_version)
1485 for (j=0; (j<n); j++) {
1486 do_symstr(fio, &(ri[j].name),bRead,symtab);
1487 if (file_version >= 63) {
1488 gmx_fio_do_int(fio,ri[j].nr);
1489 gmx_fio_do_uchar(fio, ri[j].ic);
1497 static void do_atoms(t_fileio *fio, t_atoms *atoms,gmx_bool bRead,t_symtab *symtab,
1499 gmx_groups_t *groups)
1503 gmx_fio_do_int(fio,atoms->nr);
1504 gmx_fio_do_int(fio,atoms->nres);
1505 if (file_version < 57) {
1506 gmx_fio_do_int(fio,groups->ngrpname);
1507 for(i=0; i<egcNR; i++) {
1508 groups->ngrpnr[i] = atoms->nr;
1509 snew(groups->grpnr[i],groups->ngrpnr[i]);
1513 snew(atoms->atom,atoms->nr);
1514 snew(atoms->atomname,atoms->nr);
1515 snew(atoms->atomtype,atoms->nr);
1516 snew(atoms->atomtypeB,atoms->nr);
1517 snew(atoms->resinfo,atoms->nres);
1518 if (file_version < 57) {
1519 snew(groups->grpname,groups->ngrpname);
1521 atoms->pdbinfo = NULL;
1523 for(i=0; (i<atoms->nr); i++) {
1524 do_atom(fio, &atoms->atom[i],egcNR,bRead, file_version,groups,i);
1526 do_strstr(fio, atoms->nr,atoms->atomname,bRead,symtab);
1527 if (bRead && (file_version <= 20)) {
1528 for(i=0; i<atoms->nr; i++) {
1529 atoms->atomtype[i] = put_symtab(symtab,"?");
1530 atoms->atomtypeB[i] = put_symtab(symtab,"?");
1533 do_strstr(fio, atoms->nr,atoms->atomtype,bRead,symtab);
1534 do_strstr(fio, atoms->nr,atoms->atomtypeB,bRead,symtab);
1536 do_resinfo(fio, atoms->nres,atoms->resinfo,bRead,symtab,file_version);
1538 if (file_version < 57) {
1539 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1541 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1545 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
1546 gmx_bool bRead,t_symtab *symtab,
1552 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1553 gmx_fio_do_int(fio,groups->ngrpname);
1555 snew(groups->grpname,groups->ngrpname);
1557 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1558 for(g=0; g<egcNR; g++) {
1559 gmx_fio_do_int(fio,groups->ngrpnr[g]);
1560 if (groups->ngrpnr[g] == 0) {
1562 groups->grpnr[g] = NULL;
1566 snew(groups->grpnr[g],groups->ngrpnr[g]);
1568 bDum=gmx_fio_ndo_uchar(fio, groups->grpnr[g],groups->ngrpnr[g]);
1573 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes,gmx_bool bRead,
1574 t_symtab *symtab,int file_version)
1577 gmx_bool bDum = TRUE;
1579 if (file_version > 25) {
1580 gmx_fio_do_int(fio,atomtypes->nr);
1583 snew(atomtypes->radius,j);
1584 snew(atomtypes->vol,j);
1585 snew(atomtypes->surftens,j);
1586 snew(atomtypes->atomnumber,j);
1587 snew(atomtypes->gb_radius,j);
1588 snew(atomtypes->S_hct,j);
1590 bDum=gmx_fio_ndo_real(fio,atomtypes->radius,j);
1591 bDum=gmx_fio_ndo_real(fio,atomtypes->vol,j);
1592 bDum=gmx_fio_ndo_real(fio,atomtypes->surftens,j);
1593 if(file_version >= 40)
1595 bDum=gmx_fio_ndo_int(fio,atomtypes->atomnumber,j);
1597 if(file_version >= 60)
1599 bDum=gmx_fio_ndo_real(fio,atomtypes->gb_radius,j);
1600 bDum=gmx_fio_ndo_real(fio,atomtypes->S_hct,j);
1603 /* File versions prior to 26 cannot do GBSA,
1604 * so they dont use this structure
1607 atomtypes->radius = NULL;
1608 atomtypes->vol = NULL;
1609 atomtypes->surftens = NULL;
1610 atomtypes->atomnumber = NULL;
1611 atomtypes->gb_radius = NULL;
1612 atomtypes->S_hct = NULL;
1616 static void do_symtab(t_fileio *fio, t_symtab *symtab,gmx_bool bRead)
1622 gmx_fio_do_int(fio,symtab->nr);
1625 snew(symtab->symbuf,1);
1626 symbuf = symtab->symbuf;
1627 symbuf->bufsize = nr;
1628 snew(symbuf->buf,nr);
1629 for (i=0; (i<nr); i++) {
1630 gmx_fio_do_string(fio,buf);
1631 symbuf->buf[i]=strdup(buf);
1635 symbuf = symtab->symbuf;
1636 while (symbuf!=NULL) {
1637 for (i=0; (i<symbuf->bufsize) && (i<nr); i++)
1638 gmx_fio_do_string(fio,symbuf->buf[i]);
1640 symbuf=symbuf->next;
1643 gmx_fatal(FARGS,"nr of symtab strings left: %d",nr);
1647 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
1649 int i,j,ngrid,gs,nelem;
1651 gmx_fio_do_int(fio,cmap_grid->ngrid);
1652 gmx_fio_do_int(fio,cmap_grid->grid_spacing);
1654 ngrid = cmap_grid->ngrid;
1655 gs = cmap_grid->grid_spacing;
1660 snew(cmap_grid->cmapdata,ngrid);
1662 for(i=0;i<cmap_grid->ngrid;i++)
1664 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
1668 for(i=0;i<cmap_grid->ngrid;i++)
1670 for(j=0;j<nelem;j++)
1672 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4]);
1673 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+1]);
1674 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+2]);
1675 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+3]);
1681 void tpx_make_chain_identifiers(t_atoms *atoms,t_block *mols)
1687 /* We always assign a new chain number, but save the chain id characters
1688 * for larger molecules.
1690 #define CHAIN_MIN_ATOMS 15
1694 for(m=0; m<mols->nr; m++)
1697 a1=mols->index[m+1];
1698 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
1707 for(a=a0; a<a1; a++)
1709 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
1710 atoms->resinfo[atoms->atom[a].resind].chainid = c;
1715 /* Blank out the chain id if there was only one chain */
1718 for(r=0; r<atoms->nres; r++)
1720 atoms->resinfo[r].chainid = ' ';
1725 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt,gmx_bool bRead,
1726 t_symtab *symtab, int file_version,
1727 gmx_groups_t *groups)
1731 if (file_version >= 57) {
1732 do_symstr(fio, &(molt->name),bRead,symtab);
1735 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
1737 if (bRead && gmx_debug_at) {
1738 pr_atoms(debug,0,"atoms",&molt->atoms,TRUE);
1741 if (file_version >= 57) {
1742 do_ilists(fio, molt->ilist,bRead,file_version);
1744 do_block(fio, &molt->cgs,bRead,file_version);
1745 if (bRead && gmx_debug_at) {
1746 pr_block(debug,0,"cgs",&molt->cgs,TRUE);
1750 /* This used to be in the atoms struct */
1751 do_blocka(fio, &molt->excls, bRead, file_version);
1754 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb,gmx_bool bRead,
1759 gmx_fio_do_int(fio,molb->type);
1760 gmx_fio_do_int(fio,molb->nmol);
1761 gmx_fio_do_int(fio,molb->natoms_mol);
1762 /* Position restraint coordinates */
1763 gmx_fio_do_int(fio,molb->nposres_xA);
1764 if (molb->nposres_xA > 0) {
1766 snew(molb->posres_xA,molb->nposres_xA);
1768 gmx_fio_ndo_rvec(fio,molb->posres_xA,molb->nposres_xA);
1770 gmx_fio_do_int(fio,molb->nposres_xB);
1771 if (molb->nposres_xB > 0) {
1773 snew(molb->posres_xB,molb->nposres_xB);
1775 gmx_fio_ndo_rvec(fio,molb->posres_xB,molb->nposres_xB);
1780 static t_block mtop_mols(gmx_mtop_t *mtop)
1786 for(mb=0; mb<mtop->nmolblock; mb++) {
1787 mols.nr += mtop->molblock[mb].nmol;
1789 mols.nalloc_index = mols.nr + 1;
1790 snew(mols.index,mols.nalloc_index);
1795 for(mb=0; mb<mtop->nmolblock; mb++) {
1796 for(mol=0; mol<mtop->molblock[mb].nmol; mol++) {
1797 a += mtop->molblock[mb].natoms_mol;
1806 static void add_posres_molblock(gmx_mtop_t *mtop)
1811 gmx_molblock_t *molb;
1814 il = &mtop->moltype[0].ilist[F_POSRES];
1820 for(i=0; i<il->nr; i+=2) {
1821 ip = &mtop->ffparams.iparams[il->iatoms[i]];
1822 am = max(am,il->iatoms[i+1]);
1823 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
1824 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
1825 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ]) {
1829 /* Make the posres coordinate block end at a molecule end */
1831 while(am >= mtop->mols.index[mol+1]) {
1834 molb = &mtop->molblock[0];
1835 molb->nposres_xA = mtop->mols.index[mol+1];
1836 snew(molb->posres_xA,molb->nposres_xA);
1838 molb->nposres_xB = molb->nposres_xA;
1839 snew(molb->posres_xB,molb->nposres_xB);
1841 molb->nposres_xB = 0;
1843 for(i=0; i<il->nr; i+=2) {
1844 ip = &mtop->ffparams.iparams[il->iatoms[i]];
1845 a = il->iatoms[i+1];
1846 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
1847 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
1848 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
1850 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
1851 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
1852 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
1857 static void set_disres_npair(gmx_mtop_t *mtop)
1864 ip = mtop->ffparams.iparams;
1866 for(mt=0; mt<mtop->nmoltype; mt++) {
1867 il = &mtop->moltype[mt].ilist[F_DISRES];
1871 for(i=0; i<il->nr; i+=3) {
1873 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label) {
1874 ip[a[i]].disres.npair = npair;
1882 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop,gmx_bool bRead,
1890 do_symtab(fio, &(mtop->symtab),bRead);
1892 pr_symtab(debug,0,"symtab",&mtop->symtab);
1894 do_symstr(fio, &(mtop->name),bRead,&(mtop->symtab));
1896 if (file_version >= 57) {
1897 do_ffparams(fio, &mtop->ffparams,bRead,file_version);
1899 gmx_fio_do_int(fio,mtop->nmoltype);
1904 snew(mtop->moltype,mtop->nmoltype);
1905 if (file_version < 57) {
1906 mtop->moltype[0].name = mtop->name;
1909 for(mt=0; mt<mtop->nmoltype; mt++) {
1910 do_moltype(fio, &mtop->moltype[mt],bRead,&mtop->symtab,file_version,
1914 if (file_version >= 57) {
1915 gmx_fio_do_int(fio,mtop->nmolblock);
1917 mtop->nmolblock = 1;
1920 snew(mtop->molblock,mtop->nmolblock);
1922 if (file_version >= 57) {
1923 for(mb=0; mb<mtop->nmolblock; mb++) {
1924 do_molblock(fio, &mtop->molblock[mb],bRead,file_version);
1926 gmx_fio_do_int(fio,mtop->natoms);
1928 mtop->molblock[0].type = 0;
1929 mtop->molblock[0].nmol = 1;
1930 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
1931 mtop->molblock[0].nposres_xA = 0;
1932 mtop->molblock[0].nposres_xB = 0;
1935 do_atomtypes (fio, &(mtop->atomtypes),bRead,&(mtop->symtab), file_version);
1937 pr_atomtypes(debug,0,"atomtypes",&mtop->atomtypes,TRUE);
1939 if (file_version < 57) {
1940 /* Debug statements are inside do_idef */
1941 do_idef (fio, &mtop->ffparams,&mtop->moltype[0],bRead,file_version);
1942 mtop->natoms = mtop->moltype[0].atoms.nr;
1945 if(file_version >= 65)
1947 do_cmap(fio, &mtop->ffparams.cmap_grid,bRead);
1951 mtop->ffparams.cmap_grid.ngrid = 0;
1952 mtop->ffparams.cmap_grid.grid_spacing = 0.1;
1953 mtop->ffparams.cmap_grid.cmapdata = NULL;
1956 if (file_version >= 57) {
1957 do_groups(fio, &mtop->groups,bRead,&(mtop->symtab),file_version);
1960 if (file_version < 57) {
1961 do_block(fio, &mtop->moltype[0].cgs,bRead,file_version);
1962 if (bRead && gmx_debug_at) {
1963 pr_block(debug,0,"cgs",&mtop->moltype[0].cgs,TRUE);
1965 do_block(fio, &mtop->mols,bRead,file_version);
1966 /* Add the posres coordinates to the molblock */
1967 add_posres_molblock(mtop);
1970 if (file_version >= 57) {
1971 mtop->mols = mtop_mols(mtop);
1974 pr_block(debug,0,"mols",&mtop->mols,TRUE);
1978 if (file_version < 51) {
1979 /* Here used to be the shake blocks */
1980 do_blocka(fio, &dumb,bRead,file_version);
1988 close_symtab(&(mtop->symtab));
1992 /* If TopOnlyOK is TRUE then we can read even future versions
1993 * of tpx files, provided the file_generation hasn't changed.
1994 * If it is FALSE, we need the inputrecord too, and bail out
1995 * if the file is newer than the program.
1997 * The version and generation if the topology (see top of this file)
1998 * are returned in the two last arguments.
2000 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2002 static void do_tpxheader(t_fileio *fio,gmx_bool bRead,t_tpxheader *tpx,
2003 gmx_bool TopOnlyOK, int *file_version,
2004 int *file_generation)
2013 gmx_fio_checktype(fio);
2014 gmx_fio_setdebug(fio,bDebugMode());
2016 /* NEW! XDR tpb file */
2017 precision = sizeof(real);
2019 gmx_fio_do_string(fio,buf);
2020 if (strncmp(buf,"VERSION",7))
2021 gmx_fatal(FARGS,"Can not read file %s,\n"
2022 " this file is from a Gromacs version which is older than 2.0\n"
2023 " Make a new one with grompp or use a gro or pdb file, if possible",
2024 gmx_fio_getname(fio));
2025 gmx_fio_do_int(fio,precision);
2026 bDouble = (precision == sizeof(double));
2027 if ((precision != sizeof(float)) && !bDouble)
2028 gmx_fatal(FARGS,"Unknown precision in file %s: real is %d bytes "
2029 "instead of %d or %d",
2030 gmx_fio_getname(fio),precision,sizeof(float),sizeof(double));
2031 gmx_fio_setprecision(fio,bDouble);
2032 fprintf(stderr,"Reading file %s, %s (%s precision)\n",
2033 gmx_fio_getname(fio),buf,bDouble ? "double" : "single");
2036 gmx_fio_write_string(fio,GromacsVersion());
2037 bDouble = (precision == sizeof(double));
2038 gmx_fio_setprecision(fio,bDouble);
2039 gmx_fio_do_int(fio,precision);
2041 fgen = tpx_generation;
2044 /* Check versions! */
2045 gmx_fio_do_int(fio,fver);
2048 gmx_fio_do_int(fio,fgen);
2052 if(file_version!=NULL)
2053 *file_version = fver;
2054 if(file_generation!=NULL)
2055 *file_generation = fgen;
2058 if ((fver <= tpx_incompatible_version) ||
2059 ((fver > tpx_version) && !TopOnlyOK) ||
2060 (fgen > tpx_generation))
2061 gmx_fatal(FARGS,"reading tpx file (%s) version %d with version %d program",
2062 gmx_fio_getname(fio),fver,tpx_version);
2064 do_section(fio,eitemHEADER,bRead);
2065 gmx_fio_do_int(fio,tpx->natoms);
2067 gmx_fio_do_int(fio,tpx->ngtc);
2071 gmx_fio_do_int(fio,idum);
2072 gmx_fio_do_real(fio,rdum);
2074 gmx_fio_do_real(fio,tpx->lambda);
2075 gmx_fio_do_int(fio,tpx->bIr);
2076 gmx_fio_do_int(fio,tpx->bTop);
2077 gmx_fio_do_int(fio,tpx->bX);
2078 gmx_fio_do_int(fio,tpx->bV);
2079 gmx_fio_do_int(fio,tpx->bF);
2080 gmx_fio_do_int(fio,tpx->bBox);
2082 if((fgen > tpx_generation)) {
2083 /* This can only happen if TopOnlyOK=TRUE */
2088 static int do_tpx(t_fileio *fio, gmx_bool bRead,
2089 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop,
2090 gmx_bool bXVallocated)
2095 gmx_bool TopOnlyOK,bDum=TRUE;
2096 int file_version,file_generation;
2100 gmx_bool bPeriodicMols;
2103 tpx.natoms = state->natoms;
2104 tpx.ngtc = state->ngtc;
2105 tpx.lambda = state->lambda;
2106 tpx.bIr = (ir != NULL);
2107 tpx.bTop = (mtop != NULL);
2108 tpx.bX = (state->x != NULL);
2109 tpx.bV = (state->v != NULL);
2110 tpx.bF = (f != NULL);
2114 TopOnlyOK = (ir==NULL);
2116 do_tpxheader(fio,bRead,&tpx,TopOnlyOK,&file_version,&file_generation);
2120 state->lambda = tpx.lambda;
2121 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
2125 init_state(state,0,tpx.ngtc,0,0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
2126 state->natoms = tpx.natoms;
2127 state->nalloc = tpx.natoms;
2131 init_state(state,tpx.natoms,tpx.ngtc,0,0); /* nose-hoover chains */
2135 #define do_test(fio,b,p) if (bRead && (p!=NULL) && !b) gmx_fatal(FARGS,"No %s in %s",#p,gmx_fio_getname(fio))
2137 do_test(fio,tpx.bBox,state->box);
2138 do_section(fio,eitemBOX,bRead);
2140 gmx_fio_ndo_rvec(fio,state->box,DIM);
2141 if (file_version >= 51) {
2142 gmx_fio_ndo_rvec(fio,state->box_rel,DIM);
2144 /* We initialize box_rel after reading the inputrec */
2145 clear_mat(state->box_rel);
2147 if (file_version >= 28) {
2148 gmx_fio_ndo_rvec(fio,state->boxv,DIM);
2149 if (file_version < 56) {
2151 gmx_fio_ndo_rvec(fio,mdum,DIM);
2156 if (state->ngtc > 0 && file_version >= 28) {
2158 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
2159 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
2160 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
2161 snew(dumv,state->ngtc);
2162 if (file_version < 69) {
2163 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2165 /* These used to be the Berendsen tcoupl_lambda's */
2166 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2170 /* Prior to tpx version 26, the inputrec was here.
2171 * I moved it to enable partial forward-compatibility
2172 * for analysis/viewer programs.
2174 if(file_version<26) {
2175 do_test(fio,tpx.bIr,ir);
2176 do_section(fio,eitemIR,bRead);
2179 do_inputrec(fio, ir,bRead,file_version,
2180 mtop ? &mtop->ffparams.fudgeQQ : NULL);
2182 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2185 do_inputrec(fio, &dum_ir,bRead,file_version,
2186 mtop ? &mtop->ffparams.fudgeQQ :NULL);
2188 pr_inputrec(debug,0,"inputrec",&dum_ir,FALSE);
2189 done_inputrec(&dum_ir);
2195 do_test(fio,tpx.bTop,mtop);
2196 do_section(fio,eitemTOP,bRead);
2199 do_mtop(fio,mtop,bRead, file_version);
2201 do_mtop(fio,&dum_top,bRead,file_version);
2202 done_mtop(&dum_top,TRUE);
2205 do_test(fio,tpx.bX,state->x);
2206 do_section(fio,eitemX,bRead);
2209 state->flags |= (1<<estX);
2211 gmx_fio_ndo_rvec(fio,state->x,state->natoms);
2214 do_test(fio,tpx.bV,state->v);
2215 do_section(fio,eitemV,bRead);
2218 state->flags |= (1<<estV);
2220 gmx_fio_ndo_rvec(fio,state->v,state->natoms);
2223 do_test(fio,tpx.bF,f);
2224 do_section(fio,eitemF,bRead);
2225 if (tpx.bF) gmx_fio_ndo_rvec(fio,f,state->natoms);
2227 /* Starting with tpx version 26, we have the inputrec
2228 * at the end of the file, so we can ignore it
2229 * if the file is never than the software (but still the
2230 * same generation - see comments at the top of this file.
2235 bPeriodicMols = FALSE;
2236 if (file_version >= 26) {
2237 do_test(fio,tpx.bIr,ir);
2238 do_section(fio,eitemIR,bRead);
2240 if (file_version >= 53) {
2241 /* Removed the pbc info from do_inputrec, since we always want it */
2244 bPeriodicMols = ir->bPeriodicMols;
2246 gmx_fio_do_int(fio,ePBC);
2247 gmx_fio_do_gmx_bool(fio,bPeriodicMols);
2249 if (file_generation <= tpx_generation && ir) {
2250 do_inputrec(fio, ir,bRead,file_version,mtop ? &mtop->ffparams.fudgeQQ : NULL);
2252 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2253 if (file_version < 51)
2254 set_box_rel(ir,state);
2255 if (file_version < 53) {
2257 bPeriodicMols = ir->bPeriodicMols;
2260 if (bRead && ir && file_version >= 53) {
2261 /* We need to do this after do_inputrec, since that initializes ir */
2263 ir->bPeriodicMols = bPeriodicMols;
2272 if (state->ngtc == 0)
2274 /* Reading old version without tcoupl state data: set it */
2275 init_gtc_state(state,ir->opts.ngtc,0,ir->opts.nhchainlength);
2277 if (tpx.bTop && mtop)
2279 if (file_version < 57)
2281 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
2283 ir->eDisre = edrSimple;
2287 ir->eDisre = edrNone;
2290 set_disres_npair(mtop);
2294 if (tpx.bTop && mtop)
2296 gmx_mtop_finalize(mtop);
2299 if (file_version >= 57)
2303 env = getenv("GMX_NOCHARGEGROUPS");
2306 sscanf(env,"%d",&ienv);
2307 fprintf(stderr,"\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
2312 "Will make single atomic charge groups in non-solvent%s\n",
2313 ienv > 1 ? " and solvent" : "");
2314 gmx_mtop_make_atomic_charge_groups(mtop,ienv==1);
2316 fprintf(stderr,"\n");
2324 /************************************************************
2326 * The following routines are the exported ones
2328 ************************************************************/
2330 t_fileio *open_tpx(const char *fn,const char *mode)
2332 return gmx_fio_open(fn,mode);
2335 void close_tpx(t_fileio *fio)
2340 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
2341 int *file_version, int *file_generation)
2345 fio = open_tpx(fn,"r");
2346 do_tpxheader(fio,TRUE,tpx,TopOnlyOK,file_version,file_generation);
2350 void write_tpx_state(const char *fn,
2351 t_inputrec *ir,t_state *state,gmx_mtop_t *mtop)
2355 fio = open_tpx(fn,"w");
2356 do_tpx(fio,FALSE,ir,state,NULL,mtop,FALSE);
2360 void read_tpx_state(const char *fn,
2361 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop)
2365 fio = open_tpx(fn,"r");
2366 do_tpx(fio,TRUE,ir,state,f,mtop,FALSE);
2370 int read_tpx(const char *fn,
2371 t_inputrec *ir, matrix box,int *natoms,
2372 rvec *x,rvec *v,rvec *f,gmx_mtop_t *mtop)
2380 fio = open_tpx(fn,"r");
2381 ePBC = do_tpx(fio,TRUE,ir,&state,f,mtop,TRUE);
2383 *natoms = state.natoms;
2385 copy_mat(state.box,box);
2393 int read_tpx_top(const char *fn,
2394 t_inputrec *ir, matrix box,int *natoms,
2395 rvec *x,rvec *v,rvec *f,t_topology *top)
2401 ePBC = read_tpx(fn,ir,box,natoms,x,v,f,&mtop);
2403 *top = gmx_mtop_t_to_t_topology(&mtop);
2408 gmx_bool fn2bTPX(const char *file)
2410 switch (fn2ftp(file)) {
2420 gmx_bool read_tps_conf(const char *infile,char *title,t_topology *top,int *ePBC,
2421 rvec **x,rvec **v,matrix box,gmx_bool bMass)
2424 int natoms,i,version,generation;
2425 gmx_bool bTop,bXNULL;
2427 t_topology *topconv;
2430 bTop = fn2bTPX(infile);
2433 read_tpxheader(infile,&header,TRUE,&version,&generation);
2435 snew(*x,header.natoms);
2437 snew(*v,header.natoms);
2439 *ePBC = read_tpx(infile,NULL,box,&natoms,
2440 (x==NULL) ? NULL : *x,(v==NULL) ? NULL : *v,NULL,mtop);
2441 *top = gmx_mtop_t_to_t_topology(mtop);
2443 strcpy(title,*top->name);
2444 tpx_make_chain_identifiers(&top->atoms,&top->mols);
2447 get_stx_coordnum(infile,&natoms);
2448 init_t_atoms(&top->atoms,natoms,FALSE);
2449 bXNULL = (x == NULL);
2453 read_stx_conf(infile,title,&top->atoms,*x,(v==NULL) ? NULL : *v,ePBC,box);
2459 aps = gmx_atomprop_init();
2460 for(i=0; (i<natoms); i++)
2461 if (!gmx_atomprop_query(aps,epropMass,
2462 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2463 *top->atoms.atomname[i],
2464 &(top->atoms.atom[i].m))) {
2466 fprintf(debug,"Can not find mass for atom %s %d %s, setting to 1\n",
2467 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2468 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
2469 *top->atoms.atomname[i]);
2471 gmx_atomprop_destroy(aps);
2473 top->idef.ntypes=-1;