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 = 76;
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 { 76, F_LINEAR_ANGLES },
130 { 30, F_CROSS_BOND_BONDS },
131 { 30, F_CROSS_BOND_ANGLES },
132 { 30, F_UREY_BRADLEY },
133 { 34, F_QUARTIC_ANGLES },
143 { 72, F_NPSOLVATION },
145 { 41, F_LJC_PAIRS_NB },
148 { 32, F_COUL_RECIP },
150 { 30, F_POLARIZATION },
152 { 22, F_DISRESVIOL },
156 { 26, F_DIHRESVIOL },
161 { 46, F_ECONSERVED },
167 #define NFTUPD asize(ftupd)
169 /* Needed for backward compatibility */
172 static void _do_section(t_fileio *fio,int key,gmx_bool bRead,const char *src,
178 if (gmx_fio_getftp(fio) == efTPA) {
180 gmx_fio_write_string(fio,itemstr[key]);
181 bDbg = gmx_fio_getdebug(fio);
182 gmx_fio_setdebug(fio,FALSE);
183 gmx_fio_write_string(fio,comment_str[key]);
184 gmx_fio_setdebug(fio,bDbg);
187 if (gmx_fio_getdebug(fio))
188 fprintf(stderr,"Looking for section %s (%s, %d)",
189 itemstr[key],src,line);
192 gmx_fio_do_string(fio,buf);
193 } while ((gmx_strcasecmp(buf,itemstr[key]) != 0));
195 if (gmx_strcasecmp(buf,itemstr[key]) != 0)
196 gmx_fatal(FARGS,"\nCould not find section heading %s",itemstr[key]);
197 else if (gmx_fio_getdebug(fio))
198 fprintf(stderr," and found it\n");
203 #define do_section(fio,key,bRead) _do_section(fio,key,bRead,__FILE__,__LINE__)
205 /**************************************************************
207 * Now the higer level routines that do io of the structures and arrays
209 **************************************************************/
210 static void do_pullgrp(t_fileio *fio, t_pullgrp *pgrp, gmx_bool bRead,
216 gmx_fio_do_int(fio,pgrp->nat);
218 snew(pgrp->ind,pgrp->nat);
219 bDum=gmx_fio_ndo_int(fio,pgrp->ind,pgrp->nat);
220 gmx_fio_do_int(fio,pgrp->nweight);
222 snew(pgrp->weight,pgrp->nweight);
223 bDum=gmx_fio_ndo_real(fio,pgrp->weight,pgrp->nweight);
224 gmx_fio_do_int(fio,pgrp->pbcatom);
225 gmx_fio_do_rvec(fio,pgrp->vec);
226 gmx_fio_do_rvec(fio,pgrp->init);
227 gmx_fio_do_real(fio,pgrp->rate);
228 gmx_fio_do_real(fio,pgrp->k);
229 if (file_version >= 56) {
230 gmx_fio_do_real(fio,pgrp->kB);
236 static void do_pull(t_fileio *fio, t_pull *pull,gmx_bool bRead, int file_version)
240 gmx_fio_do_int(fio,pull->ngrp);
241 gmx_fio_do_int(fio,pull->eGeom);
242 gmx_fio_do_ivec(fio,pull->dim);
243 gmx_fio_do_real(fio,pull->cyl_r1);
244 gmx_fio_do_real(fio,pull->cyl_r0);
245 gmx_fio_do_real(fio,pull->constr_tol);
246 gmx_fio_do_int(fio,pull->nstxout);
247 gmx_fio_do_int(fio,pull->nstfout);
249 snew(pull->grp,pull->ngrp+1);
250 for(g=0; g<pull->ngrp+1; g++)
251 do_pullgrp(fio,&pull->grp[g],bRead,file_version);
255 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg,gmx_bool bRead, int file_version)
260 gmx_fio_do_int(fio,rotg->eType);
261 gmx_fio_do_int(fio,rotg->bMassW);
262 gmx_fio_do_int(fio,rotg->nat);
264 snew(rotg->ind,rotg->nat);
265 gmx_fio_ndo_int(fio,rotg->ind,rotg->nat);
267 snew(rotg->x_ref,rotg->nat);
268 gmx_fio_ndo_rvec(fio,rotg->x_ref,rotg->nat);
269 gmx_fio_do_rvec(fio,rotg->vec);
270 gmx_fio_do_rvec(fio,rotg->pivot);
271 gmx_fio_do_real(fio,rotg->rate);
272 gmx_fio_do_real(fio,rotg->k);
273 gmx_fio_do_real(fio,rotg->slab_dist);
274 gmx_fio_do_real(fio,rotg->min_gaussian);
275 gmx_fio_do_real(fio,rotg->eps);
276 gmx_fio_do_int(fio,rotg->eFittype);
277 gmx_fio_do_int(fio,rotg->PotAngle_nstep);
278 gmx_fio_do_real(fio,rotg->PotAngle_step);
281 static void do_rot(t_fileio *fio, t_rot *rot,gmx_bool bRead, int file_version)
285 gmx_fio_do_int(fio,rot->ngrp);
286 gmx_fio_do_int(fio,rot->nstrout);
287 gmx_fio_do_int(fio,rot->nstsout);
289 snew(rot->grp,rot->ngrp);
290 for(g=0; g<rot->ngrp; g++)
291 do_rotgrp(fio, &rot->grp[g],bRead,file_version);
295 static void do_inputrec(t_fileio *fio, t_inputrec *ir,gmx_bool bRead,
296 int file_version, real *fudgeQQ)
298 int i,j,k,*tmp,idum=0;
303 real zerotemptime,finish_t,init_temp,finish_temp;
305 if (file_version != tpx_version)
307 /* Give a warning about features that are not accessible */
308 fprintf(stderr,"Note: tpx file_version %d, software version %d\n",
309 file_version,tpx_version);
317 if (file_version == 0)
322 /* Basic inputrec stuff */
323 gmx_fio_do_int(fio,ir->eI);
324 if (file_version >= 62) {
325 gmx_fio_do_gmx_large_int(fio, ir->nsteps);
327 gmx_fio_do_int(fio,idum);
330 if(file_version > 25) {
331 if (file_version >= 62) {
332 gmx_fio_do_gmx_large_int(fio, ir->init_step);
334 gmx_fio_do_int(fio,idum);
335 ir->init_step = idum;
341 if(file_version >= 58)
342 gmx_fio_do_int(fio,ir->simulation_part);
344 ir->simulation_part=1;
346 if (file_version >= 67) {
347 gmx_fio_do_int(fio,ir->nstcalcenergy);
349 ir->nstcalcenergy = 1;
351 if (file_version < 53) {
352 /* The pbc info has been moved out of do_inputrec,
353 * since we always want it, also without reading the inputrec.
355 gmx_fio_do_int(fio,ir->ePBC);
356 if ((file_version <= 15) && (ir->ePBC == 2))
358 if (file_version >= 45) {
359 gmx_fio_do_int(fio,ir->bPeriodicMols);
363 ir->bPeriodicMols = TRUE;
365 ir->bPeriodicMols = FALSE;
369 gmx_fio_do_int(fio,ir->ns_type);
370 gmx_fio_do_int(fio,ir->nstlist);
371 gmx_fio_do_int(fio,ir->ndelta);
372 if (file_version < 41) {
373 gmx_fio_do_int(fio,idum);
374 gmx_fio_do_int(fio,idum);
376 if (file_version >= 45)
377 gmx_fio_do_real(fio,ir->rtpi);
380 gmx_fio_do_int(fio,ir->nstcomm);
381 if (file_version > 34)
382 gmx_fio_do_int(fio,ir->comm_mode);
383 else if (ir->nstcomm < 0)
384 ir->comm_mode = ecmANGULAR;
386 ir->comm_mode = ecmLINEAR;
387 ir->nstcomm = abs(ir->nstcomm);
389 if(file_version > 25)
390 gmx_fio_do_int(fio,ir->nstcheckpoint);
394 gmx_fio_do_int(fio,ir->nstcgsteep);
397 gmx_fio_do_int(fio,ir->nbfgscorr);
401 gmx_fio_do_int(fio,ir->nstlog);
402 gmx_fio_do_int(fio,ir->nstxout);
403 gmx_fio_do_int(fio,ir->nstvout);
404 gmx_fio_do_int(fio,ir->nstfout);
405 gmx_fio_do_int(fio,ir->nstenergy);
406 gmx_fio_do_int(fio,ir->nstxtcout);
407 if (file_version >= 59) {
408 gmx_fio_do_double(fio,ir->init_t);
409 gmx_fio_do_double(fio,ir->delta_t);
411 gmx_fio_do_real(fio,rdum);
413 gmx_fio_do_real(fio,rdum);
416 gmx_fio_do_real(fio,ir->xtcprec);
417 if (file_version < 19) {
418 gmx_fio_do_int(fio,idum);
419 gmx_fio_do_int(fio,idum);
421 if(file_version < 18)
422 gmx_fio_do_int(fio,idum);
423 gmx_fio_do_real(fio,ir->rlist);
424 if (file_version >= 67) {
425 gmx_fio_do_real(fio,ir->rlistlong);
427 gmx_fio_do_int(fio,ir->coulombtype);
428 if (file_version < 32 && ir->coulombtype == eelRF)
429 ir->coulombtype = eelRF_NEC;
430 gmx_fio_do_real(fio,ir->rcoulomb_switch);
431 gmx_fio_do_real(fio,ir->rcoulomb);
432 gmx_fio_do_int(fio,ir->vdwtype);
433 gmx_fio_do_real(fio,ir->rvdw_switch);
434 gmx_fio_do_real(fio,ir->rvdw);
435 if (file_version < 67) {
436 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
438 gmx_fio_do_int(fio,ir->eDispCorr);
439 gmx_fio_do_real(fio,ir->epsilon_r);
440 if (file_version >= 37) {
441 gmx_fio_do_real(fio,ir->epsilon_rf);
443 if (EEL_RF(ir->coulombtype)) {
444 ir->epsilon_rf = ir->epsilon_r;
447 ir->epsilon_rf = 1.0;
450 if (file_version >= 29)
451 gmx_fio_do_real(fio,ir->tabext);
455 if(file_version > 25) {
456 gmx_fio_do_int(fio,ir->gb_algorithm);
457 gmx_fio_do_int(fio,ir->nstgbradii);
458 gmx_fio_do_real(fio,ir->rgbradii);
459 gmx_fio_do_real(fio,ir->gb_saltconc);
460 gmx_fio_do_int(fio,ir->implicit_solvent);
462 ir->gb_algorithm=egbSTILL;
466 ir->implicit_solvent=eisNO;
470 gmx_fio_do_real(fio,ir->gb_epsilon_solvent);
471 gmx_fio_do_real(fio,ir->gb_obc_alpha);
472 gmx_fio_do_real(fio,ir->gb_obc_beta);
473 gmx_fio_do_real(fio,ir->gb_obc_gamma);
476 gmx_fio_do_real(fio,ir->gb_dielectric_offset);
477 gmx_fio_do_int(fio,ir->sa_algorithm);
481 ir->gb_dielectric_offset = 0.009;
482 ir->sa_algorithm = esaAPPROX;
484 gmx_fio_do_real(fio,ir->sa_surface_tension);
486 /* Override sa_surface_tension if it is not changed in the mpd-file */
487 if(ir->sa_surface_tension<0)
489 if(ir->gb_algorithm==egbSTILL)
491 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
493 else if(ir->gb_algorithm==egbHCT || ir->gb_algorithm==egbOBC)
495 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
502 /* Better use sensible values than insane (0.0) ones... */
503 ir->gb_epsilon_solvent = 80;
504 ir->gb_obc_alpha = 1.0;
505 ir->gb_obc_beta = 0.8;
506 ir->gb_obc_gamma = 4.85;
507 ir->sa_surface_tension = 2.092;
511 gmx_fio_do_int(fio,ir->nkx);
512 gmx_fio_do_int(fio,ir->nky);
513 gmx_fio_do_int(fio,ir->nkz);
514 gmx_fio_do_int(fio,ir->pme_order);
515 gmx_fio_do_real(fio,ir->ewald_rtol);
517 if (file_version >=24)
518 gmx_fio_do_int(fio,ir->ewald_geometry);
520 ir->ewald_geometry=eewg3D;
522 if (file_version <=17) {
523 ir->epsilon_surface=0;
524 if (file_version==17)
525 gmx_fio_do_int(fio,idum);
528 gmx_fio_do_real(fio,ir->epsilon_surface);
530 gmx_fio_do_gmx_bool(fio,ir->bOptFFT);
532 gmx_fio_do_gmx_bool(fio,ir->bContinuation);
533 gmx_fio_do_int(fio,ir->etc);
534 /* before version 18, ir->etc was a gmx_bool (ir->btc),
535 * but the values 0 and 1 still mean no and
536 * berendsen temperature coupling, respectively.
538 if (file_version >= 71)
540 gmx_fio_do_int(fio,ir->nsttcouple);
544 ir->nsttcouple = ir->nstcalcenergy;
546 if (file_version <= 15)
548 gmx_fio_do_int(fio,idum);
550 if (file_version <=17)
552 gmx_fio_do_int(fio,ir->epct);
553 if (file_version <= 15)
557 ir->epct = epctSURFACETENSION;
559 gmx_fio_do_int(fio,idum);
562 /* we have removed the NO alternative at the beginning */
566 ir->epct=epctISOTROPIC;
570 ir->epc=epcBERENDSEN;
575 gmx_fio_do_int(fio,ir->epc);
576 gmx_fio_do_int(fio,ir->epct);
578 if (file_version >= 71)
580 gmx_fio_do_int(fio,ir->nstpcouple);
584 ir->nstpcouple = ir->nstcalcenergy;
586 gmx_fio_do_real(fio,ir->tau_p);
587 if (file_version <= 15) {
588 gmx_fio_do_rvec(fio,vdum);
589 clear_mat(ir->ref_p);
591 ir->ref_p[i][i] = vdum[i];
593 gmx_fio_do_rvec(fio,ir->ref_p[XX]);
594 gmx_fio_do_rvec(fio,ir->ref_p[YY]);
595 gmx_fio_do_rvec(fio,ir->ref_p[ZZ]);
597 if (file_version <= 15) {
598 gmx_fio_do_rvec(fio,vdum);
599 clear_mat(ir->compress);
601 ir->compress[i][i] = vdum[i];
604 gmx_fio_do_rvec(fio,ir->compress[XX]);
605 gmx_fio_do_rvec(fio,ir->compress[YY]);
606 gmx_fio_do_rvec(fio,ir->compress[ZZ]);
608 if (file_version >= 47) {
609 gmx_fio_do_int(fio,ir->refcoord_scaling);
610 gmx_fio_do_rvec(fio,ir->posres_com);
611 gmx_fio_do_rvec(fio,ir->posres_comB);
613 ir->refcoord_scaling = erscNO;
614 clear_rvec(ir->posres_com);
615 clear_rvec(ir->posres_comB);
617 if(file_version > 25)
618 gmx_fio_do_int(fio,ir->andersen_seed);
622 if(file_version < 26) {
623 gmx_fio_do_gmx_bool(fio,bSimAnn);
624 gmx_fio_do_real(fio,zerotemptime);
627 if (file_version < 37)
628 gmx_fio_do_real(fio,rdum);
630 gmx_fio_do_real(fio,ir->shake_tol);
631 if (file_version < 54)
632 gmx_fio_do_real(fio,*fudgeQQ);
633 gmx_fio_do_int(fio,ir->efep);
634 if (file_version <= 14 && ir->efep > efepNO)
636 if (file_version >= 59) {
637 gmx_fio_do_double(fio, ir->init_lambda);
638 gmx_fio_do_double(fio, ir->delta_lambda);
640 gmx_fio_do_real(fio,rdum);
641 ir->init_lambda = rdum;
642 gmx_fio_do_real(fio,rdum);
643 ir->delta_lambda = rdum;
645 if (file_version >= 64) {
646 gmx_fio_do_int(fio,ir->n_flambda);
648 snew(ir->flambda,ir->n_flambda);
650 bDum=gmx_fio_ndo_double(fio,ir->flambda,ir->n_flambda);
655 if (file_version >= 13)
656 gmx_fio_do_real(fio,ir->sc_alpha);
659 if (file_version >= 38)
660 gmx_fio_do_int(fio,ir->sc_power);
663 if (file_version >= 15)
664 gmx_fio_do_real(fio,ir->sc_sigma);
669 if (file_version >= 71)
671 ir->sc_sigma_min = ir->sc_sigma;
675 ir->sc_sigma_min = 0;
678 if (file_version >= 64) {
679 gmx_fio_do_int(fio,ir->nstdhdl);
684 if (file_version >= 73)
686 gmx_fio_do_int(fio, ir->separate_dhdl_file);
687 gmx_fio_do_int(fio, ir->dhdl_derivatives);
691 ir->separate_dhdl_file = sepdhdlfileYES;
692 ir->dhdl_derivatives = dhdlderivativesYES;
695 if (file_version >= 71)
697 gmx_fio_do_int(fio,ir->dh_hist_size);
698 gmx_fio_do_double(fio,ir->dh_hist_spacing);
702 ir->dh_hist_size = 0;
703 ir->dh_hist_spacing = 0.1;
705 if (file_version >= 57) {
706 gmx_fio_do_int(fio,ir->eDisre);
708 gmx_fio_do_int(fio,ir->eDisreWeighting);
709 if (file_version < 22) {
710 if (ir->eDisreWeighting == 0)
711 ir->eDisreWeighting = edrwEqual;
713 ir->eDisreWeighting = edrwConservative;
715 gmx_fio_do_gmx_bool(fio,ir->bDisreMixed);
716 gmx_fio_do_real(fio,ir->dr_fc);
717 gmx_fio_do_real(fio,ir->dr_tau);
718 gmx_fio_do_int(fio,ir->nstdisreout);
719 if (file_version >= 22) {
720 gmx_fio_do_real(fio,ir->orires_fc);
721 gmx_fio_do_real(fio,ir->orires_tau);
722 gmx_fio_do_int(fio,ir->nstorireout);
728 if(file_version >= 26) {
729 gmx_fio_do_real(fio,ir->dihre_fc);
730 if (file_version < 56) {
731 gmx_fio_do_real(fio,rdum);
732 gmx_fio_do_int(fio,idum);
738 gmx_fio_do_real(fio,ir->em_stepsize);
739 gmx_fio_do_real(fio,ir->em_tol);
740 if (file_version >= 22)
741 gmx_fio_do_gmx_bool(fio,ir->bShakeSOR);
743 ir->bShakeSOR = TRUE;
744 if (file_version >= 11)
745 gmx_fio_do_int(fio,ir->niter);
748 fprintf(stderr,"Note: niter not in run input file, setting it to %d\n",
751 if (file_version >= 21)
752 gmx_fio_do_real(fio,ir->fc_stepsize);
755 gmx_fio_do_int(fio,ir->eConstrAlg);
756 gmx_fio_do_int(fio,ir->nProjOrder);
757 gmx_fio_do_real(fio,ir->LincsWarnAngle);
758 if (file_version <= 14)
759 gmx_fio_do_int(fio,idum);
760 if (file_version >=26)
761 gmx_fio_do_int(fio,ir->nLincsIter);
764 fprintf(stderr,"Note: nLincsIter not in run input file, setting it to %d\n",
767 if (file_version < 33)
768 gmx_fio_do_real(fio,bd_temp);
769 gmx_fio_do_real(fio,ir->bd_fric);
770 gmx_fio_do_int(fio,ir->ld_seed);
771 if (file_version >= 33) {
773 gmx_fio_do_rvec(fio,ir->deform[i]);
776 clear_rvec(ir->deform[i]);
778 if (file_version >= 14)
779 gmx_fio_do_real(fio,ir->cos_accel);
782 gmx_fio_do_int(fio,ir->userint1);
783 gmx_fio_do_int(fio,ir->userint2);
784 gmx_fio_do_int(fio,ir->userint3);
785 gmx_fio_do_int(fio,ir->userint4);
786 gmx_fio_do_real(fio,ir->userreal1);
787 gmx_fio_do_real(fio,ir->userreal2);
788 gmx_fio_do_real(fio,ir->userreal3);
789 gmx_fio_do_real(fio,ir->userreal4);
792 if (file_version >= 75) {
793 gmx_fio_do_gmx_bool(fio,ir->bAdress);
795 if (bRead) snew(ir->adress, 1);
796 gmx_fio_do_int(fio,ir->adress->type);
797 gmx_fio_do_real(fio,ir->adress->const_wf);
798 gmx_fio_do_real(fio,ir->adress->ex_width);
799 gmx_fio_do_real(fio,ir->adress->hy_width);
800 gmx_fio_do_int(fio,ir->adress->icor);
801 gmx_fio_do_int(fio,ir->adress->site);
802 gmx_fio_do_rvec(fio,ir->adress->refs);
803 gmx_fio_do_int(fio,ir->adress->n_tf_grps);
804 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
805 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
806 gmx_fio_do_int(fio,ir->adress->do_hybridpairs);
808 if (bRead)snew(ir->adress->tf_table_index,ir->adress->n_tf_grps);
809 if (ir->adress->n_tf_grps > 0) {
810 bDum=gmx_fio_ndo_int(fio,ir->adress->tf_table_index,ir->adress->n_tf_grps);
812 if (bRead)snew(ir->adress->group_explicit,ir->adress->n_energy_grps);
813 if (ir->adress->n_energy_grps > 0) {
814 bDum=gmx_fio_ndo_int(fio, ir->adress->group_explicit,ir->adress->n_energy_grps);
822 if (file_version >= 48) {
823 gmx_fio_do_int(fio,ir->ePull);
824 if (ir->ePull != epullNO) {
827 do_pull(fio, ir->pull,bRead,file_version);
833 /* Enforced rotation */
834 if (file_version >= 74) {
835 gmx_fio_do_int(fio,ir->bRot);
836 if (ir->bRot == TRUE) {
839 do_rot(fio, ir->rot,bRead,file_version);
846 gmx_fio_do_int(fio,ir->opts.ngtc);
847 if (file_version >= 69) {
848 gmx_fio_do_int(fio,ir->opts.nhchainlength);
850 ir->opts.nhchainlength = 1;
852 gmx_fio_do_int(fio,ir->opts.ngacc);
853 gmx_fio_do_int(fio,ir->opts.ngfrz);
854 gmx_fio_do_int(fio,ir->opts.ngener);
857 snew(ir->opts.nrdf, ir->opts.ngtc);
858 snew(ir->opts.ref_t, ir->opts.ngtc);
859 snew(ir->opts.annealing, ir->opts.ngtc);
860 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
861 snew(ir->opts.anneal_time, ir->opts.ngtc);
862 snew(ir->opts.anneal_temp, ir->opts.ngtc);
863 snew(ir->opts.tau_t, ir->opts.ngtc);
864 snew(ir->opts.nFreeze,ir->opts.ngfrz);
865 snew(ir->opts.acc, ir->opts.ngacc);
866 snew(ir->opts.egp_flags,ir->opts.ngener*ir->opts.ngener);
868 if (ir->opts.ngtc > 0) {
869 if (bRead && file_version<13) {
870 snew(tmp,ir->opts.ngtc);
871 bDum=gmx_fio_ndo_int(fio,tmp, ir->opts.ngtc);
872 for(i=0; i<ir->opts.ngtc; i++)
873 ir->opts.nrdf[i] = tmp[i];
876 bDum=gmx_fio_ndo_real(fio,ir->opts.nrdf, ir->opts.ngtc);
878 bDum=gmx_fio_ndo_real(fio,ir->opts.ref_t,ir->opts.ngtc);
879 bDum=gmx_fio_ndo_real(fio,ir->opts.tau_t,ir->opts.ngtc);
880 if (file_version<33 && ir->eI==eiBD) {
881 for(i=0; i<ir->opts.ngtc; i++)
882 ir->opts.tau_t[i] = bd_temp;
885 if (ir->opts.ngfrz > 0)
886 bDum=gmx_fio_ndo_ivec(fio,ir->opts.nFreeze,ir->opts.ngfrz);
887 if (ir->opts.ngacc > 0)
888 gmx_fio_ndo_rvec(fio,ir->opts.acc,ir->opts.ngacc);
889 if (file_version >= 12)
890 bDum=gmx_fio_ndo_int(fio,ir->opts.egp_flags,
891 ir->opts.ngener*ir->opts.ngener);
893 if(bRead && file_version < 26) {
894 for(i=0;i<ir->opts.ngtc;i++) {
896 ir->opts.annealing[i] = eannSINGLE;
897 ir->opts.anneal_npoints[i] = 2;
898 snew(ir->opts.anneal_time[i],2);
899 snew(ir->opts.anneal_temp[i],2);
900 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
901 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
902 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
903 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
904 ir->opts.anneal_time[i][0] = ir->init_t;
905 ir->opts.anneal_time[i][1] = finish_t;
906 ir->opts.anneal_temp[i][0] = init_temp;
907 ir->opts.anneal_temp[i][1] = finish_temp;
909 ir->opts.annealing[i] = eannNO;
910 ir->opts.anneal_npoints[i] = 0;
914 /* file version 26 or later */
915 /* First read the lists with annealing and npoints for each group */
916 bDum=gmx_fio_ndo_int(fio,ir->opts.annealing,ir->opts.ngtc);
917 bDum=gmx_fio_ndo_int(fio,ir->opts.anneal_npoints,ir->opts.ngtc);
918 for(j=0;j<(ir->opts.ngtc);j++) {
919 k=ir->opts.anneal_npoints[j];
921 snew(ir->opts.anneal_time[j],k);
922 snew(ir->opts.anneal_temp[j],k);
924 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_time[j],k);
925 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_temp[j],k);
929 if (file_version >= 45) {
930 gmx_fio_do_int(fio,ir->nwall);
931 gmx_fio_do_int(fio,ir->wall_type);
932 if (file_version >= 50)
933 gmx_fio_do_real(fio,ir->wall_r_linpot);
935 ir->wall_r_linpot = -1;
936 gmx_fio_do_int(fio,ir->wall_atomtype[0]);
937 gmx_fio_do_int(fio,ir->wall_atomtype[1]);
938 gmx_fio_do_real(fio,ir->wall_density[0]);
939 gmx_fio_do_real(fio,ir->wall_density[1]);
940 gmx_fio_do_real(fio,ir->wall_ewald_zfac);
944 ir->wall_atomtype[0] = -1;
945 ir->wall_atomtype[1] = -1;
946 ir->wall_density[0] = 0;
947 ir->wall_density[1] = 0;
948 ir->wall_ewald_zfac = 3;
950 /* Cosine stuff for electric fields */
951 for(j=0; (j<DIM); j++) {
952 gmx_fio_do_int(fio,ir->ex[j].n);
953 gmx_fio_do_int(fio,ir->et[j].n);
955 snew(ir->ex[j].a, ir->ex[j].n);
956 snew(ir->ex[j].phi,ir->ex[j].n);
957 snew(ir->et[j].a, ir->et[j].n);
958 snew(ir->et[j].phi,ir->et[j].n);
960 bDum=gmx_fio_ndo_real(fio,ir->ex[j].a, ir->ex[j].n);
961 bDum=gmx_fio_ndo_real(fio,ir->ex[j].phi,ir->ex[j].n);
962 bDum=gmx_fio_ndo_real(fio,ir->et[j].a, ir->et[j].n);
963 bDum=gmx_fio_ndo_real(fio,ir->et[j].phi,ir->et[j].n);
967 if(file_version>=39){
968 gmx_fio_do_gmx_bool(fio,ir->bQMMM);
969 gmx_fio_do_int(fio,ir->QMMMscheme);
970 gmx_fio_do_real(fio,ir->scalefactor);
971 gmx_fio_do_int(fio,ir->opts.ngQM);
973 snew(ir->opts.QMmethod, ir->opts.ngQM);
974 snew(ir->opts.QMbasis, ir->opts.ngQM);
975 snew(ir->opts.QMcharge, ir->opts.ngQM);
976 snew(ir->opts.QMmult, ir->opts.ngQM);
977 snew(ir->opts.bSH, ir->opts.ngQM);
978 snew(ir->opts.CASorbitals, ir->opts.ngQM);
979 snew(ir->opts.CASelectrons,ir->opts.ngQM);
980 snew(ir->opts.SAon, ir->opts.ngQM);
981 snew(ir->opts.SAoff, ir->opts.ngQM);
982 snew(ir->opts.SAsteps, ir->opts.ngQM);
983 snew(ir->opts.bOPT, ir->opts.ngQM);
984 snew(ir->opts.bTS, ir->opts.ngQM);
986 if (ir->opts.ngQM > 0) {
987 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmethod,ir->opts.ngQM);
988 bDum=gmx_fio_ndo_int(fio,ir->opts.QMbasis,ir->opts.ngQM);
989 bDum=gmx_fio_ndo_int(fio,ir->opts.QMcharge,ir->opts.ngQM);
990 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmult,ir->opts.ngQM);
991 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bSH,ir->opts.ngQM);
992 bDum=gmx_fio_ndo_int(fio,ir->opts.CASorbitals,ir->opts.ngQM);
993 bDum=gmx_fio_ndo_int(fio,ir->opts.CASelectrons,ir->opts.ngQM);
994 bDum=gmx_fio_ndo_real(fio,ir->opts.SAon,ir->opts.ngQM);
995 bDum=gmx_fio_ndo_real(fio,ir->opts.SAoff,ir->opts.ngQM);
996 bDum=gmx_fio_ndo_int(fio,ir->opts.SAsteps,ir->opts.ngQM);
997 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bOPT,ir->opts.ngQM);
998 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bTS,ir->opts.ngQM);
1000 /* end of QMMM stuff */
1005 static void do_harm(t_fileio *fio, t_iparams *iparams,gmx_bool bRead)
1007 gmx_fio_do_real(fio,iparams->harmonic.rA);
1008 gmx_fio_do_real(fio,iparams->harmonic.krA);
1009 gmx_fio_do_real(fio,iparams->harmonic.rB);
1010 gmx_fio_do_real(fio,iparams->harmonic.krB);
1013 void do_iparams(t_fileio *fio, t_functype ftype,t_iparams *iparams,
1014 gmx_bool bRead, int file_version)
1021 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1029 do_harm(fio, iparams,bRead);
1030 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead) {
1031 /* Correct incorrect storage of parameters */
1032 iparams->pdihs.phiB = iparams->pdihs.phiA;
1033 iparams->pdihs.cpB = iparams->pdihs.cpA;
1036 case F_LINEAR_ANGLES:
1037 gmx_fio_do_real(fio,iparams->linangle.klinA);
1038 gmx_fio_do_real(fio,iparams->linangle.aA);
1039 gmx_fio_do_real(fio,iparams->linangle.klinB);
1040 gmx_fio_do_real(fio,iparams->linangle.aB);
1043 gmx_fio_do_real(fio,iparams->fene.bm);
1044 gmx_fio_do_real(fio,iparams->fene.kb);
1047 gmx_fio_do_real(fio,iparams->restraint.lowA);
1048 gmx_fio_do_real(fio,iparams->restraint.up1A);
1049 gmx_fio_do_real(fio,iparams->restraint.up2A);
1050 gmx_fio_do_real(fio,iparams->restraint.kA);
1051 gmx_fio_do_real(fio,iparams->restraint.lowB);
1052 gmx_fio_do_real(fio,iparams->restraint.up1B);
1053 gmx_fio_do_real(fio,iparams->restraint.up2B);
1054 gmx_fio_do_real(fio,iparams->restraint.kB);
1060 gmx_fio_do_real(fio,iparams->tab.kA);
1061 gmx_fio_do_int(fio,iparams->tab.table);
1062 gmx_fio_do_real(fio,iparams->tab.kB);
1064 case F_CROSS_BOND_BONDS:
1065 gmx_fio_do_real(fio,iparams->cross_bb.r1e);
1066 gmx_fio_do_real(fio,iparams->cross_bb.r2e);
1067 gmx_fio_do_real(fio,iparams->cross_bb.krr);
1069 case F_CROSS_BOND_ANGLES:
1070 gmx_fio_do_real(fio,iparams->cross_ba.r1e);
1071 gmx_fio_do_real(fio,iparams->cross_ba.r2e);
1072 gmx_fio_do_real(fio,iparams->cross_ba.r3e);
1073 gmx_fio_do_real(fio,iparams->cross_ba.krt);
1075 case F_UREY_BRADLEY:
1076 gmx_fio_do_real(fio,iparams->u_b.theta);
1077 gmx_fio_do_real(fio,iparams->u_b.ktheta);
1078 gmx_fio_do_real(fio,iparams->u_b.r13);
1079 gmx_fio_do_real(fio,iparams->u_b.kUB);
1081 case F_QUARTIC_ANGLES:
1082 gmx_fio_do_real(fio,iparams->qangle.theta);
1083 bDum=gmx_fio_ndo_real(fio,iparams->qangle.c,5);
1086 gmx_fio_do_real(fio,iparams->bham.a);
1087 gmx_fio_do_real(fio,iparams->bham.b);
1088 gmx_fio_do_real(fio,iparams->bham.c);
1091 gmx_fio_do_real(fio,iparams->morse.b0);
1092 gmx_fio_do_real(fio,iparams->morse.cb);
1093 gmx_fio_do_real(fio,iparams->morse.beta);
1096 gmx_fio_do_real(fio,iparams->cubic.b0);
1097 gmx_fio_do_real(fio,iparams->cubic.kb);
1098 gmx_fio_do_real(fio,iparams->cubic.kcub);
1102 case F_POLARIZATION:
1103 gmx_fio_do_real(fio,iparams->polarize.alpha);
1106 gmx_fio_do_real(fio,iparams->anharm_polarize.alpha);
1107 gmx_fio_do_real(fio,iparams->anharm_polarize.drcut);
1108 gmx_fio_do_real(fio,iparams->anharm_polarize.khyp);
1111 if (file_version < 31)
1112 gmx_fatal(FARGS,"Old tpr files with water_polarization not supported. Make a new.");
1113 gmx_fio_do_real(fio,iparams->wpol.al_x);
1114 gmx_fio_do_real(fio,iparams->wpol.al_y);
1115 gmx_fio_do_real(fio,iparams->wpol.al_z);
1116 gmx_fio_do_real(fio,iparams->wpol.rOH);
1117 gmx_fio_do_real(fio,iparams->wpol.rHH);
1118 gmx_fio_do_real(fio,iparams->wpol.rOD);
1121 gmx_fio_do_real(fio,iparams->thole.a);
1122 gmx_fio_do_real(fio,iparams->thole.alpha1);
1123 gmx_fio_do_real(fio,iparams->thole.alpha2);
1124 gmx_fio_do_real(fio,iparams->thole.rfac);
1127 gmx_fio_do_real(fio,iparams->lj.c6);
1128 gmx_fio_do_real(fio,iparams->lj.c12);
1131 gmx_fio_do_real(fio,iparams->lj14.c6A);
1132 gmx_fio_do_real(fio,iparams->lj14.c12A);
1133 gmx_fio_do_real(fio,iparams->lj14.c6B);
1134 gmx_fio_do_real(fio,iparams->lj14.c12B);
1137 gmx_fio_do_real(fio,iparams->ljc14.fqq);
1138 gmx_fio_do_real(fio,iparams->ljc14.qi);
1139 gmx_fio_do_real(fio,iparams->ljc14.qj);
1140 gmx_fio_do_real(fio,iparams->ljc14.c6);
1141 gmx_fio_do_real(fio,iparams->ljc14.c12);
1143 case F_LJC_PAIRS_NB:
1144 gmx_fio_do_real(fio,iparams->ljcnb.qi);
1145 gmx_fio_do_real(fio,iparams->ljcnb.qj);
1146 gmx_fio_do_real(fio,iparams->ljcnb.c6);
1147 gmx_fio_do_real(fio,iparams->ljcnb.c12);
1153 gmx_fio_do_real(fio,iparams->pdihs.phiA);
1154 gmx_fio_do_real(fio,iparams->pdihs.cpA);
1155 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42) {
1156 /* Read the incorrectly stored multiplicity */
1157 gmx_fio_do_real(fio,iparams->harmonic.rB);
1158 gmx_fio_do_real(fio,iparams->harmonic.krB);
1159 iparams->pdihs.phiB = iparams->pdihs.phiA;
1160 iparams->pdihs.cpB = iparams->pdihs.cpA;
1162 gmx_fio_do_real(fio,iparams->pdihs.phiB);
1163 gmx_fio_do_real(fio,iparams->pdihs.cpB);
1164 gmx_fio_do_int(fio,iparams->pdihs.mult);
1168 gmx_fio_do_int(fio,iparams->disres.label);
1169 gmx_fio_do_int(fio,iparams->disres.type);
1170 gmx_fio_do_real(fio,iparams->disres.low);
1171 gmx_fio_do_real(fio,iparams->disres.up1);
1172 gmx_fio_do_real(fio,iparams->disres.up2);
1173 gmx_fio_do_real(fio,iparams->disres.kfac);
1176 gmx_fio_do_int(fio,iparams->orires.ex);
1177 gmx_fio_do_int(fio,iparams->orires.label);
1178 gmx_fio_do_int(fio,iparams->orires.power);
1179 gmx_fio_do_real(fio,iparams->orires.c);
1180 gmx_fio_do_real(fio,iparams->orires.obs);
1181 gmx_fio_do_real(fio,iparams->orires.kfac);
1184 gmx_fio_do_int(fio,iparams->dihres.power);
1185 gmx_fio_do_int(fio,iparams->dihres.label);
1186 gmx_fio_do_real(fio,iparams->dihres.phi);
1187 gmx_fio_do_real(fio,iparams->dihres.dphi);
1188 gmx_fio_do_real(fio,iparams->dihres.kfac);
1191 gmx_fio_do_rvec(fio,iparams->posres.pos0A);
1192 gmx_fio_do_rvec(fio,iparams->posres.fcA);
1193 if (bRead && file_version < 27) {
1194 copy_rvec(iparams->posres.pos0A,iparams->posres.pos0B);
1195 copy_rvec(iparams->posres.fcA,iparams->posres.fcB);
1197 gmx_fio_do_rvec(fio,iparams->posres.pos0B);
1198 gmx_fio_do_rvec(fio,iparams->posres.fcB);
1202 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA,NR_RBDIHS);
1203 if(file_version>=25)
1204 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB,NR_RBDIHS);
1207 /* Fourier dihedrals are internally represented
1208 * as Ryckaert-Bellemans since those are faster to compute.
1210 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA, NR_RBDIHS);
1211 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB, NR_RBDIHS);
1215 gmx_fio_do_real(fio,iparams->constr.dA);
1216 gmx_fio_do_real(fio,iparams->constr.dB);
1219 gmx_fio_do_real(fio,iparams->settle.doh);
1220 gmx_fio_do_real(fio,iparams->settle.dhh);
1223 gmx_fio_do_real(fio,iparams->vsite.a);
1228 gmx_fio_do_real(fio,iparams->vsite.a);
1229 gmx_fio_do_real(fio,iparams->vsite.b);
1234 gmx_fio_do_real(fio,iparams->vsite.a);
1235 gmx_fio_do_real(fio,iparams->vsite.b);
1236 gmx_fio_do_real(fio,iparams->vsite.c);
1239 gmx_fio_do_int(fio,iparams->vsiten.n);
1240 gmx_fio_do_real(fio,iparams->vsiten.a);
1245 /* We got rid of some parameters in version 68 */
1246 if(bRead && file_version<68)
1248 gmx_fio_do_real(fio,rdum);
1249 gmx_fio_do_real(fio,rdum);
1250 gmx_fio_do_real(fio,rdum);
1251 gmx_fio_do_real(fio,rdum);
1253 gmx_fio_do_real(fio,iparams->gb.sar);
1254 gmx_fio_do_real(fio,iparams->gb.st);
1255 gmx_fio_do_real(fio,iparams->gb.pi);
1256 gmx_fio_do_real(fio,iparams->gb.gbr);
1257 gmx_fio_do_real(fio,iparams->gb.bmlt);
1260 gmx_fio_do_int(fio,iparams->cmap.cmapA);
1261 gmx_fio_do_int(fio,iparams->cmap.cmapB);
1264 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
1266 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
1269 gmx_fio_unset_comment(fio);
1272 static void do_ilist(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,int file_version,
1279 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1281 if (file_version < 44) {
1282 for(i=0; i<MAXNODES; i++)
1283 gmx_fio_do_int(fio,idum);
1285 gmx_fio_do_int(fio,ilist->nr);
1287 snew(ilist->iatoms,ilist->nr);
1288 bDum=gmx_fio_ndo_int(fio,ilist->iatoms,ilist->nr);
1290 gmx_fio_unset_comment(fio);
1293 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1294 gmx_bool bRead, int file_version)
1300 gmx_fio_do_int(fio,ffparams->atnr);
1301 if (file_version < 57) {
1302 gmx_fio_do_int(fio,idum);
1304 gmx_fio_do_int(fio,ffparams->ntypes);
1306 fprintf(debug,"ffparams->atnr = %d, ntypes = %d\n",
1307 ffparams->atnr,ffparams->ntypes);
1309 snew(ffparams->functype,ffparams->ntypes);
1310 snew(ffparams->iparams,ffparams->ntypes);
1312 /* Read/write all the function types */
1313 bDum=gmx_fio_ndo_int(fio,ffparams->functype,ffparams->ntypes);
1315 pr_ivec(debug,0,"functype",ffparams->functype,ffparams->ntypes,TRUE);
1317 if (file_version >= 66) {
1318 gmx_fio_do_double(fio,ffparams->reppow);
1320 ffparams->reppow = 12.0;
1323 if (file_version >= 57) {
1324 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1327 /* Check whether all these function types are supported by the code.
1328 * In practice the code is backwards compatible, which means that the
1329 * numbering may have to be altered from old numbering to new numbering
1331 for (i=0; (i<ffparams->ntypes); i++) {
1333 /* Loop over file versions */
1334 for (k=0; (k<NFTUPD); k++)
1335 /* Compare the read file_version to the update table */
1336 if ((file_version < ftupd[k].fvnr) &&
1337 (ffparams->functype[i] >= ftupd[k].ftype)) {
1338 ffparams->functype[i] += 1;
1340 fprintf(debug,"Incrementing function type %d to %d (due to %s)\n",
1341 i,ffparams->functype[i],
1342 interaction_function[ftupd[k].ftype].longname);
1347 do_iparams(fio, ffparams->functype[i],&ffparams->iparams[i],bRead,
1350 pr_iparams(debug,ffparams->functype[i],&ffparams->iparams[i]);
1354 static void do_ilists(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,
1357 int i,j,renum[F_NRE];
1358 gmx_bool bDum=TRUE,bClear;
1361 for(j=0; (j<F_NRE); j++) {
1364 for (k=0; k<NFTUPD; k++)
1365 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
1369 ilist[j].iatoms = NULL;
1371 do_ilist(fio, &ilist[j],bRead,file_version,j);
1374 if (bRead && gmx_debug_at)
1375 pr_ilist(debug,0,interaction_function[j].longname,
1376 functype,&ilist[j],TRUE);
1381 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams,gmx_moltype_t *molt,
1382 gmx_bool bRead, int file_version)
1384 do_ffparams(fio, ffparams,bRead,file_version);
1386 if (file_version >= 54) {
1387 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1390 do_ilists(fio, molt->ilist,bRead,file_version);
1393 static void do_block(t_fileio *fio, t_block *block,gmx_bool bRead,int file_version)
1395 int i,idum,dum_nra,*dum_a;
1398 if (file_version < 44)
1399 for(i=0; i<MAXNODES; i++)
1400 gmx_fio_do_int(fio,idum);
1401 gmx_fio_do_int(fio,block->nr);
1402 if (file_version < 51)
1403 gmx_fio_do_int(fio,dum_nra);
1405 block->nalloc_index = block->nr+1;
1406 snew(block->index,block->nalloc_index);
1408 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1410 if (file_version < 51 && dum_nra > 0) {
1411 snew(dum_a,dum_nra);
1412 bDum=gmx_fio_ndo_int(fio,dum_a,dum_nra);
1417 static void do_blocka(t_fileio *fio, t_blocka *block,gmx_bool bRead,
1423 if (file_version < 44)
1424 for(i=0; i<MAXNODES; i++)
1425 gmx_fio_do_int(fio,idum);
1426 gmx_fio_do_int(fio,block->nr);
1427 gmx_fio_do_int(fio,block->nra);
1429 block->nalloc_index = block->nr+1;
1430 snew(block->index,block->nalloc_index);
1431 block->nalloc_a = block->nra;
1432 snew(block->a,block->nalloc_a);
1434 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1435 bDum=gmx_fio_ndo_int(fio,block->a,block->nra);
1438 static void do_atom(t_fileio *fio, t_atom *atom,int ngrp,gmx_bool bRead,
1439 int file_version, gmx_groups_t *groups,int atnr)
1443 gmx_fio_do_real(fio,atom->m);
1444 gmx_fio_do_real(fio,atom->q);
1445 gmx_fio_do_real(fio,atom->mB);
1446 gmx_fio_do_real(fio,atom->qB);
1447 gmx_fio_do_ushort(fio, atom->type);
1448 gmx_fio_do_ushort(fio, atom->typeB);
1449 gmx_fio_do_int(fio,atom->ptype);
1450 gmx_fio_do_int(fio,atom->resind);
1451 if (file_version >= 52)
1452 gmx_fio_do_int(fio,atom->atomnumber);
1454 atom->atomnumber = NOTSET;
1455 if (file_version < 23)
1457 else if (file_version < 39)
1462 if (file_version < 57) {
1463 unsigned char uchar[egcNR];
1464 gmx_fio_ndo_uchar(fio,uchar,myngrp);
1465 for(i=myngrp; (i<ngrp); i++) {
1468 /* Copy the old data format to the groups struct */
1469 for(i=0; i<ngrp; i++) {
1470 groups->grpnr[i][atnr] = uchar[i];
1475 static void do_grps(t_fileio *fio, int ngrp,t_grps grps[],gmx_bool bRead,
1481 if (file_version < 23)
1483 else if (file_version < 39)
1488 for(j=0; (j<ngrp); j++) {
1490 gmx_fio_do_int(fio,grps[j].nr);
1492 snew(grps[j].nm_ind,grps[j].nr);
1493 bDum=gmx_fio_ndo_int(fio,grps[j].nm_ind,grps[j].nr);
1497 snew(grps[j].nm_ind,grps[j].nr);
1502 static void do_symstr(t_fileio *fio, char ***nm,gmx_bool bRead,t_symtab *symtab)
1507 gmx_fio_do_int(fio,ls);
1508 *nm = get_symtab_handle(symtab,ls);
1511 ls = lookup_symtab(symtab,*nm);
1512 gmx_fio_do_int(fio,ls);
1516 static void do_strstr(t_fileio *fio, int nstr,char ***nm,gmx_bool bRead,
1521 for (j=0; (j<nstr); j++)
1522 do_symstr(fio, &(nm[j]),bRead,symtab);
1525 static void do_resinfo(t_fileio *fio, int n,t_resinfo *ri,gmx_bool bRead,
1526 t_symtab *symtab, int file_version)
1530 for (j=0; (j<n); j++) {
1531 do_symstr(fio, &(ri[j].name),bRead,symtab);
1532 if (file_version >= 63) {
1533 gmx_fio_do_int(fio,ri[j].nr);
1534 gmx_fio_do_uchar(fio, ri[j].ic);
1542 static void do_atoms(t_fileio *fio, t_atoms *atoms,gmx_bool bRead,t_symtab *symtab,
1544 gmx_groups_t *groups)
1548 gmx_fio_do_int(fio,atoms->nr);
1549 gmx_fio_do_int(fio,atoms->nres);
1550 if (file_version < 57) {
1551 gmx_fio_do_int(fio,groups->ngrpname);
1552 for(i=0; i<egcNR; i++) {
1553 groups->ngrpnr[i] = atoms->nr;
1554 snew(groups->grpnr[i],groups->ngrpnr[i]);
1558 snew(atoms->atom,atoms->nr);
1559 snew(atoms->atomname,atoms->nr);
1560 snew(atoms->atomtype,atoms->nr);
1561 snew(atoms->atomtypeB,atoms->nr);
1562 snew(atoms->resinfo,atoms->nres);
1563 if (file_version < 57) {
1564 snew(groups->grpname,groups->ngrpname);
1566 atoms->pdbinfo = NULL;
1568 for(i=0; (i<atoms->nr); i++) {
1569 do_atom(fio, &atoms->atom[i],egcNR,bRead, file_version,groups,i);
1571 do_strstr(fio, atoms->nr,atoms->atomname,bRead,symtab);
1572 if (bRead && (file_version <= 20)) {
1573 for(i=0; i<atoms->nr; i++) {
1574 atoms->atomtype[i] = put_symtab(symtab,"?");
1575 atoms->atomtypeB[i] = put_symtab(symtab,"?");
1578 do_strstr(fio, atoms->nr,atoms->atomtype,bRead,symtab);
1579 do_strstr(fio, atoms->nr,atoms->atomtypeB,bRead,symtab);
1581 do_resinfo(fio, atoms->nres,atoms->resinfo,bRead,symtab,file_version);
1583 if (file_version < 57) {
1584 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1586 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1590 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
1591 gmx_bool bRead,t_symtab *symtab,
1597 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1598 gmx_fio_do_int(fio,groups->ngrpname);
1600 snew(groups->grpname,groups->ngrpname);
1602 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1603 for(g=0; g<egcNR; g++) {
1604 gmx_fio_do_int(fio,groups->ngrpnr[g]);
1605 if (groups->ngrpnr[g] == 0) {
1607 groups->grpnr[g] = NULL;
1611 snew(groups->grpnr[g],groups->ngrpnr[g]);
1613 bDum=gmx_fio_ndo_uchar(fio, groups->grpnr[g],groups->ngrpnr[g]);
1618 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes,gmx_bool bRead,
1619 t_symtab *symtab,int file_version)
1622 gmx_bool bDum = TRUE;
1624 if (file_version > 25) {
1625 gmx_fio_do_int(fio,atomtypes->nr);
1628 snew(atomtypes->radius,j);
1629 snew(atomtypes->vol,j);
1630 snew(atomtypes->surftens,j);
1631 snew(atomtypes->atomnumber,j);
1632 snew(atomtypes->gb_radius,j);
1633 snew(atomtypes->S_hct,j);
1635 bDum=gmx_fio_ndo_real(fio,atomtypes->radius,j);
1636 bDum=gmx_fio_ndo_real(fio,atomtypes->vol,j);
1637 bDum=gmx_fio_ndo_real(fio,atomtypes->surftens,j);
1638 if(file_version >= 40)
1640 bDum=gmx_fio_ndo_int(fio,atomtypes->atomnumber,j);
1642 if(file_version >= 60)
1644 bDum=gmx_fio_ndo_real(fio,atomtypes->gb_radius,j);
1645 bDum=gmx_fio_ndo_real(fio,atomtypes->S_hct,j);
1648 /* File versions prior to 26 cannot do GBSA,
1649 * so they dont use this structure
1652 atomtypes->radius = NULL;
1653 atomtypes->vol = NULL;
1654 atomtypes->surftens = NULL;
1655 atomtypes->atomnumber = NULL;
1656 atomtypes->gb_radius = NULL;
1657 atomtypes->S_hct = NULL;
1661 static void do_symtab(t_fileio *fio, t_symtab *symtab,gmx_bool bRead)
1667 gmx_fio_do_int(fio,symtab->nr);
1670 snew(symtab->symbuf,1);
1671 symbuf = symtab->symbuf;
1672 symbuf->bufsize = nr;
1673 snew(symbuf->buf,nr);
1674 for (i=0; (i<nr); i++) {
1675 gmx_fio_do_string(fio,buf);
1676 symbuf->buf[i]=strdup(buf);
1680 symbuf = symtab->symbuf;
1681 while (symbuf!=NULL) {
1682 for (i=0; (i<symbuf->bufsize) && (i<nr); i++)
1683 gmx_fio_do_string(fio,symbuf->buf[i]);
1685 symbuf=symbuf->next;
1688 gmx_fatal(FARGS,"nr of symtab strings left: %d",nr);
1692 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
1694 int i,j,ngrid,gs,nelem;
1696 gmx_fio_do_int(fio,cmap_grid->ngrid);
1697 gmx_fio_do_int(fio,cmap_grid->grid_spacing);
1699 ngrid = cmap_grid->ngrid;
1700 gs = cmap_grid->grid_spacing;
1705 snew(cmap_grid->cmapdata,ngrid);
1707 for(i=0;i<cmap_grid->ngrid;i++)
1709 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
1713 for(i=0;i<cmap_grid->ngrid;i++)
1715 for(j=0;j<nelem;j++)
1717 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4]);
1718 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+1]);
1719 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+2]);
1720 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+3]);
1726 void tpx_make_chain_identifiers(t_atoms *atoms,t_block *mols)
1732 /* We always assign a new chain number, but save the chain id characters
1733 * for larger molecules.
1735 #define CHAIN_MIN_ATOMS 15
1739 for(m=0; m<mols->nr; m++)
1742 a1=mols->index[m+1];
1743 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
1752 for(a=a0; a<a1; a++)
1754 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
1755 atoms->resinfo[atoms->atom[a].resind].chainid = c;
1760 /* Blank out the chain id if there was only one chain */
1763 for(r=0; r<atoms->nres; r++)
1765 atoms->resinfo[r].chainid = ' ';
1770 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt,gmx_bool bRead,
1771 t_symtab *symtab, int file_version,
1772 gmx_groups_t *groups)
1776 if (file_version >= 57) {
1777 do_symstr(fio, &(molt->name),bRead,symtab);
1780 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
1782 if (bRead && gmx_debug_at) {
1783 pr_atoms(debug,0,"atoms",&molt->atoms,TRUE);
1786 if (file_version >= 57) {
1787 do_ilists(fio, molt->ilist,bRead,file_version);
1789 do_block(fio, &molt->cgs,bRead,file_version);
1790 if (bRead && gmx_debug_at) {
1791 pr_block(debug,0,"cgs",&molt->cgs,TRUE);
1795 /* This used to be in the atoms struct */
1796 do_blocka(fio, &molt->excls, bRead, file_version);
1799 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb,gmx_bool bRead,
1804 gmx_fio_do_int(fio,molb->type);
1805 gmx_fio_do_int(fio,molb->nmol);
1806 gmx_fio_do_int(fio,molb->natoms_mol);
1807 /* Position restraint coordinates */
1808 gmx_fio_do_int(fio,molb->nposres_xA);
1809 if (molb->nposres_xA > 0) {
1811 snew(molb->posres_xA,molb->nposres_xA);
1813 gmx_fio_ndo_rvec(fio,molb->posres_xA,molb->nposres_xA);
1815 gmx_fio_do_int(fio,molb->nposres_xB);
1816 if (molb->nposres_xB > 0) {
1818 snew(molb->posres_xB,molb->nposres_xB);
1820 gmx_fio_ndo_rvec(fio,molb->posres_xB,molb->nposres_xB);
1825 static t_block mtop_mols(gmx_mtop_t *mtop)
1831 for(mb=0; mb<mtop->nmolblock; mb++) {
1832 mols.nr += mtop->molblock[mb].nmol;
1834 mols.nalloc_index = mols.nr + 1;
1835 snew(mols.index,mols.nalloc_index);
1840 for(mb=0; mb<mtop->nmolblock; mb++) {
1841 for(mol=0; mol<mtop->molblock[mb].nmol; mol++) {
1842 a += mtop->molblock[mb].natoms_mol;
1851 static void add_posres_molblock(gmx_mtop_t *mtop)
1856 gmx_molblock_t *molb;
1859 il = &mtop->moltype[0].ilist[F_POSRES];
1865 for(i=0; i<il->nr; i+=2) {
1866 ip = &mtop->ffparams.iparams[il->iatoms[i]];
1867 am = max(am,il->iatoms[i+1]);
1868 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
1869 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
1870 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ]) {
1874 /* Make the posres coordinate block end at a molecule end */
1876 while(am >= mtop->mols.index[mol+1]) {
1879 molb = &mtop->molblock[0];
1880 molb->nposres_xA = mtop->mols.index[mol+1];
1881 snew(molb->posres_xA,molb->nposres_xA);
1883 molb->nposres_xB = molb->nposres_xA;
1884 snew(molb->posres_xB,molb->nposres_xB);
1886 molb->nposres_xB = 0;
1888 for(i=0; i<il->nr; i+=2) {
1889 ip = &mtop->ffparams.iparams[il->iatoms[i]];
1890 a = il->iatoms[i+1];
1891 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
1892 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
1893 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
1895 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
1896 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
1897 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
1902 static void set_disres_npair(gmx_mtop_t *mtop)
1909 ip = mtop->ffparams.iparams;
1911 for(mt=0; mt<mtop->nmoltype; mt++) {
1912 il = &mtop->moltype[mt].ilist[F_DISRES];
1916 for(i=0; i<il->nr; i+=3) {
1918 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label) {
1919 ip[a[i]].disres.npair = npair;
1927 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop,gmx_bool bRead,
1935 do_symtab(fio, &(mtop->symtab),bRead);
1937 pr_symtab(debug,0,"symtab",&mtop->symtab);
1939 do_symstr(fio, &(mtop->name),bRead,&(mtop->symtab));
1941 if (file_version >= 57) {
1942 do_ffparams(fio, &mtop->ffparams,bRead,file_version);
1944 gmx_fio_do_int(fio,mtop->nmoltype);
1949 snew(mtop->moltype,mtop->nmoltype);
1950 if (file_version < 57) {
1951 mtop->moltype[0].name = mtop->name;
1954 for(mt=0; mt<mtop->nmoltype; mt++) {
1955 do_moltype(fio, &mtop->moltype[mt],bRead,&mtop->symtab,file_version,
1959 if (file_version >= 57) {
1960 gmx_fio_do_int(fio,mtop->nmolblock);
1962 mtop->nmolblock = 1;
1965 snew(mtop->molblock,mtop->nmolblock);
1967 if (file_version >= 57) {
1968 for(mb=0; mb<mtop->nmolblock; mb++) {
1969 do_molblock(fio, &mtop->molblock[mb],bRead,file_version);
1971 gmx_fio_do_int(fio,mtop->natoms);
1973 mtop->molblock[0].type = 0;
1974 mtop->molblock[0].nmol = 1;
1975 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
1976 mtop->molblock[0].nposres_xA = 0;
1977 mtop->molblock[0].nposres_xB = 0;
1980 do_atomtypes (fio, &(mtop->atomtypes),bRead,&(mtop->symtab), file_version);
1982 pr_atomtypes(debug,0,"atomtypes",&mtop->atomtypes,TRUE);
1984 if (file_version < 57) {
1985 /* Debug statements are inside do_idef */
1986 do_idef (fio, &mtop->ffparams,&mtop->moltype[0],bRead,file_version);
1987 mtop->natoms = mtop->moltype[0].atoms.nr;
1990 if(file_version >= 65)
1992 do_cmap(fio, &mtop->ffparams.cmap_grid,bRead);
1996 mtop->ffparams.cmap_grid.ngrid = 0;
1997 mtop->ffparams.cmap_grid.grid_spacing = 0.1;
1998 mtop->ffparams.cmap_grid.cmapdata = NULL;
2001 if (file_version >= 57) {
2002 do_groups(fio, &mtop->groups,bRead,&(mtop->symtab),file_version);
2005 if (file_version < 57) {
2006 do_block(fio, &mtop->moltype[0].cgs,bRead,file_version);
2007 if (bRead && gmx_debug_at) {
2008 pr_block(debug,0,"cgs",&mtop->moltype[0].cgs,TRUE);
2010 do_block(fio, &mtop->mols,bRead,file_version);
2011 /* Add the posres coordinates to the molblock */
2012 add_posres_molblock(mtop);
2015 if (file_version >= 57) {
2016 mtop->mols = mtop_mols(mtop);
2019 pr_block(debug,0,"mols",&mtop->mols,TRUE);
2023 if (file_version < 51) {
2024 /* Here used to be the shake blocks */
2025 do_blocka(fio, &dumb,bRead,file_version);
2033 close_symtab(&(mtop->symtab));
2037 /* If TopOnlyOK is TRUE then we can read even future versions
2038 * of tpx files, provided the file_generation hasn't changed.
2039 * If it is FALSE, we need the inputrecord too, and bail out
2040 * if the file is newer than the program.
2042 * The version and generation if the topology (see top of this file)
2043 * are returned in the two last arguments.
2045 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2047 static void do_tpxheader(t_fileio *fio,gmx_bool bRead,t_tpxheader *tpx,
2048 gmx_bool TopOnlyOK, int *file_version,
2049 int *file_generation)
2058 gmx_fio_checktype(fio);
2059 gmx_fio_setdebug(fio,bDebugMode());
2061 /* NEW! XDR tpb file */
2062 precision = sizeof(real);
2064 gmx_fio_do_string(fio,buf);
2065 if (strncmp(buf,"VERSION",7))
2066 gmx_fatal(FARGS,"Can not read file %s,\n"
2067 " this file is from a Gromacs version which is older than 2.0\n"
2068 " Make a new one with grompp or use a gro or pdb file, if possible",
2069 gmx_fio_getname(fio));
2070 gmx_fio_do_int(fio,precision);
2071 bDouble = (precision == sizeof(double));
2072 if ((precision != sizeof(float)) && !bDouble)
2073 gmx_fatal(FARGS,"Unknown precision in file %s: real is %d bytes "
2074 "instead of %d or %d",
2075 gmx_fio_getname(fio),precision,sizeof(float),sizeof(double));
2076 gmx_fio_setprecision(fio,bDouble);
2077 fprintf(stderr,"Reading file %s, %s (%s precision)\n",
2078 gmx_fio_getname(fio),buf,bDouble ? "double" : "single");
2081 gmx_fio_write_string(fio,GromacsVersion());
2082 bDouble = (precision == sizeof(double));
2083 gmx_fio_setprecision(fio,bDouble);
2084 gmx_fio_do_int(fio,precision);
2086 fgen = tpx_generation;
2089 /* Check versions! */
2090 gmx_fio_do_int(fio,fver);
2093 gmx_fio_do_int(fio,fgen);
2097 if(file_version!=NULL)
2098 *file_version = fver;
2099 if(file_generation!=NULL)
2100 *file_generation = fgen;
2103 if ((fver <= tpx_incompatible_version) ||
2104 ((fver > tpx_version) && !TopOnlyOK) ||
2105 (fgen > tpx_generation))
2106 gmx_fatal(FARGS,"reading tpx file (%s) version %d with version %d program",
2107 gmx_fio_getname(fio),fver,tpx_version);
2109 do_section(fio,eitemHEADER,bRead);
2110 gmx_fio_do_int(fio,tpx->natoms);
2112 gmx_fio_do_int(fio,tpx->ngtc);
2116 gmx_fio_do_int(fio,idum);
2117 gmx_fio_do_real(fio,rdum);
2119 gmx_fio_do_real(fio,tpx->lambda);
2120 gmx_fio_do_int(fio,tpx->bIr);
2121 gmx_fio_do_int(fio,tpx->bTop);
2122 gmx_fio_do_int(fio,tpx->bX);
2123 gmx_fio_do_int(fio,tpx->bV);
2124 gmx_fio_do_int(fio,tpx->bF);
2125 gmx_fio_do_int(fio,tpx->bBox);
2127 if((fgen > tpx_generation)) {
2128 /* This can only happen if TopOnlyOK=TRUE */
2133 static int do_tpx(t_fileio *fio, gmx_bool bRead,
2134 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop,
2135 gmx_bool bXVallocated)
2140 gmx_bool TopOnlyOK,bDum=TRUE;
2141 int file_version,file_generation;
2145 gmx_bool bPeriodicMols;
2148 tpx.natoms = state->natoms;
2149 tpx.ngtc = state->ngtc;
2150 tpx.lambda = state->lambda;
2151 tpx.bIr = (ir != NULL);
2152 tpx.bTop = (mtop != NULL);
2153 tpx.bX = (state->x != NULL);
2154 tpx.bV = (state->v != NULL);
2155 tpx.bF = (f != NULL);
2159 TopOnlyOK = (ir==NULL);
2161 do_tpxheader(fio,bRead,&tpx,TopOnlyOK,&file_version,&file_generation);
2165 state->lambda = tpx.lambda;
2166 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
2170 init_state(state,0,tpx.ngtc,0,0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
2171 state->natoms = tpx.natoms;
2172 state->nalloc = tpx.natoms;
2176 init_state(state,tpx.natoms,tpx.ngtc,0,0); /* nose-hoover chains */
2180 #define do_test(fio,b,p) if (bRead && (p!=NULL) && !b) gmx_fatal(FARGS,"No %s in %s",#p,gmx_fio_getname(fio))
2182 do_test(fio,tpx.bBox,state->box);
2183 do_section(fio,eitemBOX,bRead);
2185 gmx_fio_ndo_rvec(fio,state->box,DIM);
2186 if (file_version >= 51) {
2187 gmx_fio_ndo_rvec(fio,state->box_rel,DIM);
2189 /* We initialize box_rel after reading the inputrec */
2190 clear_mat(state->box_rel);
2192 if (file_version >= 28) {
2193 gmx_fio_ndo_rvec(fio,state->boxv,DIM);
2194 if (file_version < 56) {
2196 gmx_fio_ndo_rvec(fio,mdum,DIM);
2201 if (state->ngtc > 0 && file_version >= 28) {
2203 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
2204 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
2205 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
2206 snew(dumv,state->ngtc);
2207 if (file_version < 69) {
2208 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2210 /* These used to be the Berendsen tcoupl_lambda's */
2211 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2215 /* Prior to tpx version 26, the inputrec was here.
2216 * I moved it to enable partial forward-compatibility
2217 * for analysis/viewer programs.
2219 if(file_version<26) {
2220 do_test(fio,tpx.bIr,ir);
2221 do_section(fio,eitemIR,bRead);
2224 do_inputrec(fio, ir,bRead,file_version,
2225 mtop ? &mtop->ffparams.fudgeQQ : NULL);
2227 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2230 do_inputrec(fio, &dum_ir,bRead,file_version,
2231 mtop ? &mtop->ffparams.fudgeQQ :NULL);
2233 pr_inputrec(debug,0,"inputrec",&dum_ir,FALSE);
2234 done_inputrec(&dum_ir);
2240 do_test(fio,tpx.bTop,mtop);
2241 do_section(fio,eitemTOP,bRead);
2244 do_mtop(fio,mtop,bRead, file_version);
2246 do_mtop(fio,&dum_top,bRead,file_version);
2247 done_mtop(&dum_top,TRUE);
2250 do_test(fio,tpx.bX,state->x);
2251 do_section(fio,eitemX,bRead);
2254 state->flags |= (1<<estX);
2256 gmx_fio_ndo_rvec(fio,state->x,state->natoms);
2259 do_test(fio,tpx.bV,state->v);
2260 do_section(fio,eitemV,bRead);
2263 state->flags |= (1<<estV);
2265 gmx_fio_ndo_rvec(fio,state->v,state->natoms);
2268 do_test(fio,tpx.bF,f);
2269 do_section(fio,eitemF,bRead);
2270 if (tpx.bF) gmx_fio_ndo_rvec(fio,f,state->natoms);
2272 /* Starting with tpx version 26, we have the inputrec
2273 * at the end of the file, so we can ignore it
2274 * if the file is never than the software (but still the
2275 * same generation - see comments at the top of this file.
2280 bPeriodicMols = FALSE;
2281 if (file_version >= 26) {
2282 do_test(fio,tpx.bIr,ir);
2283 do_section(fio,eitemIR,bRead);
2285 if (file_version >= 53) {
2286 /* Removed the pbc info from do_inputrec, since we always want it */
2289 bPeriodicMols = ir->bPeriodicMols;
2291 gmx_fio_do_int(fio,ePBC);
2292 gmx_fio_do_gmx_bool(fio,bPeriodicMols);
2294 if (file_generation <= tpx_generation && ir) {
2295 do_inputrec(fio, ir,bRead,file_version,mtop ? &mtop->ffparams.fudgeQQ : NULL);
2297 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2298 if (file_version < 51)
2299 set_box_rel(ir,state);
2300 if (file_version < 53) {
2302 bPeriodicMols = ir->bPeriodicMols;
2305 if (bRead && ir && file_version >= 53) {
2306 /* We need to do this after do_inputrec, since that initializes ir */
2308 ir->bPeriodicMols = bPeriodicMols;
2317 if (state->ngtc == 0)
2319 /* Reading old version without tcoupl state data: set it */
2320 init_gtc_state(state,ir->opts.ngtc,0,ir->opts.nhchainlength);
2322 if (tpx.bTop && mtop)
2324 if (file_version < 57)
2326 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
2328 ir->eDisre = edrSimple;
2332 ir->eDisre = edrNone;
2335 set_disres_npair(mtop);
2339 if (tpx.bTop && mtop)
2341 gmx_mtop_finalize(mtop);
2344 if (file_version >= 57)
2348 env = getenv("GMX_NOCHARGEGROUPS");
2351 sscanf(env,"%d",&ienv);
2352 fprintf(stderr,"\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
2357 "Will make single atomic charge groups in non-solvent%s\n",
2358 ienv > 1 ? " and solvent" : "");
2359 gmx_mtop_make_atomic_charge_groups(mtop,ienv==1);
2361 fprintf(stderr,"\n");
2369 /************************************************************
2371 * The following routines are the exported ones
2373 ************************************************************/
2375 t_fileio *open_tpx(const char *fn,const char *mode)
2377 return gmx_fio_open(fn,mode);
2380 void close_tpx(t_fileio *fio)
2385 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
2386 int *file_version, int *file_generation)
2390 fio = open_tpx(fn,"r");
2391 do_tpxheader(fio,TRUE,tpx,TopOnlyOK,file_version,file_generation);
2395 void write_tpx_state(const char *fn,
2396 t_inputrec *ir,t_state *state,gmx_mtop_t *mtop)
2400 fio = open_tpx(fn,"w");
2401 do_tpx(fio,FALSE,ir,state,NULL,mtop,FALSE);
2405 void read_tpx_state(const char *fn,
2406 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop)
2410 fio = open_tpx(fn,"r");
2411 do_tpx(fio,TRUE,ir,state,f,mtop,FALSE);
2415 int read_tpx(const char *fn,
2416 t_inputrec *ir, matrix box,int *natoms,
2417 rvec *x,rvec *v,rvec *f,gmx_mtop_t *mtop)
2425 fio = open_tpx(fn,"r");
2426 ePBC = do_tpx(fio,TRUE,ir,&state,f,mtop,TRUE);
2428 *natoms = state.natoms;
2430 copy_mat(state.box,box);
2438 int read_tpx_top(const char *fn,
2439 t_inputrec *ir, matrix box,int *natoms,
2440 rvec *x,rvec *v,rvec *f,t_topology *top)
2446 ePBC = read_tpx(fn,ir,box,natoms,x,v,f,&mtop);
2448 *top = gmx_mtop_t_to_t_topology(&mtop);
2453 gmx_bool fn2bTPX(const char *file)
2455 switch (fn2ftp(file)) {
2465 gmx_bool read_tps_conf(const char *infile,char *title,t_topology *top,int *ePBC,
2466 rvec **x,rvec **v,matrix box,gmx_bool bMass)
2469 int natoms,i,version,generation;
2470 gmx_bool bTop,bXNULL;
2472 t_topology *topconv;
2475 bTop = fn2bTPX(infile);
2478 read_tpxheader(infile,&header,TRUE,&version,&generation);
2480 snew(*x,header.natoms);
2482 snew(*v,header.natoms);
2484 *ePBC = read_tpx(infile,NULL,box,&natoms,
2485 (x==NULL) ? NULL : *x,(v==NULL) ? NULL : *v,NULL,mtop);
2486 *top = gmx_mtop_t_to_t_topology(mtop);
2488 strcpy(title,*top->name);
2489 tpx_make_chain_identifiers(&top->atoms,&top->mols);
2492 get_stx_coordnum(infile,&natoms);
2493 init_t_atoms(&top->atoms,natoms,FALSE);
2494 bXNULL = (x == NULL);
2498 read_stx_conf(infile,title,&top->atoms,*x,(v==NULL) ? NULL : *v,ePBC,box);
2504 aps = gmx_atomprop_init();
2505 for(i=0; (i<natoms); i++)
2506 if (!gmx_atomprop_query(aps,epropMass,
2507 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2508 *top->atoms.atomname[i],
2509 &(top->atoms.atom[i].m))) {
2511 fprintf(debug,"Can not find mass for atom %s %d %s, setting to 1\n",
2512 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2513 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
2514 *top->atoms.atomname[i]);
2516 gmx_atomprop_destroy(aps);
2518 top->idef.ntypes=-1;