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 = 78;
78 /* This number should only be increased when you edit the TOPOLOGY section
79 * of the tpx format. This way we can maintain forward compatibility too
80 * for all analysis tools and/or external programs that only need to
81 * know the atom/residue names, charges, and bond connectivity.
83 * It first appeared in tpx version 26, when I also moved the inputrecord
84 * to the end of the tpx file, so we can just skip it if we only
87 static const int tpx_generation = 24;
89 /* This number should be the most recent backwards incompatible version
90 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
92 static const int tpx_incompatible_version = 9;
96 /* Struct used to maintain tpx compatibility when function types are added */
98 int fvnr; /* file version number in which the function type first appeared */
99 int ftype; /* function type */
103 *The entries should be ordered in:
104 * 1. ascending file version number
105 * 2. ascending function type number
107 /*static const t_ftupd ftupd[] = {
108 { 20, F_CUBICBONDS },
112 { 22, F_DISRESVIOL },
118 { 26, F_DIHRESVIOL },
119 { 30, F_CROSS_BOND_BONDS },
120 { 30, F_CROSS_BOND_ANGLES },
121 { 30, F_UREY_BRADLEY },
122 { 30, F_POLARIZATION }
126 *The entries should be ordered in:
127 * 1. ascending function type number
128 * 2. ascending file version number
130 static const t_ftupd ftupd[] = {
131 { 20, F_CUBICBONDS },
136 { 43, F_TABBONDSNC },
137 { 70, F_RESTRBONDS },
138 { 76, F_LINEAR_ANGLES },
139 { 30, F_CROSS_BOND_BONDS },
140 { 30, F_CROSS_BOND_ANGLES },
141 { 30, F_UREY_BRADLEY },
142 { 34, F_QUARTIC_ANGLES },
152 { 72, F_NPSOLVATION },
154 { 41, F_LJC_PAIRS_NB },
157 { 32, F_COUL_RECIP },
159 { 30, F_POLARIZATION },
161 { 22, F_DISRESVIOL },
165 { 26, F_DIHRESVIOL },
170 { 46, F_ECONSERVED },
176 #define NFTUPD asize(ftupd)
178 /* Needed for backward compatibility */
181 static void _do_section(t_fileio *fio,int key,gmx_bool bRead,const char *src,
187 if (gmx_fio_getftp(fio) == efTPA) {
189 gmx_fio_write_string(fio,itemstr[key]);
190 bDbg = gmx_fio_getdebug(fio);
191 gmx_fio_setdebug(fio,FALSE);
192 gmx_fio_write_string(fio,comment_str[key]);
193 gmx_fio_setdebug(fio,bDbg);
196 if (gmx_fio_getdebug(fio))
197 fprintf(stderr,"Looking for section %s (%s, %d)",
198 itemstr[key],src,line);
201 gmx_fio_do_string(fio,buf);
202 } while ((gmx_strcasecmp(buf,itemstr[key]) != 0));
204 if (gmx_strcasecmp(buf,itemstr[key]) != 0)
205 gmx_fatal(FARGS,"\nCould not find section heading %s",itemstr[key]);
206 else if (gmx_fio_getdebug(fio))
207 fprintf(stderr," and found it\n");
212 #define do_section(fio,key,bRead) _do_section(fio,key,bRead,__FILE__,__LINE__)
214 /**************************************************************
216 * Now the higer level routines that do io of the structures and arrays
218 **************************************************************/
219 static void do_pullgrp(t_fileio *fio, t_pullgrp *pgrp, gmx_bool bRead,
225 gmx_fio_do_int(fio,pgrp->nat);
227 snew(pgrp->ind,pgrp->nat);
228 bDum=gmx_fio_ndo_int(fio,pgrp->ind,pgrp->nat);
229 gmx_fio_do_int(fio,pgrp->nweight);
231 snew(pgrp->weight,pgrp->nweight);
232 bDum=gmx_fio_ndo_real(fio,pgrp->weight,pgrp->nweight);
233 gmx_fio_do_int(fio,pgrp->pbcatom);
234 gmx_fio_do_rvec(fio,pgrp->vec);
235 gmx_fio_do_rvec(fio,pgrp->init);
236 gmx_fio_do_real(fio,pgrp->rate);
237 gmx_fio_do_real(fio,pgrp->k);
238 if (file_version >= 56) {
239 gmx_fio_do_real(fio,pgrp->kB);
245 static void do_pull(t_fileio *fio, t_pull *pull,gmx_bool bRead, int file_version)
249 gmx_fio_do_int(fio,pull->ngrp);
250 gmx_fio_do_int(fio,pull->eGeom);
251 gmx_fio_do_ivec(fio,pull->dim);
252 gmx_fio_do_real(fio,pull->cyl_r1);
253 gmx_fio_do_real(fio,pull->cyl_r0);
254 gmx_fio_do_real(fio,pull->constr_tol);
255 gmx_fio_do_int(fio,pull->nstxout);
256 gmx_fio_do_int(fio,pull->nstfout);
258 snew(pull->grp,pull->ngrp+1);
259 for(g=0; g<pull->ngrp+1; g++)
260 do_pullgrp(fio,&pull->grp[g],bRead,file_version);
264 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg,gmx_bool bRead, int file_version)
269 gmx_fio_do_int(fio,rotg->eType);
270 gmx_fio_do_int(fio,rotg->bMassW);
271 gmx_fio_do_int(fio,rotg->nat);
273 snew(rotg->ind,rotg->nat);
274 gmx_fio_ndo_int(fio,rotg->ind,rotg->nat);
276 snew(rotg->x_ref,rotg->nat);
277 gmx_fio_ndo_rvec(fio,rotg->x_ref,rotg->nat);
278 gmx_fio_do_rvec(fio,rotg->vec);
279 gmx_fio_do_rvec(fio,rotg->pivot);
280 gmx_fio_do_real(fio,rotg->rate);
281 gmx_fio_do_real(fio,rotg->k);
282 gmx_fio_do_real(fio,rotg->slab_dist);
283 gmx_fio_do_real(fio,rotg->min_gaussian);
284 gmx_fio_do_real(fio,rotg->eps);
285 gmx_fio_do_int(fio,rotg->eFittype);
286 gmx_fio_do_int(fio,rotg->PotAngle_nstep);
287 gmx_fio_do_real(fio,rotg->PotAngle_step);
290 static void do_rot(t_fileio *fio, t_rot *rot,gmx_bool bRead, int file_version)
294 gmx_fio_do_int(fio,rot->ngrp);
295 gmx_fio_do_int(fio,rot->nstrout);
296 gmx_fio_do_int(fio,rot->nstsout);
298 snew(rot->grp,rot->ngrp);
299 for(g=0; g<rot->ngrp; g++)
300 do_rotgrp(fio, &rot->grp[g],bRead,file_version);
304 static void do_inputrec(t_fileio *fio, t_inputrec *ir,gmx_bool bRead,
305 int file_version, real *fudgeQQ)
307 int i,j,k,*tmp,idum=0;
312 real zerotemptime,finish_t,init_temp,finish_temp;
314 if (file_version != tpx_version)
316 /* Give a warning about features that are not accessible */
317 fprintf(stderr,"Note: file tpx version %d, software tpx version %d\n",
318 file_version,tpx_version);
326 if (file_version == 0)
331 /* Basic inputrec stuff */
332 gmx_fio_do_int(fio,ir->eI);
333 if (file_version >= 62) {
334 gmx_fio_do_gmx_large_int(fio, ir->nsteps);
336 gmx_fio_do_int(fio,idum);
339 if(file_version > 25) {
340 if (file_version >= 62) {
341 gmx_fio_do_gmx_large_int(fio, ir->init_step);
343 gmx_fio_do_int(fio,idum);
344 ir->init_step = idum;
350 if(file_version >= 58)
351 gmx_fio_do_int(fio,ir->simulation_part);
353 ir->simulation_part=1;
355 if (file_version >= 67) {
356 gmx_fio_do_int(fio,ir->nstcalcenergy);
358 ir->nstcalcenergy = 1;
360 if (file_version < 53) {
361 /* The pbc info has been moved out of do_inputrec,
362 * since we always want it, also without reading the inputrec.
364 gmx_fio_do_int(fio,ir->ePBC);
365 if ((file_version <= 15) && (ir->ePBC == 2))
367 if (file_version >= 45) {
368 gmx_fio_do_int(fio,ir->bPeriodicMols);
372 ir->bPeriodicMols = TRUE;
374 ir->bPeriodicMols = FALSE;
378 gmx_fio_do_int(fio,ir->ns_type);
379 gmx_fio_do_int(fio,ir->nstlist);
380 gmx_fio_do_int(fio,ir->ndelta);
381 if (file_version < 41) {
382 gmx_fio_do_int(fio,idum);
383 gmx_fio_do_int(fio,idum);
385 if (file_version >= 45)
386 gmx_fio_do_real(fio,ir->rtpi);
389 gmx_fio_do_int(fio,ir->nstcomm);
390 if (file_version > 34)
391 gmx_fio_do_int(fio,ir->comm_mode);
392 else if (ir->nstcomm < 0)
393 ir->comm_mode = ecmANGULAR;
395 ir->comm_mode = ecmLINEAR;
396 ir->nstcomm = abs(ir->nstcomm);
398 if(file_version > 25)
399 gmx_fio_do_int(fio,ir->nstcheckpoint);
403 gmx_fio_do_int(fio,ir->nstcgsteep);
406 gmx_fio_do_int(fio,ir->nbfgscorr);
410 gmx_fio_do_int(fio,ir->nstlog);
411 gmx_fio_do_int(fio,ir->nstxout);
412 gmx_fio_do_int(fio,ir->nstvout);
413 gmx_fio_do_int(fio,ir->nstfout);
414 gmx_fio_do_int(fio,ir->nstenergy);
415 gmx_fio_do_int(fio,ir->nstxtcout);
416 if (file_version >= 59) {
417 gmx_fio_do_double(fio,ir->init_t);
418 gmx_fio_do_double(fio,ir->delta_t);
420 gmx_fio_do_real(fio,rdum);
422 gmx_fio_do_real(fio,rdum);
425 gmx_fio_do_real(fio,ir->xtcprec);
426 if (file_version < 19) {
427 gmx_fio_do_int(fio,idum);
428 gmx_fio_do_int(fio,idum);
430 if(file_version < 18)
431 gmx_fio_do_int(fio,idum);
432 gmx_fio_do_real(fio,ir->rlist);
433 if (file_version >= 67) {
434 gmx_fio_do_real(fio,ir->rlistlong);
436 gmx_fio_do_int(fio,ir->coulombtype);
437 if (file_version < 32 && ir->coulombtype == eelRF)
438 ir->coulombtype = eelRF_NEC;
439 gmx_fio_do_real(fio,ir->rcoulomb_switch);
440 gmx_fio_do_real(fio,ir->rcoulomb);
441 gmx_fio_do_int(fio,ir->vdwtype);
442 gmx_fio_do_real(fio,ir->rvdw_switch);
443 gmx_fio_do_real(fio,ir->rvdw);
444 if (file_version < 67) {
445 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
447 gmx_fio_do_int(fio,ir->eDispCorr);
448 gmx_fio_do_real(fio,ir->epsilon_r);
449 if (file_version >= 37) {
450 gmx_fio_do_real(fio,ir->epsilon_rf);
452 if (EEL_RF(ir->coulombtype)) {
453 ir->epsilon_rf = ir->epsilon_r;
456 ir->epsilon_rf = 1.0;
459 if (file_version >= 29)
460 gmx_fio_do_real(fio,ir->tabext);
464 if(file_version > 25) {
465 gmx_fio_do_int(fio,ir->gb_algorithm);
466 gmx_fio_do_int(fio,ir->nstgbradii);
467 gmx_fio_do_real(fio,ir->rgbradii);
468 gmx_fio_do_real(fio,ir->gb_saltconc);
469 gmx_fio_do_int(fio,ir->implicit_solvent);
471 ir->gb_algorithm=egbSTILL;
475 ir->implicit_solvent=eisNO;
479 gmx_fio_do_real(fio,ir->gb_epsilon_solvent);
480 gmx_fio_do_real(fio,ir->gb_obc_alpha);
481 gmx_fio_do_real(fio,ir->gb_obc_beta);
482 gmx_fio_do_real(fio,ir->gb_obc_gamma);
485 gmx_fio_do_real(fio,ir->gb_dielectric_offset);
486 gmx_fio_do_int(fio,ir->sa_algorithm);
490 ir->gb_dielectric_offset = 0.009;
491 ir->sa_algorithm = esaAPPROX;
493 gmx_fio_do_real(fio,ir->sa_surface_tension);
495 /* Override sa_surface_tension if it is not changed in the mpd-file */
496 if(ir->sa_surface_tension<0)
498 if(ir->gb_algorithm==egbSTILL)
500 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
502 else if(ir->gb_algorithm==egbHCT || ir->gb_algorithm==egbOBC)
504 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
511 /* Better use sensible values than insane (0.0) ones... */
512 ir->gb_epsilon_solvent = 80;
513 ir->gb_obc_alpha = 1.0;
514 ir->gb_obc_beta = 0.8;
515 ir->gb_obc_gamma = 4.85;
516 ir->sa_surface_tension = 2.092;
520 gmx_fio_do_int(fio,ir->nkx);
521 gmx_fio_do_int(fio,ir->nky);
522 gmx_fio_do_int(fio,ir->nkz);
523 gmx_fio_do_int(fio,ir->pme_order);
524 gmx_fio_do_real(fio,ir->ewald_rtol);
526 if (file_version >=24)
527 gmx_fio_do_int(fio,ir->ewald_geometry);
529 ir->ewald_geometry=eewg3D;
531 if (file_version <=17) {
532 ir->epsilon_surface=0;
533 if (file_version==17)
534 gmx_fio_do_int(fio,idum);
537 gmx_fio_do_real(fio,ir->epsilon_surface);
539 gmx_fio_do_gmx_bool(fio,ir->bOptFFT);
541 gmx_fio_do_gmx_bool(fio,ir->bContinuation);
542 gmx_fio_do_int(fio,ir->etc);
543 /* before version 18, ir->etc was a gmx_bool (ir->btc),
544 * but the values 0 and 1 still mean no and
545 * berendsen temperature coupling, respectively.
547 if (file_version >= 71)
549 gmx_fio_do_int(fio,ir->nsttcouple);
553 ir->nsttcouple = ir->nstcalcenergy;
555 if (file_version <= 15)
557 gmx_fio_do_int(fio,idum);
559 if (file_version <=17)
561 gmx_fio_do_int(fio,ir->epct);
562 if (file_version <= 15)
566 ir->epct = epctSURFACETENSION;
568 gmx_fio_do_int(fio,idum);
571 /* we have removed the NO alternative at the beginning */
575 ir->epct=epctISOTROPIC;
579 ir->epc=epcBERENDSEN;
584 gmx_fio_do_int(fio,ir->epc);
585 gmx_fio_do_int(fio,ir->epct);
587 if (file_version >= 71)
589 gmx_fio_do_int(fio,ir->nstpcouple);
593 ir->nstpcouple = ir->nstcalcenergy;
595 gmx_fio_do_real(fio,ir->tau_p);
596 if (file_version <= 15) {
597 gmx_fio_do_rvec(fio,vdum);
598 clear_mat(ir->ref_p);
600 ir->ref_p[i][i] = vdum[i];
602 gmx_fio_do_rvec(fio,ir->ref_p[XX]);
603 gmx_fio_do_rvec(fio,ir->ref_p[YY]);
604 gmx_fio_do_rvec(fio,ir->ref_p[ZZ]);
606 if (file_version <= 15) {
607 gmx_fio_do_rvec(fio,vdum);
608 clear_mat(ir->compress);
610 ir->compress[i][i] = vdum[i];
613 gmx_fio_do_rvec(fio,ir->compress[XX]);
614 gmx_fio_do_rvec(fio,ir->compress[YY]);
615 gmx_fio_do_rvec(fio,ir->compress[ZZ]);
617 if (file_version >= 47) {
618 gmx_fio_do_int(fio,ir->refcoord_scaling);
619 gmx_fio_do_rvec(fio,ir->posres_com);
620 gmx_fio_do_rvec(fio,ir->posres_comB);
622 ir->refcoord_scaling = erscNO;
623 clear_rvec(ir->posres_com);
624 clear_rvec(ir->posres_comB);
626 if(file_version > 25)
627 gmx_fio_do_int(fio,ir->andersen_seed);
631 if(file_version < 26) {
632 gmx_fio_do_gmx_bool(fio,bSimAnn);
633 gmx_fio_do_real(fio,zerotemptime);
636 if (file_version < 37)
637 gmx_fio_do_real(fio,rdum);
639 gmx_fio_do_real(fio,ir->shake_tol);
640 if (file_version < 54)
641 gmx_fio_do_real(fio,*fudgeQQ);
642 gmx_fio_do_int(fio,ir->efep);
643 if (file_version <= 14 && ir->efep > efepNO)
645 if (file_version >= 59) {
646 gmx_fio_do_double(fio, ir->init_lambda);
647 gmx_fio_do_double(fio, ir->delta_lambda);
649 gmx_fio_do_real(fio,rdum);
650 ir->init_lambda = rdum;
651 gmx_fio_do_real(fio,rdum);
652 ir->delta_lambda = rdum;
654 if (file_version >= 64) {
655 gmx_fio_do_int(fio,ir->n_flambda);
657 snew(ir->flambda,ir->n_flambda);
659 bDum=gmx_fio_ndo_double(fio,ir->flambda,ir->n_flambda);
664 if (file_version >= 13)
665 gmx_fio_do_real(fio,ir->sc_alpha);
668 if (file_version >= 38)
669 gmx_fio_do_int(fio,ir->sc_power);
672 if (file_version >= 15)
673 gmx_fio_do_real(fio,ir->sc_sigma);
678 if (file_version >= 71)
680 ir->sc_sigma_min = ir->sc_sigma;
684 ir->sc_sigma_min = 0;
687 if (file_version >= 64) {
688 gmx_fio_do_int(fio,ir->nstdhdl);
693 if (file_version >= 73)
695 gmx_fio_do_int(fio, ir->separate_dhdl_file);
696 gmx_fio_do_int(fio, ir->dhdl_derivatives);
700 ir->separate_dhdl_file = sepdhdlfileYES;
701 ir->dhdl_derivatives = dhdlderivativesYES;
704 if (file_version >= 71)
706 gmx_fio_do_int(fio,ir->dh_hist_size);
707 gmx_fio_do_double(fio,ir->dh_hist_spacing);
711 ir->dh_hist_size = 0;
712 ir->dh_hist_spacing = 0.1;
714 if (file_version >= 57) {
715 gmx_fio_do_int(fio,ir->eDisre);
717 gmx_fio_do_int(fio,ir->eDisreWeighting);
718 if (file_version < 22) {
719 if (ir->eDisreWeighting == 0)
720 ir->eDisreWeighting = edrwEqual;
722 ir->eDisreWeighting = edrwConservative;
724 gmx_fio_do_gmx_bool(fio,ir->bDisreMixed);
725 gmx_fio_do_real(fio,ir->dr_fc);
726 gmx_fio_do_real(fio,ir->dr_tau);
727 gmx_fio_do_int(fio,ir->nstdisreout);
728 if (file_version >= 22) {
729 gmx_fio_do_real(fio,ir->orires_fc);
730 gmx_fio_do_real(fio,ir->orires_tau);
731 gmx_fio_do_int(fio,ir->nstorireout);
737 if(file_version >= 26) {
738 gmx_fio_do_real(fio,ir->dihre_fc);
739 if (file_version < 56) {
740 gmx_fio_do_real(fio,rdum);
741 gmx_fio_do_int(fio,idum);
747 gmx_fio_do_real(fio,ir->em_stepsize);
748 gmx_fio_do_real(fio,ir->em_tol);
749 if (file_version >= 22)
750 gmx_fio_do_gmx_bool(fio,ir->bShakeSOR);
752 ir->bShakeSOR = TRUE;
753 if (file_version >= 11)
754 gmx_fio_do_int(fio,ir->niter);
757 fprintf(stderr,"Note: niter not in run input file, setting it to %d\n",
760 if (file_version >= 21)
761 gmx_fio_do_real(fio,ir->fc_stepsize);
764 gmx_fio_do_int(fio,ir->eConstrAlg);
765 gmx_fio_do_int(fio,ir->nProjOrder);
766 gmx_fio_do_real(fio,ir->LincsWarnAngle);
767 if (file_version <= 14)
768 gmx_fio_do_int(fio,idum);
769 if (file_version >=26)
770 gmx_fio_do_int(fio,ir->nLincsIter);
773 fprintf(stderr,"Note: nLincsIter not in run input file, setting it to %d\n",
776 if (file_version < 33)
777 gmx_fio_do_real(fio,bd_temp);
778 gmx_fio_do_real(fio,ir->bd_fric);
779 gmx_fio_do_int(fio,ir->ld_seed);
780 if (file_version >= 33) {
782 gmx_fio_do_rvec(fio,ir->deform[i]);
785 clear_rvec(ir->deform[i]);
787 if (file_version >= 14)
788 gmx_fio_do_real(fio,ir->cos_accel);
791 gmx_fio_do_int(fio,ir->userint1);
792 gmx_fio_do_int(fio,ir->userint2);
793 gmx_fio_do_int(fio,ir->userint3);
794 gmx_fio_do_int(fio,ir->userint4);
795 gmx_fio_do_real(fio,ir->userreal1);
796 gmx_fio_do_real(fio,ir->userreal2);
797 gmx_fio_do_real(fio,ir->userreal3);
798 gmx_fio_do_real(fio,ir->userreal4);
801 if (file_version >= 75) {
802 gmx_fio_do_gmx_bool(fio,ir->bAdress);
804 if (bRead) snew(ir->adress, 1);
805 gmx_fio_do_int(fio,ir->adress->type);
806 gmx_fio_do_real(fio,ir->adress->const_wf);
807 gmx_fio_do_real(fio,ir->adress->ex_width);
808 gmx_fio_do_real(fio,ir->adress->hy_width);
809 gmx_fio_do_int(fio,ir->adress->icor);
810 gmx_fio_do_int(fio,ir->adress->site);
811 gmx_fio_do_rvec(fio,ir->adress->refs);
812 gmx_fio_do_int(fio,ir->adress->n_tf_grps);
813 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
814 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
815 gmx_fio_do_int(fio,ir->adress->do_hybridpairs);
817 if (bRead)snew(ir->adress->tf_table_index,ir->adress->n_tf_grps);
818 if (ir->adress->n_tf_grps > 0) {
819 bDum=gmx_fio_ndo_int(fio,ir->adress->tf_table_index,ir->adress->n_tf_grps);
821 if (bRead)snew(ir->adress->group_explicit,ir->adress->n_energy_grps);
822 if (ir->adress->n_energy_grps > 0) {
823 bDum=gmx_fio_ndo_int(fio, ir->adress->group_explicit,ir->adress->n_energy_grps);
831 if (file_version >= 48) {
832 gmx_fio_do_int(fio,ir->ePull);
833 if (ir->ePull != epullNO) {
836 do_pull(fio, ir->pull,bRead,file_version);
842 /* Enforced rotation */
843 if (file_version >= 74) {
844 gmx_fio_do_int(fio,ir->bRot);
845 if (ir->bRot == TRUE) {
848 do_rot(fio, ir->rot,bRead,file_version);
855 gmx_fio_do_int(fio,ir->opts.ngtc);
856 if (file_version >= 69) {
857 gmx_fio_do_int(fio,ir->opts.nhchainlength);
859 ir->opts.nhchainlength = 1;
861 gmx_fio_do_int(fio,ir->opts.ngacc);
862 gmx_fio_do_int(fio,ir->opts.ngfrz);
863 gmx_fio_do_int(fio,ir->opts.ngener);
866 snew(ir->opts.nrdf, ir->opts.ngtc);
867 snew(ir->opts.ref_t, ir->opts.ngtc);
868 snew(ir->opts.annealing, ir->opts.ngtc);
869 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
870 snew(ir->opts.anneal_time, ir->opts.ngtc);
871 snew(ir->opts.anneal_temp, ir->opts.ngtc);
872 snew(ir->opts.tau_t, ir->opts.ngtc);
873 snew(ir->opts.nFreeze,ir->opts.ngfrz);
874 snew(ir->opts.acc, ir->opts.ngacc);
875 snew(ir->opts.egp_flags,ir->opts.ngener*ir->opts.ngener);
877 if (ir->opts.ngtc > 0) {
878 if (bRead && file_version<13) {
879 snew(tmp,ir->opts.ngtc);
880 bDum=gmx_fio_ndo_int(fio,tmp, ir->opts.ngtc);
881 for(i=0; i<ir->opts.ngtc; i++)
882 ir->opts.nrdf[i] = tmp[i];
885 bDum=gmx_fio_ndo_real(fio,ir->opts.nrdf, ir->opts.ngtc);
887 bDum=gmx_fio_ndo_real(fio,ir->opts.ref_t,ir->opts.ngtc);
888 bDum=gmx_fio_ndo_real(fio,ir->opts.tau_t,ir->opts.ngtc);
889 if (file_version<33 && ir->eI==eiBD) {
890 for(i=0; i<ir->opts.ngtc; i++)
891 ir->opts.tau_t[i] = bd_temp;
894 if (ir->opts.ngfrz > 0)
895 bDum=gmx_fio_ndo_ivec(fio,ir->opts.nFreeze,ir->opts.ngfrz);
896 if (ir->opts.ngacc > 0)
897 gmx_fio_ndo_rvec(fio,ir->opts.acc,ir->opts.ngacc);
898 if (file_version >= 12)
899 bDum=gmx_fio_ndo_int(fio,ir->opts.egp_flags,
900 ir->opts.ngener*ir->opts.ngener);
902 if(bRead && file_version < 26) {
903 for(i=0;i<ir->opts.ngtc;i++) {
905 ir->opts.annealing[i] = eannSINGLE;
906 ir->opts.anneal_npoints[i] = 2;
907 snew(ir->opts.anneal_time[i],2);
908 snew(ir->opts.anneal_temp[i],2);
909 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
910 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
911 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
912 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
913 ir->opts.anneal_time[i][0] = ir->init_t;
914 ir->opts.anneal_time[i][1] = finish_t;
915 ir->opts.anneal_temp[i][0] = init_temp;
916 ir->opts.anneal_temp[i][1] = finish_temp;
918 ir->opts.annealing[i] = eannNO;
919 ir->opts.anneal_npoints[i] = 0;
923 /* file version 26 or later */
924 /* First read the lists with annealing and npoints for each group */
925 bDum=gmx_fio_ndo_int(fio,ir->opts.annealing,ir->opts.ngtc);
926 bDum=gmx_fio_ndo_int(fio,ir->opts.anneal_npoints,ir->opts.ngtc);
927 for(j=0;j<(ir->opts.ngtc);j++) {
928 k=ir->opts.anneal_npoints[j];
930 snew(ir->opts.anneal_time[j],k);
931 snew(ir->opts.anneal_temp[j],k);
933 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_time[j],k);
934 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_temp[j],k);
938 if (file_version >= 45) {
939 gmx_fio_do_int(fio,ir->nwall);
940 gmx_fio_do_int(fio,ir->wall_type);
941 if (file_version >= 50)
942 gmx_fio_do_real(fio,ir->wall_r_linpot);
944 ir->wall_r_linpot = -1;
945 gmx_fio_do_int(fio,ir->wall_atomtype[0]);
946 gmx_fio_do_int(fio,ir->wall_atomtype[1]);
947 gmx_fio_do_real(fio,ir->wall_density[0]);
948 gmx_fio_do_real(fio,ir->wall_density[1]);
949 gmx_fio_do_real(fio,ir->wall_ewald_zfac);
953 ir->wall_atomtype[0] = -1;
954 ir->wall_atomtype[1] = -1;
955 ir->wall_density[0] = 0;
956 ir->wall_density[1] = 0;
957 ir->wall_ewald_zfac = 3;
959 /* Cosine stuff for electric fields */
960 for(j=0; (j<DIM); j++) {
961 gmx_fio_do_int(fio,ir->ex[j].n);
962 gmx_fio_do_int(fio,ir->et[j].n);
964 snew(ir->ex[j].a, ir->ex[j].n);
965 snew(ir->ex[j].phi,ir->ex[j].n);
966 snew(ir->et[j].a, ir->et[j].n);
967 snew(ir->et[j].phi,ir->et[j].n);
969 bDum=gmx_fio_ndo_real(fio,ir->ex[j].a, ir->ex[j].n);
970 bDum=gmx_fio_ndo_real(fio,ir->ex[j].phi,ir->ex[j].n);
971 bDum=gmx_fio_ndo_real(fio,ir->et[j].a, ir->et[j].n);
972 bDum=gmx_fio_ndo_real(fio,ir->et[j].phi,ir->et[j].n);
976 if(file_version>=39){
977 gmx_fio_do_gmx_bool(fio,ir->bQMMM);
978 gmx_fio_do_int(fio,ir->QMMMscheme);
979 gmx_fio_do_real(fio,ir->scalefactor);
980 gmx_fio_do_int(fio,ir->opts.ngQM);
982 snew(ir->opts.QMmethod, ir->opts.ngQM);
983 snew(ir->opts.QMbasis, ir->opts.ngQM);
984 snew(ir->opts.QMcharge, ir->opts.ngQM);
985 snew(ir->opts.QMmult, ir->opts.ngQM);
986 snew(ir->opts.bSH, ir->opts.ngQM);
987 snew(ir->opts.CASorbitals, ir->opts.ngQM);
988 snew(ir->opts.CASelectrons,ir->opts.ngQM);
989 snew(ir->opts.SAon, ir->opts.ngQM);
990 snew(ir->opts.SAoff, ir->opts.ngQM);
991 snew(ir->opts.SAsteps, ir->opts.ngQM);
992 snew(ir->opts.bOPT, ir->opts.ngQM);
993 snew(ir->opts.bTS, ir->opts.ngQM);
995 if (ir->opts.ngQM > 0) {
996 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmethod,ir->opts.ngQM);
997 bDum=gmx_fio_ndo_int(fio,ir->opts.QMbasis,ir->opts.ngQM);
998 bDum=gmx_fio_ndo_int(fio,ir->opts.QMcharge,ir->opts.ngQM);
999 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmult,ir->opts.ngQM);
1000 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bSH,ir->opts.ngQM);
1001 bDum=gmx_fio_ndo_int(fio,ir->opts.CASorbitals,ir->opts.ngQM);
1002 bDum=gmx_fio_ndo_int(fio,ir->opts.CASelectrons,ir->opts.ngQM);
1003 bDum=gmx_fio_ndo_real(fio,ir->opts.SAon,ir->opts.ngQM);
1004 bDum=gmx_fio_ndo_real(fio,ir->opts.SAoff,ir->opts.ngQM);
1005 bDum=gmx_fio_ndo_int(fio,ir->opts.SAsteps,ir->opts.ngQM);
1006 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bOPT,ir->opts.ngQM);
1007 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bTS,ir->opts.ngQM);
1009 /* end of QMMM stuff */
1014 static void do_harm(t_fileio *fio, t_iparams *iparams,gmx_bool bRead)
1016 gmx_fio_do_real(fio,iparams->harmonic.rA);
1017 gmx_fio_do_real(fio,iparams->harmonic.krA);
1018 gmx_fio_do_real(fio,iparams->harmonic.rB);
1019 gmx_fio_do_real(fio,iparams->harmonic.krB);
1022 void do_iparams(t_fileio *fio, t_functype ftype,t_iparams *iparams,
1023 gmx_bool bRead, int file_version)
1030 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1038 do_harm(fio, iparams,bRead);
1039 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead) {
1040 /* Correct incorrect storage of parameters */
1041 iparams->pdihs.phiB = iparams->pdihs.phiA;
1042 iparams->pdihs.cpB = iparams->pdihs.cpA;
1045 case F_LINEAR_ANGLES:
1046 gmx_fio_do_real(fio,iparams->linangle.klinA);
1047 gmx_fio_do_real(fio,iparams->linangle.aA);
1048 gmx_fio_do_real(fio,iparams->linangle.klinB);
1049 gmx_fio_do_real(fio,iparams->linangle.aB);
1052 gmx_fio_do_real(fio,iparams->fene.bm);
1053 gmx_fio_do_real(fio,iparams->fene.kb);
1056 gmx_fio_do_real(fio,iparams->restraint.lowA);
1057 gmx_fio_do_real(fio,iparams->restraint.up1A);
1058 gmx_fio_do_real(fio,iparams->restraint.up2A);
1059 gmx_fio_do_real(fio,iparams->restraint.kA);
1060 gmx_fio_do_real(fio,iparams->restraint.lowB);
1061 gmx_fio_do_real(fio,iparams->restraint.up1B);
1062 gmx_fio_do_real(fio,iparams->restraint.up2B);
1063 gmx_fio_do_real(fio,iparams->restraint.kB);
1069 gmx_fio_do_real(fio,iparams->tab.kA);
1070 gmx_fio_do_int(fio,iparams->tab.table);
1071 gmx_fio_do_real(fio,iparams->tab.kB);
1073 case F_CROSS_BOND_BONDS:
1074 gmx_fio_do_real(fio,iparams->cross_bb.r1e);
1075 gmx_fio_do_real(fio,iparams->cross_bb.r2e);
1076 gmx_fio_do_real(fio,iparams->cross_bb.krr);
1078 case F_CROSS_BOND_ANGLES:
1079 gmx_fio_do_real(fio,iparams->cross_ba.r1e);
1080 gmx_fio_do_real(fio,iparams->cross_ba.r2e);
1081 gmx_fio_do_real(fio,iparams->cross_ba.r3e);
1082 gmx_fio_do_real(fio,iparams->cross_ba.krt);
1084 case F_UREY_BRADLEY:
1085 gmx_fio_do_real(fio,iparams->u_b.theta);
1086 gmx_fio_do_real(fio,iparams->u_b.ktheta);
1087 gmx_fio_do_real(fio,iparams->u_b.r13);
1088 gmx_fio_do_real(fio,iparams->u_b.kUB);
1090 case F_QUARTIC_ANGLES:
1091 gmx_fio_do_real(fio,iparams->qangle.theta);
1092 bDum=gmx_fio_ndo_real(fio,iparams->qangle.c,5);
1095 gmx_fio_do_real(fio,iparams->bham.a);
1096 gmx_fio_do_real(fio,iparams->bham.b);
1097 gmx_fio_do_real(fio,iparams->bham.c);
1100 gmx_fio_do_real(fio,iparams->morse.b0);
1101 gmx_fio_do_real(fio,iparams->morse.cb);
1102 gmx_fio_do_real(fio,iparams->morse.beta);
1105 gmx_fio_do_real(fio,iparams->cubic.b0);
1106 gmx_fio_do_real(fio,iparams->cubic.kb);
1107 gmx_fio_do_real(fio,iparams->cubic.kcub);
1111 case F_POLARIZATION:
1112 gmx_fio_do_real(fio,iparams->polarize.alpha);
1115 gmx_fio_do_real(fio,iparams->anharm_polarize.alpha);
1116 gmx_fio_do_real(fio,iparams->anharm_polarize.drcut);
1117 gmx_fio_do_real(fio,iparams->anharm_polarize.khyp);
1120 if (file_version < 31)
1121 gmx_fatal(FARGS,"Old tpr files with water_polarization not supported. Make a new.");
1122 gmx_fio_do_real(fio,iparams->wpol.al_x);
1123 gmx_fio_do_real(fio,iparams->wpol.al_y);
1124 gmx_fio_do_real(fio,iparams->wpol.al_z);
1125 gmx_fio_do_real(fio,iparams->wpol.rOH);
1126 gmx_fio_do_real(fio,iparams->wpol.rHH);
1127 gmx_fio_do_real(fio,iparams->wpol.rOD);
1130 gmx_fio_do_real(fio,iparams->thole.a);
1131 gmx_fio_do_real(fio,iparams->thole.alpha1);
1132 gmx_fio_do_real(fio,iparams->thole.alpha2);
1133 gmx_fio_do_real(fio,iparams->thole.rfac);
1136 gmx_fio_do_real(fio,iparams->lj.c6);
1137 gmx_fio_do_real(fio,iparams->lj.c12);
1140 gmx_fio_do_real(fio,iparams->lj14.c6A);
1141 gmx_fio_do_real(fio,iparams->lj14.c12A);
1142 gmx_fio_do_real(fio,iparams->lj14.c6B);
1143 gmx_fio_do_real(fio,iparams->lj14.c12B);
1146 gmx_fio_do_real(fio,iparams->ljc14.fqq);
1147 gmx_fio_do_real(fio,iparams->ljc14.qi);
1148 gmx_fio_do_real(fio,iparams->ljc14.qj);
1149 gmx_fio_do_real(fio,iparams->ljc14.c6);
1150 gmx_fio_do_real(fio,iparams->ljc14.c12);
1152 case F_LJC_PAIRS_NB:
1153 gmx_fio_do_real(fio,iparams->ljcnb.qi);
1154 gmx_fio_do_real(fio,iparams->ljcnb.qj);
1155 gmx_fio_do_real(fio,iparams->ljcnb.c6);
1156 gmx_fio_do_real(fio,iparams->ljcnb.c12);
1162 gmx_fio_do_real(fio,iparams->pdihs.phiA);
1163 gmx_fio_do_real(fio,iparams->pdihs.cpA);
1164 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42) {
1165 /* Read the incorrectly stored multiplicity */
1166 gmx_fio_do_real(fio,iparams->harmonic.rB);
1167 gmx_fio_do_real(fio,iparams->harmonic.krB);
1168 iparams->pdihs.phiB = iparams->pdihs.phiA;
1169 iparams->pdihs.cpB = iparams->pdihs.cpA;
1171 gmx_fio_do_real(fio,iparams->pdihs.phiB);
1172 gmx_fio_do_real(fio,iparams->pdihs.cpB);
1173 gmx_fio_do_int(fio,iparams->pdihs.mult);
1177 gmx_fio_do_int(fio,iparams->disres.label);
1178 gmx_fio_do_int(fio,iparams->disres.type);
1179 gmx_fio_do_real(fio,iparams->disres.low);
1180 gmx_fio_do_real(fio,iparams->disres.up1);
1181 gmx_fio_do_real(fio,iparams->disres.up2);
1182 gmx_fio_do_real(fio,iparams->disres.kfac);
1185 gmx_fio_do_int(fio,iparams->orires.ex);
1186 gmx_fio_do_int(fio,iparams->orires.label);
1187 gmx_fio_do_int(fio,iparams->orires.power);
1188 gmx_fio_do_real(fio,iparams->orires.c);
1189 gmx_fio_do_real(fio,iparams->orires.obs);
1190 gmx_fio_do_real(fio,iparams->orires.kfac);
1193 gmx_fio_do_int(fio,iparams->dihres.power);
1194 gmx_fio_do_int(fio,iparams->dihres.label);
1195 gmx_fio_do_real(fio,iparams->dihres.phi);
1196 gmx_fio_do_real(fio,iparams->dihres.dphi);
1197 gmx_fio_do_real(fio,iparams->dihres.kfac);
1200 gmx_fio_do_rvec(fio,iparams->posres.pos0A);
1201 gmx_fio_do_rvec(fio,iparams->posres.fcA);
1202 if (bRead && file_version < 27) {
1203 copy_rvec(iparams->posres.pos0A,iparams->posres.pos0B);
1204 copy_rvec(iparams->posres.fcA,iparams->posres.fcB);
1206 gmx_fio_do_rvec(fio,iparams->posres.pos0B);
1207 gmx_fio_do_rvec(fio,iparams->posres.fcB);
1211 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA,NR_RBDIHS);
1212 if(file_version>=25)
1213 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB,NR_RBDIHS);
1216 /* Fourier dihedrals are internally represented
1217 * as Ryckaert-Bellemans since those are faster to compute.
1219 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA, NR_RBDIHS);
1220 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB, NR_RBDIHS);
1224 gmx_fio_do_real(fio,iparams->constr.dA);
1225 gmx_fio_do_real(fio,iparams->constr.dB);
1228 gmx_fio_do_real(fio,iparams->settle.doh);
1229 gmx_fio_do_real(fio,iparams->settle.dhh);
1232 gmx_fio_do_real(fio,iparams->vsite.a);
1237 gmx_fio_do_real(fio,iparams->vsite.a);
1238 gmx_fio_do_real(fio,iparams->vsite.b);
1243 gmx_fio_do_real(fio,iparams->vsite.a);
1244 gmx_fio_do_real(fio,iparams->vsite.b);
1245 gmx_fio_do_real(fio,iparams->vsite.c);
1248 gmx_fio_do_int(fio,iparams->vsiten.n);
1249 gmx_fio_do_real(fio,iparams->vsiten.a);
1254 /* We got rid of some parameters in version 68 */
1255 if(bRead && file_version<68)
1257 gmx_fio_do_real(fio,rdum);
1258 gmx_fio_do_real(fio,rdum);
1259 gmx_fio_do_real(fio,rdum);
1260 gmx_fio_do_real(fio,rdum);
1262 gmx_fio_do_real(fio,iparams->gb.sar);
1263 gmx_fio_do_real(fio,iparams->gb.st);
1264 gmx_fio_do_real(fio,iparams->gb.pi);
1265 gmx_fio_do_real(fio,iparams->gb.gbr);
1266 gmx_fio_do_real(fio,iparams->gb.bmlt);
1269 gmx_fio_do_int(fio,iparams->cmap.cmapA);
1270 gmx_fio_do_int(fio,iparams->cmap.cmapB);
1273 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
1275 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
1278 gmx_fio_unset_comment(fio);
1281 static void do_ilist(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,int file_version,
1288 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1290 if (file_version < 44) {
1291 for(i=0; i<MAXNODES; i++)
1292 gmx_fio_do_int(fio,idum);
1294 gmx_fio_do_int(fio,ilist->nr);
1296 snew(ilist->iatoms,ilist->nr);
1297 bDum=gmx_fio_ndo_int(fio,ilist->iatoms,ilist->nr);
1299 gmx_fio_unset_comment(fio);
1302 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1303 gmx_bool bRead, int file_version)
1309 gmx_fio_do_int(fio,ffparams->atnr);
1310 if (file_version < 57) {
1311 gmx_fio_do_int(fio,idum);
1313 gmx_fio_do_int(fio,ffparams->ntypes);
1315 fprintf(debug,"ffparams->atnr = %d, ntypes = %d\n",
1316 ffparams->atnr,ffparams->ntypes);
1318 snew(ffparams->functype,ffparams->ntypes);
1319 snew(ffparams->iparams,ffparams->ntypes);
1321 /* Read/write all the function types */
1322 bDum=gmx_fio_ndo_int(fio,ffparams->functype,ffparams->ntypes);
1324 pr_ivec(debug,0,"functype",ffparams->functype,ffparams->ntypes,TRUE);
1326 if (file_version >= 66) {
1327 gmx_fio_do_double(fio,ffparams->reppow);
1329 ffparams->reppow = 12.0;
1332 if (file_version >= 57) {
1333 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1336 /* Check whether all these function types are supported by the code.
1337 * In practice the code is backwards compatible, which means that the
1338 * numbering may have to be altered from old numbering to new numbering
1340 for (i=0; (i<ffparams->ntypes); i++) {
1342 /* Loop over file versions */
1343 for (k=0; (k<NFTUPD); k++)
1344 /* Compare the read file_version to the update table */
1345 if ((file_version < ftupd[k].fvnr) &&
1346 (ffparams->functype[i] >= ftupd[k].ftype)) {
1347 ffparams->functype[i] += 1;
1349 fprintf(debug,"Incrementing function type %d to %d (due to %s)\n",
1350 i,ffparams->functype[i],
1351 interaction_function[ftupd[k].ftype].longname);
1356 do_iparams(fio, ffparams->functype[i],&ffparams->iparams[i],bRead,
1359 pr_iparams(debug,ffparams->functype[i],&ffparams->iparams[i]);
1363 static void add_settle_atoms(t_ilist *ilist)
1367 /* Settle used to only store the first atom: add the other two */
1368 srenew(ilist->iatoms,2*ilist->nr);
1369 for(i=ilist->nr/2-1; i>=0; i--)
1371 ilist->iatoms[4*i+0] = ilist->iatoms[2*i+0];
1372 ilist->iatoms[4*i+1] = ilist->iatoms[2*i+1];
1373 ilist->iatoms[4*i+2] = ilist->iatoms[2*i+1] + 1;
1374 ilist->iatoms[4*i+3] = ilist->iatoms[2*i+1] + 2;
1376 ilist->nr = 2*ilist->nr;
1379 static void do_ilists(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,
1382 int i,j,renum[F_NRE];
1383 gmx_bool bDum=TRUE,bClear;
1386 for(j=0; (j<F_NRE); j++) {
1389 for (k=0; k<NFTUPD; k++)
1390 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
1394 ilist[j].iatoms = NULL;
1396 do_ilist(fio, &ilist[j],bRead,file_version,j);
1397 if (file_version < 78 && j == F_SETTLE && ilist[j].nr > 0)
1399 add_settle_atoms(&ilist[j]);
1403 if (bRead && gmx_debug_at)
1404 pr_ilist(debug,0,interaction_function[j].longname,
1405 functype,&ilist[j],TRUE);
1410 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams,gmx_moltype_t *molt,
1411 gmx_bool bRead, int file_version)
1413 do_ffparams(fio, ffparams,bRead,file_version);
1415 if (file_version >= 54) {
1416 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1419 do_ilists(fio, molt->ilist,bRead,file_version);
1422 static void do_block(t_fileio *fio, t_block *block,gmx_bool bRead,int file_version)
1424 int i,idum,dum_nra,*dum_a;
1427 if (file_version < 44)
1428 for(i=0; i<MAXNODES; i++)
1429 gmx_fio_do_int(fio,idum);
1430 gmx_fio_do_int(fio,block->nr);
1431 if (file_version < 51)
1432 gmx_fio_do_int(fio,dum_nra);
1434 block->nalloc_index = block->nr+1;
1435 snew(block->index,block->nalloc_index);
1437 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1439 if (file_version < 51 && dum_nra > 0) {
1440 snew(dum_a,dum_nra);
1441 bDum=gmx_fio_ndo_int(fio,dum_a,dum_nra);
1446 static void do_blocka(t_fileio *fio, t_blocka *block,gmx_bool bRead,
1452 if (file_version < 44)
1453 for(i=0; i<MAXNODES; i++)
1454 gmx_fio_do_int(fio,idum);
1455 gmx_fio_do_int(fio,block->nr);
1456 gmx_fio_do_int(fio,block->nra);
1458 block->nalloc_index = block->nr+1;
1459 snew(block->index,block->nalloc_index);
1460 block->nalloc_a = block->nra;
1461 snew(block->a,block->nalloc_a);
1463 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1464 bDum=gmx_fio_ndo_int(fio,block->a,block->nra);
1467 static void do_atom(t_fileio *fio, t_atom *atom,int ngrp,gmx_bool bRead,
1468 int file_version, gmx_groups_t *groups,int atnr)
1472 gmx_fio_do_real(fio,atom->m);
1473 gmx_fio_do_real(fio,atom->q);
1474 gmx_fio_do_real(fio,atom->mB);
1475 gmx_fio_do_real(fio,atom->qB);
1476 gmx_fio_do_ushort(fio, atom->type);
1477 gmx_fio_do_ushort(fio, atom->typeB);
1478 gmx_fio_do_int(fio,atom->ptype);
1479 gmx_fio_do_int(fio,atom->resind);
1480 if (file_version >= 52)
1481 gmx_fio_do_int(fio,atom->atomnumber);
1483 atom->atomnumber = NOTSET;
1484 if (file_version < 23)
1486 else if (file_version < 39)
1491 if (file_version < 57) {
1492 unsigned char uchar[egcNR];
1493 gmx_fio_ndo_uchar(fio,uchar,myngrp);
1494 for(i=myngrp; (i<ngrp); i++) {
1497 /* Copy the old data format to the groups struct */
1498 for(i=0; i<ngrp; i++) {
1499 groups->grpnr[i][atnr] = uchar[i];
1504 static void do_grps(t_fileio *fio, int ngrp,t_grps grps[],gmx_bool bRead,
1510 if (file_version < 23)
1512 else if (file_version < 39)
1517 for(j=0; (j<ngrp); j++) {
1519 gmx_fio_do_int(fio,grps[j].nr);
1521 snew(grps[j].nm_ind,grps[j].nr);
1522 bDum=gmx_fio_ndo_int(fio,grps[j].nm_ind,grps[j].nr);
1526 snew(grps[j].nm_ind,grps[j].nr);
1531 static void do_symstr(t_fileio *fio, char ***nm,gmx_bool bRead,t_symtab *symtab)
1536 gmx_fio_do_int(fio,ls);
1537 *nm = get_symtab_handle(symtab,ls);
1540 ls = lookup_symtab(symtab,*nm);
1541 gmx_fio_do_int(fio,ls);
1545 static void do_strstr(t_fileio *fio, int nstr,char ***nm,gmx_bool bRead,
1550 for (j=0; (j<nstr); j++)
1551 do_symstr(fio, &(nm[j]),bRead,symtab);
1554 static void do_resinfo(t_fileio *fio, int n,t_resinfo *ri,gmx_bool bRead,
1555 t_symtab *symtab, int file_version)
1559 for (j=0; (j<n); j++) {
1560 do_symstr(fio, &(ri[j].name),bRead,symtab);
1561 if (file_version >= 63) {
1562 gmx_fio_do_int(fio,ri[j].nr);
1563 gmx_fio_do_uchar(fio, ri[j].ic);
1571 static void do_atoms(t_fileio *fio, t_atoms *atoms,gmx_bool bRead,t_symtab *symtab,
1573 gmx_groups_t *groups)
1577 gmx_fio_do_int(fio,atoms->nr);
1578 gmx_fio_do_int(fio,atoms->nres);
1579 if (file_version < 57) {
1580 gmx_fio_do_int(fio,groups->ngrpname);
1581 for(i=0; i<egcNR; i++) {
1582 groups->ngrpnr[i] = atoms->nr;
1583 snew(groups->grpnr[i],groups->ngrpnr[i]);
1587 snew(atoms->atom,atoms->nr);
1588 snew(atoms->atomname,atoms->nr);
1589 snew(atoms->atomtype,atoms->nr);
1590 snew(atoms->atomtypeB,atoms->nr);
1591 snew(atoms->resinfo,atoms->nres);
1592 if (file_version < 57) {
1593 snew(groups->grpname,groups->ngrpname);
1595 atoms->pdbinfo = NULL;
1597 for(i=0; (i<atoms->nr); i++) {
1598 do_atom(fio, &atoms->atom[i],egcNR,bRead, file_version,groups,i);
1600 do_strstr(fio, atoms->nr,atoms->atomname,bRead,symtab);
1601 if (bRead && (file_version <= 20)) {
1602 for(i=0; i<atoms->nr; i++) {
1603 atoms->atomtype[i] = put_symtab(symtab,"?");
1604 atoms->atomtypeB[i] = put_symtab(symtab,"?");
1607 do_strstr(fio, atoms->nr,atoms->atomtype,bRead,symtab);
1608 do_strstr(fio, atoms->nr,atoms->atomtypeB,bRead,symtab);
1610 do_resinfo(fio, atoms->nres,atoms->resinfo,bRead,symtab,file_version);
1612 if (file_version < 57) {
1613 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1615 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1619 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
1620 gmx_bool bRead,t_symtab *symtab,
1626 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1627 gmx_fio_do_int(fio,groups->ngrpname);
1629 snew(groups->grpname,groups->ngrpname);
1631 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1632 for(g=0; g<egcNR; g++) {
1633 gmx_fio_do_int(fio,groups->ngrpnr[g]);
1634 if (groups->ngrpnr[g] == 0) {
1636 groups->grpnr[g] = NULL;
1640 snew(groups->grpnr[g],groups->ngrpnr[g]);
1642 bDum=gmx_fio_ndo_uchar(fio, groups->grpnr[g],groups->ngrpnr[g]);
1647 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes,gmx_bool bRead,
1648 t_symtab *symtab,int file_version)
1651 gmx_bool bDum = TRUE;
1653 if (file_version > 25) {
1654 gmx_fio_do_int(fio,atomtypes->nr);
1657 snew(atomtypes->radius,j);
1658 snew(atomtypes->vol,j);
1659 snew(atomtypes->surftens,j);
1660 snew(atomtypes->atomnumber,j);
1661 snew(atomtypes->gb_radius,j);
1662 snew(atomtypes->S_hct,j);
1664 bDum=gmx_fio_ndo_real(fio,atomtypes->radius,j);
1665 bDum=gmx_fio_ndo_real(fio,atomtypes->vol,j);
1666 bDum=gmx_fio_ndo_real(fio,atomtypes->surftens,j);
1667 if(file_version >= 40)
1669 bDum=gmx_fio_ndo_int(fio,atomtypes->atomnumber,j);
1671 if(file_version >= 60)
1673 bDum=gmx_fio_ndo_real(fio,atomtypes->gb_radius,j);
1674 bDum=gmx_fio_ndo_real(fio,atomtypes->S_hct,j);
1677 /* File versions prior to 26 cannot do GBSA,
1678 * so they dont use this structure
1681 atomtypes->radius = NULL;
1682 atomtypes->vol = NULL;
1683 atomtypes->surftens = NULL;
1684 atomtypes->atomnumber = NULL;
1685 atomtypes->gb_radius = NULL;
1686 atomtypes->S_hct = NULL;
1690 static void do_symtab(t_fileio *fio, t_symtab *symtab,gmx_bool bRead)
1696 gmx_fio_do_int(fio,symtab->nr);
1699 snew(symtab->symbuf,1);
1700 symbuf = symtab->symbuf;
1701 symbuf->bufsize = nr;
1702 snew(symbuf->buf,nr);
1703 for (i=0; (i<nr); i++) {
1704 gmx_fio_do_string(fio,buf);
1705 symbuf->buf[i]=strdup(buf);
1709 symbuf = symtab->symbuf;
1710 while (symbuf!=NULL) {
1711 for (i=0; (i<symbuf->bufsize) && (i<nr); i++)
1712 gmx_fio_do_string(fio,symbuf->buf[i]);
1714 symbuf=symbuf->next;
1717 gmx_fatal(FARGS,"nr of symtab strings left: %d",nr);
1721 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
1723 int i,j,ngrid,gs,nelem;
1725 gmx_fio_do_int(fio,cmap_grid->ngrid);
1726 gmx_fio_do_int(fio,cmap_grid->grid_spacing);
1728 ngrid = cmap_grid->ngrid;
1729 gs = cmap_grid->grid_spacing;
1734 snew(cmap_grid->cmapdata,ngrid);
1736 for(i=0;i<cmap_grid->ngrid;i++)
1738 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
1742 for(i=0;i<cmap_grid->ngrid;i++)
1744 for(j=0;j<nelem;j++)
1746 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4]);
1747 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+1]);
1748 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+2]);
1749 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+3]);
1755 void tpx_make_chain_identifiers(t_atoms *atoms,t_block *mols)
1761 /* We always assign a new chain number, but save the chain id characters
1762 * for larger molecules.
1764 #define CHAIN_MIN_ATOMS 15
1768 for(m=0; m<mols->nr; m++)
1771 a1=mols->index[m+1];
1772 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
1781 for(a=a0; a<a1; a++)
1783 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
1784 atoms->resinfo[atoms->atom[a].resind].chainid = c;
1789 /* Blank out the chain id if there was only one chain */
1792 for(r=0; r<atoms->nres; r++)
1794 atoms->resinfo[r].chainid = ' ';
1799 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt,gmx_bool bRead,
1800 t_symtab *symtab, int file_version,
1801 gmx_groups_t *groups)
1805 if (file_version >= 57) {
1806 do_symstr(fio, &(molt->name),bRead,symtab);
1809 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
1811 if (bRead && gmx_debug_at) {
1812 pr_atoms(debug,0,"atoms",&molt->atoms,TRUE);
1815 if (file_version >= 57) {
1816 do_ilists(fio, molt->ilist,bRead,file_version);
1818 do_block(fio, &molt->cgs,bRead,file_version);
1819 if (bRead && gmx_debug_at) {
1820 pr_block(debug,0,"cgs",&molt->cgs,TRUE);
1824 /* This used to be in the atoms struct */
1825 do_blocka(fio, &molt->excls, bRead, file_version);
1828 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb,gmx_bool bRead,
1833 gmx_fio_do_int(fio,molb->type);
1834 gmx_fio_do_int(fio,molb->nmol);
1835 gmx_fio_do_int(fio,molb->natoms_mol);
1836 /* Position restraint coordinates */
1837 gmx_fio_do_int(fio,molb->nposres_xA);
1838 if (molb->nposres_xA > 0) {
1840 snew(molb->posres_xA,molb->nposres_xA);
1842 gmx_fio_ndo_rvec(fio,molb->posres_xA,molb->nposres_xA);
1844 gmx_fio_do_int(fio,molb->nposres_xB);
1845 if (molb->nposres_xB > 0) {
1847 snew(molb->posres_xB,molb->nposres_xB);
1849 gmx_fio_ndo_rvec(fio,molb->posres_xB,molb->nposres_xB);
1854 static t_block mtop_mols(gmx_mtop_t *mtop)
1860 for(mb=0; mb<mtop->nmolblock; mb++) {
1861 mols.nr += mtop->molblock[mb].nmol;
1863 mols.nalloc_index = mols.nr + 1;
1864 snew(mols.index,mols.nalloc_index);
1869 for(mb=0; mb<mtop->nmolblock; mb++) {
1870 for(mol=0; mol<mtop->molblock[mb].nmol; mol++) {
1871 a += mtop->molblock[mb].natoms_mol;
1880 static void add_posres_molblock(gmx_mtop_t *mtop)
1885 gmx_molblock_t *molb;
1888 il = &mtop->moltype[0].ilist[F_POSRES];
1894 for(i=0; i<il->nr; i+=2) {
1895 ip = &mtop->ffparams.iparams[il->iatoms[i]];
1896 am = max(am,il->iatoms[i+1]);
1897 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
1898 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
1899 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ]) {
1903 /* Make the posres coordinate block end at a molecule end */
1905 while(am >= mtop->mols.index[mol+1]) {
1908 molb = &mtop->molblock[0];
1909 molb->nposres_xA = mtop->mols.index[mol+1];
1910 snew(molb->posres_xA,molb->nposres_xA);
1912 molb->nposres_xB = molb->nposres_xA;
1913 snew(molb->posres_xB,molb->nposres_xB);
1915 molb->nposres_xB = 0;
1917 for(i=0; i<il->nr; i+=2) {
1918 ip = &mtop->ffparams.iparams[il->iatoms[i]];
1919 a = il->iatoms[i+1];
1920 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
1921 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
1922 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
1924 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
1925 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
1926 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
1931 static void set_disres_npair(gmx_mtop_t *mtop)
1938 ip = mtop->ffparams.iparams;
1940 for(mt=0; mt<mtop->nmoltype; mt++) {
1941 il = &mtop->moltype[mt].ilist[F_DISRES];
1945 for(i=0; i<il->nr; i+=3) {
1947 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label) {
1948 ip[a[i]].disres.npair = npair;
1956 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop,gmx_bool bRead,
1964 do_symtab(fio, &(mtop->symtab),bRead);
1966 pr_symtab(debug,0,"symtab",&mtop->symtab);
1968 do_symstr(fio, &(mtop->name),bRead,&(mtop->symtab));
1970 if (file_version >= 57) {
1971 do_ffparams(fio, &mtop->ffparams,bRead,file_version);
1973 gmx_fio_do_int(fio,mtop->nmoltype);
1978 snew(mtop->moltype,mtop->nmoltype);
1979 if (file_version < 57) {
1980 mtop->moltype[0].name = mtop->name;
1983 for(mt=0; mt<mtop->nmoltype; mt++) {
1984 do_moltype(fio, &mtop->moltype[mt],bRead,&mtop->symtab,file_version,
1988 if (file_version >= 57) {
1989 gmx_fio_do_int(fio,mtop->nmolblock);
1991 mtop->nmolblock = 1;
1994 snew(mtop->molblock,mtop->nmolblock);
1996 if (file_version >= 57) {
1997 for(mb=0; mb<mtop->nmolblock; mb++) {
1998 do_molblock(fio, &mtop->molblock[mb],bRead,file_version);
2000 gmx_fio_do_int(fio,mtop->natoms);
2002 mtop->molblock[0].type = 0;
2003 mtop->molblock[0].nmol = 1;
2004 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
2005 mtop->molblock[0].nposres_xA = 0;
2006 mtop->molblock[0].nposres_xB = 0;
2009 do_atomtypes (fio, &(mtop->atomtypes),bRead,&(mtop->symtab), file_version);
2011 pr_atomtypes(debug,0,"atomtypes",&mtop->atomtypes,TRUE);
2013 if (file_version < 57) {
2014 /* Debug statements are inside do_idef */
2015 do_idef (fio, &mtop->ffparams,&mtop->moltype[0],bRead,file_version);
2016 mtop->natoms = mtop->moltype[0].atoms.nr;
2019 if(file_version >= 65)
2021 do_cmap(fio, &mtop->ffparams.cmap_grid,bRead);
2025 mtop->ffparams.cmap_grid.ngrid = 0;
2026 mtop->ffparams.cmap_grid.grid_spacing = 0;
2027 mtop->ffparams.cmap_grid.cmapdata = NULL;
2030 if (file_version >= 57) {
2031 do_groups(fio, &mtop->groups,bRead,&(mtop->symtab),file_version);
2034 if (file_version < 57) {
2035 do_block(fio, &mtop->moltype[0].cgs,bRead,file_version);
2036 if (bRead && gmx_debug_at) {
2037 pr_block(debug,0,"cgs",&mtop->moltype[0].cgs,TRUE);
2039 do_block(fio, &mtop->mols,bRead,file_version);
2040 /* Add the posres coordinates to the molblock */
2041 add_posres_molblock(mtop);
2044 if (file_version >= 57) {
2045 mtop->mols = mtop_mols(mtop);
2048 pr_block(debug,0,"mols",&mtop->mols,TRUE);
2052 if (file_version < 51) {
2053 /* Here used to be the shake blocks */
2054 do_blocka(fio, &dumb,bRead,file_version);
2062 close_symtab(&(mtop->symtab));
2066 /* If TopOnlyOK is TRUE then we can read even future versions
2067 * of tpx files, provided the file_generation hasn't changed.
2068 * If it is FALSE, we need the inputrecord too, and bail out
2069 * if the file is newer than the program.
2071 * The version and generation if the topology (see top of this file)
2072 * are returned in the two last arguments.
2074 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2076 static void do_tpxheader(t_fileio *fio,gmx_bool bRead,t_tpxheader *tpx,
2077 gmx_bool TopOnlyOK, int *file_version,
2078 int *file_generation)
2081 char file_tag[STRLEN];
2088 gmx_fio_checktype(fio);
2089 gmx_fio_setdebug(fio,bDebugMode());
2091 /* NEW! XDR tpb file */
2092 precision = sizeof(real);
2094 gmx_fio_do_string(fio,buf);
2095 if (strncmp(buf,"VERSION",7))
2096 gmx_fatal(FARGS,"Can not read file %s,\n"
2097 " this file is from a Gromacs version which is older than 2.0\n"
2098 " Make a new one with grompp or use a gro or pdb file, if possible",
2099 gmx_fio_getname(fio));
2100 gmx_fio_do_int(fio,precision);
2101 bDouble = (precision == sizeof(double));
2102 if ((precision != sizeof(float)) && !bDouble)
2103 gmx_fatal(FARGS,"Unknown precision in file %s: real is %d bytes "
2104 "instead of %d or %d",
2105 gmx_fio_getname(fio),precision,sizeof(float),sizeof(double));
2106 gmx_fio_setprecision(fio,bDouble);
2107 fprintf(stderr,"Reading file %s, %s (%s precision)\n",
2108 gmx_fio_getname(fio),buf,bDouble ? "double" : "single");
2111 gmx_fio_write_string(fio,GromacsVersion());
2112 bDouble = (precision == sizeof(double));
2113 gmx_fio_setprecision(fio,bDouble);
2114 gmx_fio_do_int(fio,precision);
2116 sprintf(file_tag,"%s",tpx_tag);
2117 fgen = tpx_generation;
2120 /* Check versions! */
2121 gmx_fio_do_int(fio,fver);
2125 gmx_fio_do_string(fio,file_tag);
2131 /* Versions before 77 don't have the tag, set it to release */
2132 sprintf(file_tag,"%s",TPX_TAG_RELEASE);
2135 if (strcmp(file_tag,tpx_tag) != 0)
2137 fprintf(stderr,"Note: file tpx tag '%s', software tpx tag '%s'\n",
2140 /* We only support reading tpx files with the same tag as the code
2141 * or tpx files with the release tag and with lower version number.
2143 if (!(strcmp(file_tag,TPX_TAG_RELEASE) == 0 && fver < tpx_version))
2145 gmx_fatal(FARGS,"tpx tag/version mismatch: reading tpx file (%s) version %d, tag '%s' with program for tpx version %d, tag '%s'",
2146 gmx_fio_getname(fio),fver,file_tag,
2147 tpx_version,tpx_tag);
2154 gmx_fio_do_int(fio,fgen);
2161 if (file_version != NULL)
2163 *file_version = fver;
2165 if (file_generation != NULL)
2167 *file_generation = fgen;
2171 if ((fver <= tpx_incompatible_version) ||
2172 ((fver > tpx_version) && !TopOnlyOK) ||
2173 (fgen > tpx_generation))
2174 gmx_fatal(FARGS,"reading tpx file (%s) version %d with version %d program",
2175 gmx_fio_getname(fio),fver,tpx_version);
2177 do_section(fio,eitemHEADER,bRead);
2178 gmx_fio_do_int(fio,tpx->natoms);
2180 gmx_fio_do_int(fio,tpx->ngtc);
2184 gmx_fio_do_int(fio,idum);
2185 gmx_fio_do_real(fio,rdum);
2187 gmx_fio_do_real(fio,tpx->lambda);
2188 gmx_fio_do_int(fio,tpx->bIr);
2189 gmx_fio_do_int(fio,tpx->bTop);
2190 gmx_fio_do_int(fio,tpx->bX);
2191 gmx_fio_do_int(fio,tpx->bV);
2192 gmx_fio_do_int(fio,tpx->bF);
2193 gmx_fio_do_int(fio,tpx->bBox);
2195 if((fgen > tpx_generation)) {
2196 /* This can only happen if TopOnlyOK=TRUE */
2201 static int do_tpx(t_fileio *fio, gmx_bool bRead,
2202 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop,
2203 gmx_bool bXVallocated)
2208 gmx_bool TopOnlyOK,bDum=TRUE;
2209 int file_version,file_generation;
2213 gmx_bool bPeriodicMols;
2216 tpx.natoms = state->natoms;
2217 tpx.ngtc = state->ngtc;
2218 tpx.lambda = state->lambda;
2219 tpx.bIr = (ir != NULL);
2220 tpx.bTop = (mtop != NULL);
2221 tpx.bX = (state->x != NULL);
2222 tpx.bV = (state->v != NULL);
2223 tpx.bF = (f != NULL);
2227 TopOnlyOK = (ir==NULL);
2229 do_tpxheader(fio,bRead,&tpx,TopOnlyOK,&file_version,&file_generation);
2233 state->lambda = tpx.lambda;
2234 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
2238 init_state(state,0,tpx.ngtc,0,0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
2239 state->natoms = tpx.natoms;
2240 state->nalloc = tpx.natoms;
2244 init_state(state,tpx.natoms,tpx.ngtc,0,0); /* nose-hoover chains */
2248 #define do_test(fio,b,p) if (bRead && (p!=NULL) && !b) gmx_fatal(FARGS,"No %s in %s",#p,gmx_fio_getname(fio))
2250 do_test(fio,tpx.bBox,state->box);
2251 do_section(fio,eitemBOX,bRead);
2253 gmx_fio_ndo_rvec(fio,state->box,DIM);
2254 if (file_version >= 51) {
2255 gmx_fio_ndo_rvec(fio,state->box_rel,DIM);
2257 /* We initialize box_rel after reading the inputrec */
2258 clear_mat(state->box_rel);
2260 if (file_version >= 28) {
2261 gmx_fio_ndo_rvec(fio,state->boxv,DIM);
2262 if (file_version < 56) {
2264 gmx_fio_ndo_rvec(fio,mdum,DIM);
2269 if (state->ngtc > 0 && file_version >= 28) {
2271 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
2272 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
2273 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
2274 snew(dumv,state->ngtc);
2275 if (file_version < 69) {
2276 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2278 /* These used to be the Berendsen tcoupl_lambda's */
2279 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2283 /* Prior to tpx version 26, the inputrec was here.
2284 * I moved it to enable partial forward-compatibility
2285 * for analysis/viewer programs.
2287 if(file_version<26) {
2288 do_test(fio,tpx.bIr,ir);
2289 do_section(fio,eitemIR,bRead);
2292 do_inputrec(fio, ir,bRead,file_version,
2293 mtop ? &mtop->ffparams.fudgeQQ : NULL);
2295 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2298 do_inputrec(fio, &dum_ir,bRead,file_version,
2299 mtop ? &mtop->ffparams.fudgeQQ :NULL);
2301 pr_inputrec(debug,0,"inputrec",&dum_ir,FALSE);
2302 done_inputrec(&dum_ir);
2308 do_test(fio,tpx.bTop,mtop);
2309 do_section(fio,eitemTOP,bRead);
2312 do_mtop(fio,mtop,bRead, file_version);
2314 do_mtop(fio,&dum_top,bRead,file_version);
2315 done_mtop(&dum_top,TRUE);
2318 do_test(fio,tpx.bX,state->x);
2319 do_section(fio,eitemX,bRead);
2322 state->flags |= (1<<estX);
2324 gmx_fio_ndo_rvec(fio,state->x,state->natoms);
2327 do_test(fio,tpx.bV,state->v);
2328 do_section(fio,eitemV,bRead);
2331 state->flags |= (1<<estV);
2333 gmx_fio_ndo_rvec(fio,state->v,state->natoms);
2336 do_test(fio,tpx.bF,f);
2337 do_section(fio,eitemF,bRead);
2338 if (tpx.bF) gmx_fio_ndo_rvec(fio,f,state->natoms);
2340 /* Starting with tpx version 26, we have the inputrec
2341 * at the end of the file, so we can ignore it
2342 * if the file is never than the software (but still the
2343 * same generation - see comments at the top of this file.
2348 bPeriodicMols = FALSE;
2349 if (file_version >= 26) {
2350 do_test(fio,tpx.bIr,ir);
2351 do_section(fio,eitemIR,bRead);
2353 if (file_version >= 53) {
2354 /* Removed the pbc info from do_inputrec, since we always want it */
2357 bPeriodicMols = ir->bPeriodicMols;
2359 gmx_fio_do_int(fio,ePBC);
2360 gmx_fio_do_gmx_bool(fio,bPeriodicMols);
2362 if (file_generation <= tpx_generation && ir) {
2363 do_inputrec(fio, ir,bRead,file_version,mtop ? &mtop->ffparams.fudgeQQ : NULL);
2365 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2366 if (file_version < 51)
2367 set_box_rel(ir,state);
2368 if (file_version < 53) {
2370 bPeriodicMols = ir->bPeriodicMols;
2373 if (bRead && ir && file_version >= 53) {
2374 /* We need to do this after do_inputrec, since that initializes ir */
2376 ir->bPeriodicMols = bPeriodicMols;
2385 if (state->ngtc == 0)
2387 /* Reading old version without tcoupl state data: set it */
2388 init_gtc_state(state,ir->opts.ngtc,0,ir->opts.nhchainlength);
2390 if (tpx.bTop && mtop)
2392 if (file_version < 57)
2394 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
2396 ir->eDisre = edrSimple;
2400 ir->eDisre = edrNone;
2403 set_disres_npair(mtop);
2407 if (tpx.bTop && mtop)
2409 gmx_mtop_finalize(mtop);
2412 if (file_version >= 57)
2416 env = getenv("GMX_NOCHARGEGROUPS");
2419 sscanf(env,"%d",&ienv);
2420 fprintf(stderr,"\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
2425 "Will make single atomic charge groups in non-solvent%s\n",
2426 ienv > 1 ? " and solvent" : "");
2427 gmx_mtop_make_atomic_charge_groups(mtop,ienv==1);
2429 fprintf(stderr,"\n");
2437 /************************************************************
2439 * The following routines are the exported ones
2441 ************************************************************/
2443 t_fileio *open_tpx(const char *fn,const char *mode)
2445 return gmx_fio_open(fn,mode);
2448 void close_tpx(t_fileio *fio)
2453 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
2454 int *file_version, int *file_generation)
2458 fio = open_tpx(fn,"r");
2459 do_tpxheader(fio,TRUE,tpx,TopOnlyOK,file_version,file_generation);
2463 void write_tpx_state(const char *fn,
2464 t_inputrec *ir,t_state *state,gmx_mtop_t *mtop)
2468 fio = open_tpx(fn,"w");
2469 do_tpx(fio,FALSE,ir,state,NULL,mtop,FALSE);
2473 void read_tpx_state(const char *fn,
2474 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop)
2478 fio = open_tpx(fn,"r");
2479 do_tpx(fio,TRUE,ir,state,f,mtop,FALSE);
2483 int read_tpx(const char *fn,
2484 t_inputrec *ir, matrix box,int *natoms,
2485 rvec *x,rvec *v,rvec *f,gmx_mtop_t *mtop)
2493 fio = open_tpx(fn,"r");
2494 ePBC = do_tpx(fio,TRUE,ir,&state,f,mtop,TRUE);
2496 *natoms = state.natoms;
2498 copy_mat(state.box,box);
2506 int read_tpx_top(const char *fn,
2507 t_inputrec *ir, matrix box,int *natoms,
2508 rvec *x,rvec *v,rvec *f,t_topology *top)
2514 ePBC = read_tpx(fn,ir,box,natoms,x,v,f,&mtop);
2516 *top = gmx_mtop_t_to_t_topology(&mtop);
2521 gmx_bool fn2bTPX(const char *file)
2523 switch (fn2ftp(file)) {
2533 gmx_bool read_tps_conf(const char *infile,char *title,t_topology *top,int *ePBC,
2534 rvec **x,rvec **v,matrix box,gmx_bool bMass)
2537 int natoms,i,version,generation;
2538 gmx_bool bTop,bXNULL=FALSE;
2540 t_topology *topconv;
2543 bTop = fn2bTPX(infile);
2546 read_tpxheader(infile,&header,TRUE,&version,&generation);
2548 snew(*x,header.natoms);
2550 snew(*v,header.natoms);
2552 *ePBC = read_tpx(infile,NULL,box,&natoms,
2553 (x==NULL) ? NULL : *x,(v==NULL) ? NULL : *v,NULL,mtop);
2554 *top = gmx_mtop_t_to_t_topology(mtop);
2556 strcpy(title,*top->name);
2557 tpx_make_chain_identifiers(&top->atoms,&top->mols);
2560 get_stx_coordnum(infile,&natoms);
2561 init_t_atoms(&top->atoms,natoms,(fn2ftp(infile) == efPDB));
2570 read_stx_conf(infile,title,&top->atoms,*x,(v==NULL) ? NULL : *v,ePBC,box);
2577 aps = gmx_atomprop_init();
2578 for(i=0; (i<natoms); i++)
2579 if (!gmx_atomprop_query(aps,epropMass,
2580 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2581 *top->atoms.atomname[i],
2582 &(top->atoms.atom[i].m))) {
2584 fprintf(debug,"Can not find mass for atom %s %d %s, setting to 1\n",
2585 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2586 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
2587 *top->atoms.atomname[i]);
2589 gmx_atomprop_destroy(aps);
2591 top->idef.ntypes=-1;