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 = 75;
69 /* This number should only be increased when you edit the TOPOLOGY section
70 * of the tpx format. This way we can maintain forward compatibility too
71 * for all analysis tools and/or external programs that only need to
72 * know the atom/residue names, charges, and bond connectivity.
74 * It first appeared in tpx version 26, when I also moved the inputrecord
75 * to the end of the tpx file, so we can just skip it if we only
78 static const int tpx_generation = 23;
80 /* This number should be the most recent backwards incompatible version
81 * I.e., if this number is 9, we cannot read tpx version 9 with this code.
83 static const int tpx_incompatible_version = 9;
87 /* Struct used to maintain tpx compatibility when function types are added */
89 int fvnr; /* file version number in which the function type first appeared */
90 int ftype; /* function type */
94 *The entries should be ordered in:
95 * 1. ascending file version number
96 * 2. ascending function type number
98 /*static const t_ftupd ftupd[] = {
103 { 22, F_DISRESVIOL },
109 { 26, F_DIHRESVIOL },
110 { 30, F_CROSS_BOND_BONDS },
111 { 30, F_CROSS_BOND_ANGLES },
112 { 30, F_UREY_BRADLEY },
113 { 30, F_POLARIZATION }
117 *The entries should be ordered in:
118 * 1. ascending function type number
119 * 2. ascending file version number
121 static const t_ftupd ftupd[] = {
122 { 20, F_CUBICBONDS },
127 { 43, F_TABBONDSNC },
128 { 70, F_RESTRBONDS },
129 { 30, F_CROSS_BOND_BONDS },
130 { 30, F_CROSS_BOND_ANGLES },
131 { 30, F_UREY_BRADLEY },
132 { 34, F_QUARTIC_ANGLES },
142 { 72, F_NPSOLVATION },
144 { 41, F_LJC_PAIRS_NB },
147 { 32, F_COUL_RECIP },
149 { 30, F_POLARIZATION },
151 { 22, F_DISRESVIOL },
155 { 26, F_DIHRESVIOL },
160 { 46, F_ECONSERVED },
165 #define NFTUPD asize(ftupd)
167 /* Needed for backward compatibility */
170 static void _do_section(t_fileio *fio,int key,gmx_bool bRead,const char *src,
176 if (gmx_fio_getftp(fio) == efTPA) {
178 gmx_fio_write_string(fio,itemstr[key]);
179 bDbg = gmx_fio_getdebug(fio);
180 gmx_fio_setdebug(fio,FALSE);
181 gmx_fio_write_string(fio,comment_str[key]);
182 gmx_fio_setdebug(fio,bDbg);
185 if (gmx_fio_getdebug(fio))
186 fprintf(stderr,"Looking for section %s (%s, %d)",
187 itemstr[key],src,line);
190 gmx_fio_do_string(fio,buf);
191 } while ((gmx_strcasecmp(buf,itemstr[key]) != 0));
193 if (gmx_strcasecmp(buf,itemstr[key]) != 0)
194 gmx_fatal(FARGS,"\nCould not find section heading %s",itemstr[key]);
195 else if (gmx_fio_getdebug(fio))
196 fprintf(stderr," and found it\n");
201 #define do_section(fio,key,bRead) _do_section(fio,key,bRead,__FILE__,__LINE__)
203 /**************************************************************
205 * Now the higer level routines that do io of the structures and arrays
207 **************************************************************/
208 static void do_pullgrp(t_fileio *fio, t_pullgrp *pgrp, gmx_bool bRead,
214 gmx_fio_do_int(fio,pgrp->nat);
216 snew(pgrp->ind,pgrp->nat);
217 bDum=gmx_fio_ndo_int(fio,pgrp->ind,pgrp->nat);
218 gmx_fio_do_int(fio,pgrp->nweight);
220 snew(pgrp->weight,pgrp->nweight);
221 bDum=gmx_fio_ndo_real(fio,pgrp->weight,pgrp->nweight);
222 gmx_fio_do_int(fio,pgrp->pbcatom);
223 gmx_fio_do_rvec(fio,pgrp->vec);
224 gmx_fio_do_rvec(fio,pgrp->init);
225 gmx_fio_do_real(fio,pgrp->rate);
226 gmx_fio_do_real(fio,pgrp->k);
227 if (file_version >= 56) {
228 gmx_fio_do_real(fio,pgrp->kB);
234 static void do_pull(t_fileio *fio, t_pull *pull,gmx_bool bRead, int file_version)
238 gmx_fio_do_int(fio,pull->ngrp);
239 gmx_fio_do_int(fio,pull->eGeom);
240 gmx_fio_do_ivec(fio,pull->dim);
241 gmx_fio_do_real(fio,pull->cyl_r1);
242 gmx_fio_do_real(fio,pull->cyl_r0);
243 gmx_fio_do_real(fio,pull->constr_tol);
244 gmx_fio_do_int(fio,pull->nstxout);
245 gmx_fio_do_int(fio,pull->nstfout);
247 snew(pull->grp,pull->ngrp+1);
248 for(g=0; g<pull->ngrp+1; g++)
249 do_pullgrp(fio,&pull->grp[g],bRead,file_version);
253 static void do_rotgrp(t_fileio *fio, t_rotgrp *rotg,gmx_bool bRead, int file_version)
258 gmx_fio_do_int(fio,rotg->eType);
259 gmx_fio_do_int(fio,rotg->bMassW);
260 gmx_fio_do_int(fio,rotg->nat);
262 snew(rotg->ind,rotg->nat);
263 gmx_fio_ndo_int(fio,rotg->ind,rotg->nat);
265 snew(rotg->x_ref,rotg->nat);
266 gmx_fio_ndo_rvec(fio,rotg->x_ref,rotg->nat);
267 gmx_fio_do_rvec(fio,rotg->vec);
268 gmx_fio_do_rvec(fio,rotg->pivot);
269 gmx_fio_do_real(fio,rotg->rate);
270 gmx_fio_do_real(fio,rotg->k);
271 gmx_fio_do_real(fio,rotg->slab_dist);
272 gmx_fio_do_real(fio,rotg->min_gaussian);
273 gmx_fio_do_real(fio,rotg->eps);
274 gmx_fio_do_int(fio,rotg->eFittype);
275 gmx_fio_do_int(fio,rotg->PotAngle_nstep);
276 gmx_fio_do_real(fio,rotg->PotAngle_step);
279 static void do_rot(t_fileio *fio, t_rot *rot,gmx_bool bRead, int file_version)
283 gmx_fio_do_int(fio,rot->ngrp);
284 gmx_fio_do_int(fio,rot->nstrout);
285 gmx_fio_do_int(fio,rot->nstsout);
287 snew(rot->grp,rot->ngrp);
288 for(g=0; g<rot->ngrp; g++)
289 do_rotgrp(fio, &rot->grp[g],bRead,file_version);
293 static void do_inputrec(t_fileio *fio, t_inputrec *ir,gmx_bool bRead,
294 int file_version, real *fudgeQQ)
296 int i,j,k,*tmp,idum=0;
301 real zerotemptime,finish_t,init_temp,finish_temp;
303 if (file_version != tpx_version)
305 /* Give a warning about features that are not accessible */
306 fprintf(stderr,"Note: tpx file_version %d, software version %d\n",
307 file_version,tpx_version);
315 if (file_version == 0)
320 /* Basic inputrec stuff */
321 gmx_fio_do_int(fio,ir->eI);
322 if (file_version >= 62) {
323 gmx_fio_do_gmx_large_int(fio, ir->nsteps);
325 gmx_fio_do_int(fio,idum);
328 if(file_version > 25) {
329 if (file_version >= 62) {
330 gmx_fio_do_gmx_large_int(fio, ir->init_step);
332 gmx_fio_do_int(fio,idum);
333 ir->init_step = idum;
339 if(file_version >= 58)
340 gmx_fio_do_int(fio,ir->simulation_part);
342 ir->simulation_part=1;
344 if (file_version >= 67) {
345 gmx_fio_do_int(fio,ir->nstcalcenergy);
347 ir->nstcalcenergy = 1;
349 if (file_version < 53) {
350 /* The pbc info has been moved out of do_inputrec,
351 * since we always want it, also without reading the inputrec.
353 gmx_fio_do_int(fio,ir->ePBC);
354 if ((file_version <= 15) && (ir->ePBC == 2))
356 if (file_version >= 45) {
357 gmx_fio_do_int(fio,ir->bPeriodicMols);
361 ir->bPeriodicMols = TRUE;
363 ir->bPeriodicMols = FALSE;
367 gmx_fio_do_int(fio,ir->ns_type);
368 gmx_fio_do_int(fio,ir->nstlist);
369 gmx_fio_do_int(fio,ir->ndelta);
370 if (file_version < 41) {
371 gmx_fio_do_int(fio,idum);
372 gmx_fio_do_int(fio,idum);
374 if (file_version >= 45)
375 gmx_fio_do_real(fio,ir->rtpi);
378 gmx_fio_do_int(fio,ir->nstcomm);
379 if (file_version > 34)
380 gmx_fio_do_int(fio,ir->comm_mode);
381 else if (ir->nstcomm < 0)
382 ir->comm_mode = ecmANGULAR;
384 ir->comm_mode = ecmLINEAR;
385 ir->nstcomm = abs(ir->nstcomm);
387 if(file_version > 25)
388 gmx_fio_do_int(fio,ir->nstcheckpoint);
392 gmx_fio_do_int(fio,ir->nstcgsteep);
395 gmx_fio_do_int(fio,ir->nbfgscorr);
399 gmx_fio_do_int(fio,ir->nstlog);
400 gmx_fio_do_int(fio,ir->nstxout);
401 gmx_fio_do_int(fio,ir->nstvout);
402 gmx_fio_do_int(fio,ir->nstfout);
403 gmx_fio_do_int(fio,ir->nstenergy);
404 gmx_fio_do_int(fio,ir->nstxtcout);
405 if (file_version >= 59) {
406 gmx_fio_do_double(fio,ir->init_t);
407 gmx_fio_do_double(fio,ir->delta_t);
409 gmx_fio_do_real(fio,rdum);
411 gmx_fio_do_real(fio,rdum);
414 gmx_fio_do_real(fio,ir->xtcprec);
415 if (file_version < 19) {
416 gmx_fio_do_int(fio,idum);
417 gmx_fio_do_int(fio,idum);
419 if(file_version < 18)
420 gmx_fio_do_int(fio,idum);
421 gmx_fio_do_real(fio,ir->rlist);
422 if (file_version >= 67) {
423 gmx_fio_do_real(fio,ir->rlistlong);
425 gmx_fio_do_int(fio,ir->coulombtype);
426 if (file_version < 32 && ir->coulombtype == eelRF)
427 ir->coulombtype = eelRF_NEC;
428 gmx_fio_do_real(fio,ir->rcoulomb_switch);
429 gmx_fio_do_real(fio,ir->rcoulomb);
430 gmx_fio_do_int(fio,ir->vdwtype);
431 gmx_fio_do_real(fio,ir->rvdw_switch);
432 gmx_fio_do_real(fio,ir->rvdw);
433 if (file_version < 67) {
434 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
436 gmx_fio_do_int(fio,ir->eDispCorr);
437 gmx_fio_do_real(fio,ir->epsilon_r);
438 if (file_version >= 37) {
439 gmx_fio_do_real(fio,ir->epsilon_rf);
441 if (EEL_RF(ir->coulombtype)) {
442 ir->epsilon_rf = ir->epsilon_r;
445 ir->epsilon_rf = 1.0;
448 if (file_version >= 29)
449 gmx_fio_do_real(fio,ir->tabext);
453 if(file_version > 25) {
454 gmx_fio_do_int(fio,ir->gb_algorithm);
455 gmx_fio_do_int(fio,ir->nstgbradii);
456 gmx_fio_do_real(fio,ir->rgbradii);
457 gmx_fio_do_real(fio,ir->gb_saltconc);
458 gmx_fio_do_int(fio,ir->implicit_solvent);
460 ir->gb_algorithm=egbSTILL;
464 ir->implicit_solvent=eisNO;
468 gmx_fio_do_real(fio,ir->gb_epsilon_solvent);
469 gmx_fio_do_real(fio,ir->gb_obc_alpha);
470 gmx_fio_do_real(fio,ir->gb_obc_beta);
471 gmx_fio_do_real(fio,ir->gb_obc_gamma);
474 gmx_fio_do_real(fio,ir->gb_dielectric_offset);
475 gmx_fio_do_int(fio,ir->sa_algorithm);
479 ir->gb_dielectric_offset = 0.009;
480 ir->sa_algorithm = esaAPPROX;
482 gmx_fio_do_real(fio,ir->sa_surface_tension);
484 /* Override sa_surface_tension if it is not changed in the mpd-file */
485 if(ir->sa_surface_tension<0)
487 if(ir->gb_algorithm==egbSTILL)
489 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
491 else if(ir->gb_algorithm==egbHCT || ir->gb_algorithm==egbOBC)
493 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
500 /* Better use sensible values than insane (0.0) ones... */
501 ir->gb_epsilon_solvent = 80;
502 ir->gb_obc_alpha = 1.0;
503 ir->gb_obc_beta = 0.8;
504 ir->gb_obc_gamma = 4.85;
505 ir->sa_surface_tension = 2.092;
509 gmx_fio_do_int(fio,ir->nkx);
510 gmx_fio_do_int(fio,ir->nky);
511 gmx_fio_do_int(fio,ir->nkz);
512 gmx_fio_do_int(fio,ir->pme_order);
513 gmx_fio_do_real(fio,ir->ewald_rtol);
515 if (file_version >=24)
516 gmx_fio_do_int(fio,ir->ewald_geometry);
518 ir->ewald_geometry=eewg3D;
520 if (file_version <=17) {
521 ir->epsilon_surface=0;
522 if (file_version==17)
523 gmx_fio_do_int(fio,idum);
526 gmx_fio_do_real(fio,ir->epsilon_surface);
528 gmx_fio_do_gmx_bool(fio,ir->bOptFFT);
530 gmx_fio_do_gmx_bool(fio,ir->bContinuation);
531 gmx_fio_do_int(fio,ir->etc);
532 /* before version 18, ir->etc was a gmx_bool (ir->btc),
533 * but the values 0 and 1 still mean no and
534 * berendsen temperature coupling, respectively.
536 if (file_version >= 71)
538 gmx_fio_do_int(fio,ir->nsttcouple);
542 ir->nsttcouple = ir->nstcalcenergy;
544 if (file_version <= 15)
546 gmx_fio_do_int(fio,idum);
548 if (file_version <=17)
550 gmx_fio_do_int(fio,ir->epct);
551 if (file_version <= 15)
555 ir->epct = epctSURFACETENSION;
557 gmx_fio_do_int(fio,idum);
560 /* we have removed the NO alternative at the beginning */
564 ir->epct=epctISOTROPIC;
568 ir->epc=epcBERENDSEN;
573 gmx_fio_do_int(fio,ir->epc);
574 gmx_fio_do_int(fio,ir->epct);
576 if (file_version >= 71)
578 gmx_fio_do_int(fio,ir->nstpcouple);
582 ir->nstpcouple = ir->nstcalcenergy;
584 gmx_fio_do_real(fio,ir->tau_p);
585 if (file_version <= 15) {
586 gmx_fio_do_rvec(fio,vdum);
587 clear_mat(ir->ref_p);
589 ir->ref_p[i][i] = vdum[i];
591 gmx_fio_do_rvec(fio,ir->ref_p[XX]);
592 gmx_fio_do_rvec(fio,ir->ref_p[YY]);
593 gmx_fio_do_rvec(fio,ir->ref_p[ZZ]);
595 if (file_version <= 15) {
596 gmx_fio_do_rvec(fio,vdum);
597 clear_mat(ir->compress);
599 ir->compress[i][i] = vdum[i];
602 gmx_fio_do_rvec(fio,ir->compress[XX]);
603 gmx_fio_do_rvec(fio,ir->compress[YY]);
604 gmx_fio_do_rvec(fio,ir->compress[ZZ]);
606 if (file_version >= 47) {
607 gmx_fio_do_int(fio,ir->refcoord_scaling);
608 gmx_fio_do_rvec(fio,ir->posres_com);
609 gmx_fio_do_rvec(fio,ir->posres_comB);
611 ir->refcoord_scaling = erscNO;
612 clear_rvec(ir->posres_com);
613 clear_rvec(ir->posres_comB);
615 if(file_version > 25)
616 gmx_fio_do_int(fio,ir->andersen_seed);
620 if(file_version < 26) {
621 gmx_fio_do_gmx_bool(fio,bSimAnn);
622 gmx_fio_do_real(fio,zerotemptime);
625 if (file_version < 37)
626 gmx_fio_do_real(fio,rdum);
628 gmx_fio_do_real(fio,ir->shake_tol);
629 if (file_version < 54)
630 gmx_fio_do_real(fio,*fudgeQQ);
631 gmx_fio_do_int(fio,ir->efep);
632 if (file_version <= 14 && ir->efep > efepNO)
634 if (file_version >= 59) {
635 gmx_fio_do_double(fio, ir->init_lambda);
636 gmx_fio_do_double(fio, ir->delta_lambda);
638 gmx_fio_do_real(fio,rdum);
639 ir->init_lambda = rdum;
640 gmx_fio_do_real(fio,rdum);
641 ir->delta_lambda = rdum;
643 if (file_version >= 64) {
644 gmx_fio_do_int(fio,ir->n_flambda);
646 snew(ir->flambda,ir->n_flambda);
648 bDum=gmx_fio_ndo_double(fio,ir->flambda,ir->n_flambda);
653 if (file_version >= 13)
654 gmx_fio_do_real(fio,ir->sc_alpha);
657 if (file_version >= 38)
658 gmx_fio_do_int(fio,ir->sc_power);
661 if (file_version >= 15)
662 gmx_fio_do_real(fio,ir->sc_sigma);
667 if (file_version >= 71)
669 ir->sc_sigma_min = ir->sc_sigma;
673 ir->sc_sigma_min = 0;
676 if (file_version >= 64) {
677 gmx_fio_do_int(fio,ir->nstdhdl);
682 if (file_version >= 73)
684 gmx_fio_do_int(fio, ir->separate_dhdl_file);
685 gmx_fio_do_int(fio, ir->dhdl_derivatives);
689 ir->separate_dhdl_file = sepdhdlfileYES;
690 ir->dhdl_derivatives = dhdlderivativesYES;
693 if (file_version >= 71)
695 gmx_fio_do_int(fio,ir->dh_hist_size);
696 gmx_fio_do_double(fio,ir->dh_hist_spacing);
700 ir->dh_hist_size = 0;
701 ir->dh_hist_spacing = 0.1;
703 if (file_version >= 57) {
704 gmx_fio_do_int(fio,ir->eDisre);
706 gmx_fio_do_int(fio,ir->eDisreWeighting);
707 if (file_version < 22) {
708 if (ir->eDisreWeighting == 0)
709 ir->eDisreWeighting = edrwEqual;
711 ir->eDisreWeighting = edrwConservative;
713 gmx_fio_do_gmx_bool(fio,ir->bDisreMixed);
714 gmx_fio_do_real(fio,ir->dr_fc);
715 gmx_fio_do_real(fio,ir->dr_tau);
716 gmx_fio_do_int(fio,ir->nstdisreout);
717 if (file_version >= 22) {
718 gmx_fio_do_real(fio,ir->orires_fc);
719 gmx_fio_do_real(fio,ir->orires_tau);
720 gmx_fio_do_int(fio,ir->nstorireout);
726 if(file_version >= 26) {
727 gmx_fio_do_real(fio,ir->dihre_fc);
728 if (file_version < 56) {
729 gmx_fio_do_real(fio,rdum);
730 gmx_fio_do_int(fio,idum);
736 gmx_fio_do_real(fio,ir->em_stepsize);
737 gmx_fio_do_real(fio,ir->em_tol);
738 if (file_version >= 22)
739 gmx_fio_do_gmx_bool(fio,ir->bShakeSOR);
741 ir->bShakeSOR = TRUE;
742 if (file_version >= 11)
743 gmx_fio_do_int(fio,ir->niter);
746 fprintf(stderr,"Note: niter not in run input file, setting it to %d\n",
749 if (file_version >= 21)
750 gmx_fio_do_real(fio,ir->fc_stepsize);
753 gmx_fio_do_int(fio,ir->eConstrAlg);
754 gmx_fio_do_int(fio,ir->nProjOrder);
755 gmx_fio_do_real(fio,ir->LincsWarnAngle);
756 if (file_version <= 14)
757 gmx_fio_do_int(fio,idum);
758 if (file_version >=26)
759 gmx_fio_do_int(fio,ir->nLincsIter);
762 fprintf(stderr,"Note: nLincsIter not in run input file, setting it to %d\n",
765 if (file_version < 33)
766 gmx_fio_do_real(fio,bd_temp);
767 gmx_fio_do_real(fio,ir->bd_fric);
768 gmx_fio_do_int(fio,ir->ld_seed);
769 if (file_version >= 33) {
771 gmx_fio_do_rvec(fio,ir->deform[i]);
774 clear_rvec(ir->deform[i]);
776 if (file_version >= 14)
777 gmx_fio_do_real(fio,ir->cos_accel);
780 gmx_fio_do_int(fio,ir->userint1);
781 gmx_fio_do_int(fio,ir->userint2);
782 gmx_fio_do_int(fio,ir->userint3);
783 gmx_fio_do_int(fio,ir->userint4);
784 gmx_fio_do_real(fio,ir->userreal1);
785 gmx_fio_do_real(fio,ir->userreal2);
786 gmx_fio_do_real(fio,ir->userreal3);
787 gmx_fio_do_real(fio,ir->userreal4);
790 if (file_version >= 75) {
791 gmx_fio_do_gmx_bool(fio,ir->bAdress);
793 if (bRead) snew(ir->adress, 1);
794 gmx_fio_do_int(fio,ir->adress->type);
795 gmx_fio_do_real(fio,ir->adress->const_wf);
796 gmx_fio_do_real(fio,ir->adress->ex_width);
797 gmx_fio_do_real(fio,ir->adress->hy_width);
798 gmx_fio_do_int(fio,ir->adress->icor);
799 gmx_fio_do_int(fio,ir->adress->site);
800 gmx_fio_do_rvec(fio,ir->adress->refs);
801 gmx_fio_do_int(fio,ir->adress->n_tf_grps);
802 gmx_fio_do_real(fio, ir->adress->ex_forcecap);
803 gmx_fio_do_int(fio, ir->adress->n_energy_grps);
804 gmx_fio_do_int(fio,ir->adress->do_hybridpairs);
806 if (bRead)snew(ir->adress->tf_table_index,ir->adress->n_tf_grps);
807 if (ir->adress->n_tf_grps > 0) {
808 bDum=gmx_fio_ndo_int(fio,ir->adress->tf_table_index,ir->adress->n_tf_grps);
810 if (bRead)snew(ir->adress->group_explicit,ir->adress->n_energy_grps);
811 if (ir->adress->n_energy_grps > 0) {
812 bDum=gmx_fio_ndo_int(fio, ir->adress->group_explicit,ir->adress->n_energy_grps);
820 if (file_version >= 48) {
821 gmx_fio_do_int(fio,ir->ePull);
822 if (ir->ePull != epullNO) {
825 do_pull(fio, ir->pull,bRead,file_version);
831 /* Enforced rotation */
832 if (file_version >= 74) {
833 gmx_fio_do_int(fio,ir->bRot);
834 if (ir->bRot == TRUE) {
837 do_rot(fio, ir->rot,bRead,file_version);
844 gmx_fio_do_int(fio,ir->opts.ngtc);
845 if (file_version >= 69) {
846 gmx_fio_do_int(fio,ir->opts.nhchainlength);
848 ir->opts.nhchainlength = 1;
850 gmx_fio_do_int(fio,ir->opts.ngacc);
851 gmx_fio_do_int(fio,ir->opts.ngfrz);
852 gmx_fio_do_int(fio,ir->opts.ngener);
855 snew(ir->opts.nrdf, ir->opts.ngtc);
856 snew(ir->opts.ref_t, ir->opts.ngtc);
857 snew(ir->opts.annealing, ir->opts.ngtc);
858 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
859 snew(ir->opts.anneal_time, ir->opts.ngtc);
860 snew(ir->opts.anneal_temp, ir->opts.ngtc);
861 snew(ir->opts.tau_t, ir->opts.ngtc);
862 snew(ir->opts.nFreeze,ir->opts.ngfrz);
863 snew(ir->opts.acc, ir->opts.ngacc);
864 snew(ir->opts.egp_flags,ir->opts.ngener*ir->opts.ngener);
866 if (ir->opts.ngtc > 0) {
867 if (bRead && file_version<13) {
868 snew(tmp,ir->opts.ngtc);
869 bDum=gmx_fio_ndo_int(fio,tmp, ir->opts.ngtc);
870 for(i=0; i<ir->opts.ngtc; i++)
871 ir->opts.nrdf[i] = tmp[i];
874 bDum=gmx_fio_ndo_real(fio,ir->opts.nrdf, ir->opts.ngtc);
876 bDum=gmx_fio_ndo_real(fio,ir->opts.ref_t,ir->opts.ngtc);
877 bDum=gmx_fio_ndo_real(fio,ir->opts.tau_t,ir->opts.ngtc);
878 if (file_version<33 && ir->eI==eiBD) {
879 for(i=0; i<ir->opts.ngtc; i++)
880 ir->opts.tau_t[i] = bd_temp;
883 if (ir->opts.ngfrz > 0)
884 bDum=gmx_fio_ndo_ivec(fio,ir->opts.nFreeze,ir->opts.ngfrz);
885 if (ir->opts.ngacc > 0)
886 gmx_fio_ndo_rvec(fio,ir->opts.acc,ir->opts.ngacc);
887 if (file_version >= 12)
888 bDum=gmx_fio_ndo_int(fio,ir->opts.egp_flags,
889 ir->opts.ngener*ir->opts.ngener);
891 if(bRead && file_version < 26) {
892 for(i=0;i<ir->opts.ngtc;i++) {
894 ir->opts.annealing[i] = eannSINGLE;
895 ir->opts.anneal_npoints[i] = 2;
896 snew(ir->opts.anneal_time[i],2);
897 snew(ir->opts.anneal_temp[i],2);
898 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
899 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
900 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
901 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
902 ir->opts.anneal_time[i][0] = ir->init_t;
903 ir->opts.anneal_time[i][1] = finish_t;
904 ir->opts.anneal_temp[i][0] = init_temp;
905 ir->opts.anneal_temp[i][1] = finish_temp;
907 ir->opts.annealing[i] = eannNO;
908 ir->opts.anneal_npoints[i] = 0;
912 /* file version 26 or later */
913 /* First read the lists with annealing and npoints for each group */
914 bDum=gmx_fio_ndo_int(fio,ir->opts.annealing,ir->opts.ngtc);
915 bDum=gmx_fio_ndo_int(fio,ir->opts.anneal_npoints,ir->opts.ngtc);
916 for(j=0;j<(ir->opts.ngtc);j++) {
917 k=ir->opts.anneal_npoints[j];
919 snew(ir->opts.anneal_time[j],k);
920 snew(ir->opts.anneal_temp[j],k);
922 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_time[j],k);
923 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_temp[j],k);
927 if (file_version >= 45) {
928 gmx_fio_do_int(fio,ir->nwall);
929 gmx_fio_do_int(fio,ir->wall_type);
930 if (file_version >= 50)
931 gmx_fio_do_real(fio,ir->wall_r_linpot);
933 ir->wall_r_linpot = -1;
934 gmx_fio_do_int(fio,ir->wall_atomtype[0]);
935 gmx_fio_do_int(fio,ir->wall_atomtype[1]);
936 gmx_fio_do_real(fio,ir->wall_density[0]);
937 gmx_fio_do_real(fio,ir->wall_density[1]);
938 gmx_fio_do_real(fio,ir->wall_ewald_zfac);
942 ir->wall_atomtype[0] = -1;
943 ir->wall_atomtype[1] = -1;
944 ir->wall_density[0] = 0;
945 ir->wall_density[1] = 0;
946 ir->wall_ewald_zfac = 3;
948 /* Cosine stuff for electric fields */
949 for(j=0; (j<DIM); j++) {
950 gmx_fio_do_int(fio,ir->ex[j].n);
951 gmx_fio_do_int(fio,ir->et[j].n);
953 snew(ir->ex[j].a, ir->ex[j].n);
954 snew(ir->ex[j].phi,ir->ex[j].n);
955 snew(ir->et[j].a, ir->et[j].n);
956 snew(ir->et[j].phi,ir->et[j].n);
958 bDum=gmx_fio_ndo_real(fio,ir->ex[j].a, ir->ex[j].n);
959 bDum=gmx_fio_ndo_real(fio,ir->ex[j].phi,ir->ex[j].n);
960 bDum=gmx_fio_ndo_real(fio,ir->et[j].a, ir->et[j].n);
961 bDum=gmx_fio_ndo_real(fio,ir->et[j].phi,ir->et[j].n);
965 if(file_version>=39){
966 gmx_fio_do_gmx_bool(fio,ir->bQMMM);
967 gmx_fio_do_int(fio,ir->QMMMscheme);
968 gmx_fio_do_real(fio,ir->scalefactor);
969 gmx_fio_do_int(fio,ir->opts.ngQM);
971 snew(ir->opts.QMmethod, ir->opts.ngQM);
972 snew(ir->opts.QMbasis, ir->opts.ngQM);
973 snew(ir->opts.QMcharge, ir->opts.ngQM);
974 snew(ir->opts.QMmult, ir->opts.ngQM);
975 snew(ir->opts.bSH, ir->opts.ngQM);
976 snew(ir->opts.CASorbitals, ir->opts.ngQM);
977 snew(ir->opts.CASelectrons,ir->opts.ngQM);
978 snew(ir->opts.SAon, ir->opts.ngQM);
979 snew(ir->opts.SAoff, ir->opts.ngQM);
980 snew(ir->opts.SAsteps, ir->opts.ngQM);
981 snew(ir->opts.bOPT, ir->opts.ngQM);
982 snew(ir->opts.bTS, ir->opts.ngQM);
984 if (ir->opts.ngQM > 0) {
985 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmethod,ir->opts.ngQM);
986 bDum=gmx_fio_ndo_int(fio,ir->opts.QMbasis,ir->opts.ngQM);
987 bDum=gmx_fio_ndo_int(fio,ir->opts.QMcharge,ir->opts.ngQM);
988 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmult,ir->opts.ngQM);
989 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bSH,ir->opts.ngQM);
990 bDum=gmx_fio_ndo_int(fio,ir->opts.CASorbitals,ir->opts.ngQM);
991 bDum=gmx_fio_ndo_int(fio,ir->opts.CASelectrons,ir->opts.ngQM);
992 bDum=gmx_fio_ndo_real(fio,ir->opts.SAon,ir->opts.ngQM);
993 bDum=gmx_fio_ndo_real(fio,ir->opts.SAoff,ir->opts.ngQM);
994 bDum=gmx_fio_ndo_int(fio,ir->opts.SAsteps,ir->opts.ngQM);
995 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bOPT,ir->opts.ngQM);
996 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bTS,ir->opts.ngQM);
998 /* end of QMMM stuff */
1003 static void do_harm(t_fileio *fio, t_iparams *iparams,gmx_bool bRead)
1005 gmx_fio_do_real(fio,iparams->harmonic.rA);
1006 gmx_fio_do_real(fio,iparams->harmonic.krA);
1007 gmx_fio_do_real(fio,iparams->harmonic.rB);
1008 gmx_fio_do_real(fio,iparams->harmonic.krB);
1011 void do_iparams(t_fileio *fio, t_functype ftype,t_iparams *iparams,
1012 gmx_bool bRead, int file_version)
1019 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1027 do_harm(fio, iparams,bRead);
1028 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead) {
1029 /* Correct incorrect storage of parameters */
1030 iparams->pdihs.phiB = iparams->pdihs.phiA;
1031 iparams->pdihs.cpB = iparams->pdihs.cpA;
1035 gmx_fio_do_real(fio,iparams->fene.bm);
1036 gmx_fio_do_real(fio,iparams->fene.kb);
1039 gmx_fio_do_real(fio,iparams->restraint.lowA);
1040 gmx_fio_do_real(fio,iparams->restraint.up1A);
1041 gmx_fio_do_real(fio,iparams->restraint.up2A);
1042 gmx_fio_do_real(fio,iparams->restraint.kA);
1043 gmx_fio_do_real(fio,iparams->restraint.lowB);
1044 gmx_fio_do_real(fio,iparams->restraint.up1B);
1045 gmx_fio_do_real(fio,iparams->restraint.up2B);
1046 gmx_fio_do_real(fio,iparams->restraint.kB);
1052 gmx_fio_do_real(fio,iparams->tab.kA);
1053 gmx_fio_do_int(fio,iparams->tab.table);
1054 gmx_fio_do_real(fio,iparams->tab.kB);
1056 case F_CROSS_BOND_BONDS:
1057 gmx_fio_do_real(fio,iparams->cross_bb.r1e);
1058 gmx_fio_do_real(fio,iparams->cross_bb.r2e);
1059 gmx_fio_do_real(fio,iparams->cross_bb.krr);
1061 case F_CROSS_BOND_ANGLES:
1062 gmx_fio_do_real(fio,iparams->cross_ba.r1e);
1063 gmx_fio_do_real(fio,iparams->cross_ba.r2e);
1064 gmx_fio_do_real(fio,iparams->cross_ba.r3e);
1065 gmx_fio_do_real(fio,iparams->cross_ba.krt);
1067 case F_UREY_BRADLEY:
1068 gmx_fio_do_real(fio,iparams->u_b.theta);
1069 gmx_fio_do_real(fio,iparams->u_b.ktheta);
1070 gmx_fio_do_real(fio,iparams->u_b.r13);
1071 gmx_fio_do_real(fio,iparams->u_b.kUB);
1073 case F_QUARTIC_ANGLES:
1074 gmx_fio_do_real(fio,iparams->qangle.theta);
1075 bDum=gmx_fio_ndo_real(fio,iparams->qangle.c,5);
1078 gmx_fio_do_real(fio,iparams->bham.a);
1079 gmx_fio_do_real(fio,iparams->bham.b);
1080 gmx_fio_do_real(fio,iparams->bham.c);
1083 gmx_fio_do_real(fio,iparams->morse.b0);
1084 gmx_fio_do_real(fio,iparams->morse.cb);
1085 gmx_fio_do_real(fio,iparams->morse.beta);
1088 gmx_fio_do_real(fio,iparams->cubic.b0);
1089 gmx_fio_do_real(fio,iparams->cubic.kb);
1090 gmx_fio_do_real(fio,iparams->cubic.kcub);
1094 case F_POLARIZATION:
1095 gmx_fio_do_real(fio,iparams->polarize.alpha);
1098 if (file_version < 31)
1099 gmx_fatal(FARGS,"Old tpr files with water_polarization not supported. Make a new.");
1100 gmx_fio_do_real(fio,iparams->wpol.al_x);
1101 gmx_fio_do_real(fio,iparams->wpol.al_y);
1102 gmx_fio_do_real(fio,iparams->wpol.al_z);
1103 gmx_fio_do_real(fio,iparams->wpol.rOH);
1104 gmx_fio_do_real(fio,iparams->wpol.rHH);
1105 gmx_fio_do_real(fio,iparams->wpol.rOD);
1108 gmx_fio_do_real(fio,iparams->thole.a);
1109 gmx_fio_do_real(fio,iparams->thole.alpha1);
1110 gmx_fio_do_real(fio,iparams->thole.alpha2);
1111 gmx_fio_do_real(fio,iparams->thole.rfac);
1114 gmx_fio_do_real(fio,iparams->lj.c6);
1115 gmx_fio_do_real(fio,iparams->lj.c12);
1118 gmx_fio_do_real(fio,iparams->lj14.c6A);
1119 gmx_fio_do_real(fio,iparams->lj14.c12A);
1120 gmx_fio_do_real(fio,iparams->lj14.c6B);
1121 gmx_fio_do_real(fio,iparams->lj14.c12B);
1124 gmx_fio_do_real(fio,iparams->ljc14.fqq);
1125 gmx_fio_do_real(fio,iparams->ljc14.qi);
1126 gmx_fio_do_real(fio,iparams->ljc14.qj);
1127 gmx_fio_do_real(fio,iparams->ljc14.c6);
1128 gmx_fio_do_real(fio,iparams->ljc14.c12);
1130 case F_LJC_PAIRS_NB:
1131 gmx_fio_do_real(fio,iparams->ljcnb.qi);
1132 gmx_fio_do_real(fio,iparams->ljcnb.qj);
1133 gmx_fio_do_real(fio,iparams->ljcnb.c6);
1134 gmx_fio_do_real(fio,iparams->ljcnb.c12);
1140 gmx_fio_do_real(fio,iparams->pdihs.phiA);
1141 gmx_fio_do_real(fio,iparams->pdihs.cpA);
1142 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42) {
1143 /* Read the incorrectly stored multiplicity */
1144 gmx_fio_do_real(fio,iparams->harmonic.rB);
1145 gmx_fio_do_real(fio,iparams->harmonic.krB);
1146 iparams->pdihs.phiB = iparams->pdihs.phiA;
1147 iparams->pdihs.cpB = iparams->pdihs.cpA;
1149 gmx_fio_do_real(fio,iparams->pdihs.phiB);
1150 gmx_fio_do_real(fio,iparams->pdihs.cpB);
1151 gmx_fio_do_int(fio,iparams->pdihs.mult);
1155 gmx_fio_do_int(fio,iparams->disres.label);
1156 gmx_fio_do_int(fio,iparams->disres.type);
1157 gmx_fio_do_real(fio,iparams->disres.low);
1158 gmx_fio_do_real(fio,iparams->disres.up1);
1159 gmx_fio_do_real(fio,iparams->disres.up2);
1160 gmx_fio_do_real(fio,iparams->disres.kfac);
1163 gmx_fio_do_int(fio,iparams->orires.ex);
1164 gmx_fio_do_int(fio,iparams->orires.label);
1165 gmx_fio_do_int(fio,iparams->orires.power);
1166 gmx_fio_do_real(fio,iparams->orires.c);
1167 gmx_fio_do_real(fio,iparams->orires.obs);
1168 gmx_fio_do_real(fio,iparams->orires.kfac);
1171 gmx_fio_do_int(fio,iparams->dihres.power);
1172 gmx_fio_do_int(fio,iparams->dihres.label);
1173 gmx_fio_do_real(fio,iparams->dihres.phi);
1174 gmx_fio_do_real(fio,iparams->dihres.dphi);
1175 gmx_fio_do_real(fio,iparams->dihres.kfac);
1178 gmx_fio_do_rvec(fio,iparams->posres.pos0A);
1179 gmx_fio_do_rvec(fio,iparams->posres.fcA);
1180 if (bRead && file_version < 27) {
1181 copy_rvec(iparams->posres.pos0A,iparams->posres.pos0B);
1182 copy_rvec(iparams->posres.fcA,iparams->posres.fcB);
1184 gmx_fio_do_rvec(fio,iparams->posres.pos0B);
1185 gmx_fio_do_rvec(fio,iparams->posres.fcB);
1189 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA,NR_RBDIHS);
1190 if(file_version>=25)
1191 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB,NR_RBDIHS);
1194 /* Fourier dihedrals are internally represented
1195 * as Ryckaert-Bellemans since those are faster to compute.
1197 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA, NR_RBDIHS);
1198 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB, NR_RBDIHS);
1202 gmx_fio_do_real(fio,iparams->constr.dA);
1203 gmx_fio_do_real(fio,iparams->constr.dB);
1206 gmx_fio_do_real(fio,iparams->settle.doh);
1207 gmx_fio_do_real(fio,iparams->settle.dhh);
1210 gmx_fio_do_real(fio,iparams->vsite.a);
1215 gmx_fio_do_real(fio,iparams->vsite.a);
1216 gmx_fio_do_real(fio,iparams->vsite.b);
1221 gmx_fio_do_real(fio,iparams->vsite.a);
1222 gmx_fio_do_real(fio,iparams->vsite.b);
1223 gmx_fio_do_real(fio,iparams->vsite.c);
1226 gmx_fio_do_int(fio,iparams->vsiten.n);
1227 gmx_fio_do_real(fio,iparams->vsiten.a);
1232 /* We got rid of some parameters in version 68 */
1233 if(bRead && file_version<68)
1235 gmx_fio_do_real(fio,rdum);
1236 gmx_fio_do_real(fio,rdum);
1237 gmx_fio_do_real(fio,rdum);
1238 gmx_fio_do_real(fio,rdum);
1240 gmx_fio_do_real(fio,iparams->gb.sar);
1241 gmx_fio_do_real(fio,iparams->gb.st);
1242 gmx_fio_do_real(fio,iparams->gb.pi);
1243 gmx_fio_do_real(fio,iparams->gb.gbr);
1244 gmx_fio_do_real(fio,iparams->gb.bmlt);
1247 gmx_fio_do_int(fio,iparams->cmap.cmapA);
1248 gmx_fio_do_int(fio,iparams->cmap.cmapB);
1251 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
1253 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
1256 gmx_fio_unset_comment(fio);
1259 static void do_ilist(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,int file_version,
1266 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1268 if (file_version < 44) {
1269 for(i=0; i<MAXNODES; i++)
1270 gmx_fio_do_int(fio,idum);
1272 gmx_fio_do_int(fio,ilist->nr);
1274 snew(ilist->iatoms,ilist->nr);
1275 bDum=gmx_fio_ndo_int(fio,ilist->iatoms,ilist->nr);
1277 gmx_fio_unset_comment(fio);
1280 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1281 gmx_bool bRead, int file_version)
1287 gmx_fio_do_int(fio,ffparams->atnr);
1288 if (file_version < 57) {
1289 gmx_fio_do_int(fio,idum);
1291 gmx_fio_do_int(fio,ffparams->ntypes);
1293 fprintf(debug,"ffparams->atnr = %d, ntypes = %d\n",
1294 ffparams->atnr,ffparams->ntypes);
1296 snew(ffparams->functype,ffparams->ntypes);
1297 snew(ffparams->iparams,ffparams->ntypes);
1299 /* Read/write all the function types */
1300 bDum=gmx_fio_ndo_int(fio,ffparams->functype,ffparams->ntypes);
1302 pr_ivec(debug,0,"functype",ffparams->functype,ffparams->ntypes,TRUE);
1304 if (file_version >= 66) {
1305 gmx_fio_do_double(fio,ffparams->reppow);
1307 ffparams->reppow = 12.0;
1310 if (file_version >= 57) {
1311 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1314 /* Check whether all these function types are supported by the code.
1315 * In practice the code is backwards compatible, which means that the
1316 * numbering may have to be altered from old numbering to new numbering
1318 for (i=0; (i<ffparams->ntypes); i++) {
1320 /* Loop over file versions */
1321 for (k=0; (k<NFTUPD); k++)
1322 /* Compare the read file_version to the update table */
1323 if ((file_version < ftupd[k].fvnr) &&
1324 (ffparams->functype[i] >= ftupd[k].ftype)) {
1325 ffparams->functype[i] += 1;
1327 fprintf(debug,"Incrementing function type %d to %d (due to %s)\n",
1328 i,ffparams->functype[i],
1329 interaction_function[ftupd[k].ftype].longname);
1334 do_iparams(fio, ffparams->functype[i],&ffparams->iparams[i],bRead,
1337 pr_iparams(debug,ffparams->functype[i],&ffparams->iparams[i]);
1341 static void do_ilists(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,
1344 int i,j,renum[F_NRE];
1345 gmx_bool bDum=TRUE,bClear;
1348 for(j=0; (j<F_NRE); j++) {
1351 for (k=0; k<NFTUPD; k++)
1352 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
1356 ilist[j].iatoms = NULL;
1358 do_ilist(fio, &ilist[j],bRead,file_version,j);
1361 if (bRead && gmx_debug_at)
1362 pr_ilist(debug,0,interaction_function[j].longname,
1363 functype,&ilist[j],TRUE);
1368 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams,gmx_moltype_t *molt,
1369 gmx_bool bRead, int file_version)
1371 do_ffparams(fio, ffparams,bRead,file_version);
1373 if (file_version >= 54) {
1374 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1377 do_ilists(fio, molt->ilist,bRead,file_version);
1380 static void do_block(t_fileio *fio, t_block *block,gmx_bool bRead,int file_version)
1382 int i,idum,dum_nra,*dum_a;
1385 if (file_version < 44)
1386 for(i=0; i<MAXNODES; i++)
1387 gmx_fio_do_int(fio,idum);
1388 gmx_fio_do_int(fio,block->nr);
1389 if (file_version < 51)
1390 gmx_fio_do_int(fio,dum_nra);
1392 block->nalloc_index = block->nr+1;
1393 snew(block->index,block->nalloc_index);
1395 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1397 if (file_version < 51 && dum_nra > 0) {
1398 snew(dum_a,dum_nra);
1399 bDum=gmx_fio_ndo_int(fio,dum_a,dum_nra);
1404 static void do_blocka(t_fileio *fio, t_blocka *block,gmx_bool bRead,
1410 if (file_version < 44)
1411 for(i=0; i<MAXNODES; i++)
1412 gmx_fio_do_int(fio,idum);
1413 gmx_fio_do_int(fio,block->nr);
1414 gmx_fio_do_int(fio,block->nra);
1416 block->nalloc_index = block->nr+1;
1417 snew(block->index,block->nalloc_index);
1418 block->nalloc_a = block->nra;
1419 snew(block->a,block->nalloc_a);
1421 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1422 bDum=gmx_fio_ndo_int(fio,block->a,block->nra);
1425 static void do_atom(t_fileio *fio, t_atom *atom,int ngrp,gmx_bool bRead,
1426 int file_version, gmx_groups_t *groups,int atnr)
1430 gmx_fio_do_real(fio,atom->m);
1431 gmx_fio_do_real(fio,atom->q);
1432 gmx_fio_do_real(fio,atom->mB);
1433 gmx_fio_do_real(fio,atom->qB);
1434 gmx_fio_do_ushort(fio, atom->type);
1435 gmx_fio_do_ushort(fio, atom->typeB);
1436 gmx_fio_do_int(fio,atom->ptype);
1437 gmx_fio_do_int(fio,atom->resind);
1438 if (file_version >= 52)
1439 gmx_fio_do_int(fio,atom->atomnumber);
1441 atom->atomnumber = NOTSET;
1442 if (file_version < 23)
1444 else if (file_version < 39)
1449 if (file_version < 57) {
1450 unsigned char uchar[egcNR];
1451 gmx_fio_ndo_uchar(fio,uchar,myngrp);
1452 for(i=myngrp; (i<ngrp); i++) {
1455 /* Copy the old data format to the groups struct */
1456 for(i=0; i<ngrp; i++) {
1457 groups->grpnr[i][atnr] = uchar[i];
1462 static void do_grps(t_fileio *fio, int ngrp,t_grps grps[],gmx_bool bRead,
1468 if (file_version < 23)
1470 else if (file_version < 39)
1475 for(j=0; (j<ngrp); j++) {
1477 gmx_fio_do_int(fio,grps[j].nr);
1479 snew(grps[j].nm_ind,grps[j].nr);
1480 bDum=gmx_fio_ndo_int(fio,grps[j].nm_ind,grps[j].nr);
1484 snew(grps[j].nm_ind,grps[j].nr);
1489 static void do_symstr(t_fileio *fio, char ***nm,gmx_bool bRead,t_symtab *symtab)
1494 gmx_fio_do_int(fio,ls);
1495 *nm = get_symtab_handle(symtab,ls);
1498 ls = lookup_symtab(symtab,*nm);
1499 gmx_fio_do_int(fio,ls);
1503 static void do_strstr(t_fileio *fio, int nstr,char ***nm,gmx_bool bRead,
1508 for (j=0; (j<nstr); j++)
1509 do_symstr(fio, &(nm[j]),bRead,symtab);
1512 static void do_resinfo(t_fileio *fio, int n,t_resinfo *ri,gmx_bool bRead,
1513 t_symtab *symtab, int file_version)
1517 for (j=0; (j<n); j++) {
1518 do_symstr(fio, &(ri[j].name),bRead,symtab);
1519 if (file_version >= 63) {
1520 gmx_fio_do_int(fio,ri[j].nr);
1521 gmx_fio_do_uchar(fio, ri[j].ic);
1529 static void do_atoms(t_fileio *fio, t_atoms *atoms,gmx_bool bRead,t_symtab *symtab,
1531 gmx_groups_t *groups)
1535 gmx_fio_do_int(fio,atoms->nr);
1536 gmx_fio_do_int(fio,atoms->nres);
1537 if (file_version < 57) {
1538 gmx_fio_do_int(fio,groups->ngrpname);
1539 for(i=0; i<egcNR; i++) {
1540 groups->ngrpnr[i] = atoms->nr;
1541 snew(groups->grpnr[i],groups->ngrpnr[i]);
1545 snew(atoms->atom,atoms->nr);
1546 snew(atoms->atomname,atoms->nr);
1547 snew(atoms->atomtype,atoms->nr);
1548 snew(atoms->atomtypeB,atoms->nr);
1549 snew(atoms->resinfo,atoms->nres);
1550 if (file_version < 57) {
1551 snew(groups->grpname,groups->ngrpname);
1553 atoms->pdbinfo = NULL;
1555 for(i=0; (i<atoms->nr); i++) {
1556 do_atom(fio, &atoms->atom[i],egcNR,bRead, file_version,groups,i);
1558 do_strstr(fio, atoms->nr,atoms->atomname,bRead,symtab);
1559 if (bRead && (file_version <= 20)) {
1560 for(i=0; i<atoms->nr; i++) {
1561 atoms->atomtype[i] = put_symtab(symtab,"?");
1562 atoms->atomtypeB[i] = put_symtab(symtab,"?");
1565 do_strstr(fio, atoms->nr,atoms->atomtype,bRead,symtab);
1566 do_strstr(fio, atoms->nr,atoms->atomtypeB,bRead,symtab);
1568 do_resinfo(fio, atoms->nres,atoms->resinfo,bRead,symtab,file_version);
1570 if (file_version < 57) {
1571 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1573 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1577 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
1578 gmx_bool bRead,t_symtab *symtab,
1584 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1585 gmx_fio_do_int(fio,groups->ngrpname);
1587 snew(groups->grpname,groups->ngrpname);
1589 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1590 for(g=0; g<egcNR; g++) {
1591 gmx_fio_do_int(fio,groups->ngrpnr[g]);
1592 if (groups->ngrpnr[g] == 0) {
1594 groups->grpnr[g] = NULL;
1598 snew(groups->grpnr[g],groups->ngrpnr[g]);
1600 bDum=gmx_fio_ndo_uchar(fio, groups->grpnr[g],groups->ngrpnr[g]);
1605 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes,gmx_bool bRead,
1606 t_symtab *symtab,int file_version)
1609 gmx_bool bDum = TRUE;
1611 if (file_version > 25) {
1612 gmx_fio_do_int(fio,atomtypes->nr);
1615 snew(atomtypes->radius,j);
1616 snew(atomtypes->vol,j);
1617 snew(atomtypes->surftens,j);
1618 snew(atomtypes->atomnumber,j);
1619 snew(atomtypes->gb_radius,j);
1620 snew(atomtypes->S_hct,j);
1622 bDum=gmx_fio_ndo_real(fio,atomtypes->radius,j);
1623 bDum=gmx_fio_ndo_real(fio,atomtypes->vol,j);
1624 bDum=gmx_fio_ndo_real(fio,atomtypes->surftens,j);
1625 if(file_version >= 40)
1627 bDum=gmx_fio_ndo_int(fio,atomtypes->atomnumber,j);
1629 if(file_version >= 60)
1631 bDum=gmx_fio_ndo_real(fio,atomtypes->gb_radius,j);
1632 bDum=gmx_fio_ndo_real(fio,atomtypes->S_hct,j);
1635 /* File versions prior to 26 cannot do GBSA,
1636 * so they dont use this structure
1639 atomtypes->radius = NULL;
1640 atomtypes->vol = NULL;
1641 atomtypes->surftens = NULL;
1642 atomtypes->atomnumber = NULL;
1643 atomtypes->gb_radius = NULL;
1644 atomtypes->S_hct = NULL;
1648 static void do_symtab(t_fileio *fio, t_symtab *symtab,gmx_bool bRead)
1654 gmx_fio_do_int(fio,symtab->nr);
1657 snew(symtab->symbuf,1);
1658 symbuf = symtab->symbuf;
1659 symbuf->bufsize = nr;
1660 snew(symbuf->buf,nr);
1661 for (i=0; (i<nr); i++) {
1662 gmx_fio_do_string(fio,buf);
1663 symbuf->buf[i]=strdup(buf);
1667 symbuf = symtab->symbuf;
1668 while (symbuf!=NULL) {
1669 for (i=0; (i<symbuf->bufsize) && (i<nr); i++)
1670 gmx_fio_do_string(fio,symbuf->buf[i]);
1672 symbuf=symbuf->next;
1675 gmx_fatal(FARGS,"nr of symtab strings left: %d",nr);
1679 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
1681 int i,j,ngrid,gs,nelem;
1683 gmx_fio_do_int(fio,cmap_grid->ngrid);
1684 gmx_fio_do_int(fio,cmap_grid->grid_spacing);
1686 ngrid = cmap_grid->ngrid;
1687 gs = cmap_grid->grid_spacing;
1692 snew(cmap_grid->cmapdata,ngrid);
1694 for(i=0;i<cmap_grid->ngrid;i++)
1696 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
1700 for(i=0;i<cmap_grid->ngrid;i++)
1702 for(j=0;j<nelem;j++)
1704 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4]);
1705 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+1]);
1706 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+2]);
1707 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+3]);
1713 void tpx_make_chain_identifiers(t_atoms *atoms,t_block *mols)
1719 /* We always assign a new chain number, but save the chain id characters
1720 * for larger molecules.
1722 #define CHAIN_MIN_ATOMS 15
1726 for(m=0; m<mols->nr; m++)
1729 a1=mols->index[m+1];
1730 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
1739 for(a=a0; a<a1; a++)
1741 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
1742 atoms->resinfo[atoms->atom[a].resind].chainid = c;
1747 /* Blank out the chain id if there was only one chain */
1750 for(r=0; r<atoms->nres; r++)
1752 atoms->resinfo[r].chainid = ' ';
1757 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt,gmx_bool bRead,
1758 t_symtab *symtab, int file_version,
1759 gmx_groups_t *groups)
1763 if (file_version >= 57) {
1764 do_symstr(fio, &(molt->name),bRead,symtab);
1767 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
1769 if (bRead && gmx_debug_at) {
1770 pr_atoms(debug,0,"atoms",&molt->atoms,TRUE);
1773 if (file_version >= 57) {
1774 do_ilists(fio, molt->ilist,bRead,file_version);
1776 do_block(fio, &molt->cgs,bRead,file_version);
1777 if (bRead && gmx_debug_at) {
1778 pr_block(debug,0,"cgs",&molt->cgs,TRUE);
1782 /* This used to be in the atoms struct */
1783 do_blocka(fio, &molt->excls, bRead, file_version);
1786 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb,gmx_bool bRead,
1791 gmx_fio_do_int(fio,molb->type);
1792 gmx_fio_do_int(fio,molb->nmol);
1793 gmx_fio_do_int(fio,molb->natoms_mol);
1794 /* Position restraint coordinates */
1795 gmx_fio_do_int(fio,molb->nposres_xA);
1796 if (molb->nposres_xA > 0) {
1798 snew(molb->posres_xA,molb->nposres_xA);
1800 gmx_fio_ndo_rvec(fio,molb->posres_xA,molb->nposres_xA);
1802 gmx_fio_do_int(fio,molb->nposres_xB);
1803 if (molb->nposres_xB > 0) {
1805 snew(molb->posres_xB,molb->nposres_xB);
1807 gmx_fio_ndo_rvec(fio,molb->posres_xB,molb->nposres_xB);
1812 static t_block mtop_mols(gmx_mtop_t *mtop)
1818 for(mb=0; mb<mtop->nmolblock; mb++) {
1819 mols.nr += mtop->molblock[mb].nmol;
1821 mols.nalloc_index = mols.nr + 1;
1822 snew(mols.index,mols.nalloc_index);
1827 for(mb=0; mb<mtop->nmolblock; mb++) {
1828 for(mol=0; mol<mtop->molblock[mb].nmol; mol++) {
1829 a += mtop->molblock[mb].natoms_mol;
1838 static void add_posres_molblock(gmx_mtop_t *mtop)
1843 gmx_molblock_t *molb;
1846 il = &mtop->moltype[0].ilist[F_POSRES];
1852 for(i=0; i<il->nr; i+=2) {
1853 ip = &mtop->ffparams.iparams[il->iatoms[i]];
1854 am = max(am,il->iatoms[i+1]);
1855 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
1856 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
1857 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ]) {
1861 /* Make the posres coordinate block end at a molecule end */
1863 while(am >= mtop->mols.index[mol+1]) {
1866 molb = &mtop->molblock[0];
1867 molb->nposres_xA = mtop->mols.index[mol+1];
1868 snew(molb->posres_xA,molb->nposres_xA);
1870 molb->nposres_xB = molb->nposres_xA;
1871 snew(molb->posres_xB,molb->nposres_xB);
1873 molb->nposres_xB = 0;
1875 for(i=0; i<il->nr; i+=2) {
1876 ip = &mtop->ffparams.iparams[il->iatoms[i]];
1877 a = il->iatoms[i+1];
1878 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
1879 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
1880 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
1882 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
1883 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
1884 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
1889 static void set_disres_npair(gmx_mtop_t *mtop)
1896 ip = mtop->ffparams.iparams;
1898 for(mt=0; mt<mtop->nmoltype; mt++) {
1899 il = &mtop->moltype[mt].ilist[F_DISRES];
1903 for(i=0; i<il->nr; i+=3) {
1905 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label) {
1906 ip[a[i]].disres.npair = npair;
1914 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop,gmx_bool bRead,
1922 do_symtab(fio, &(mtop->symtab),bRead);
1924 pr_symtab(debug,0,"symtab",&mtop->symtab);
1926 do_symstr(fio, &(mtop->name),bRead,&(mtop->symtab));
1928 if (file_version >= 57) {
1929 do_ffparams(fio, &mtop->ffparams,bRead,file_version);
1931 gmx_fio_do_int(fio,mtop->nmoltype);
1936 snew(mtop->moltype,mtop->nmoltype);
1937 if (file_version < 57) {
1938 mtop->moltype[0].name = mtop->name;
1941 for(mt=0; mt<mtop->nmoltype; mt++) {
1942 do_moltype(fio, &mtop->moltype[mt],bRead,&mtop->symtab,file_version,
1946 if (file_version >= 57) {
1947 gmx_fio_do_int(fio,mtop->nmolblock);
1949 mtop->nmolblock = 1;
1952 snew(mtop->molblock,mtop->nmolblock);
1954 if (file_version >= 57) {
1955 for(mb=0; mb<mtop->nmolblock; mb++) {
1956 do_molblock(fio, &mtop->molblock[mb],bRead,file_version);
1958 gmx_fio_do_int(fio,mtop->natoms);
1960 mtop->molblock[0].type = 0;
1961 mtop->molblock[0].nmol = 1;
1962 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
1963 mtop->molblock[0].nposres_xA = 0;
1964 mtop->molblock[0].nposres_xB = 0;
1967 do_atomtypes (fio, &(mtop->atomtypes),bRead,&(mtop->symtab), file_version);
1969 pr_atomtypes(debug,0,"atomtypes",&mtop->atomtypes,TRUE);
1971 if (file_version < 57) {
1972 /* Debug statements are inside do_idef */
1973 do_idef (fio, &mtop->ffparams,&mtop->moltype[0],bRead,file_version);
1974 mtop->natoms = mtop->moltype[0].atoms.nr;
1977 if(file_version >= 65)
1979 do_cmap(fio, &mtop->ffparams.cmap_grid,bRead);
1983 mtop->ffparams.cmap_grid.ngrid = 0;
1984 mtop->ffparams.cmap_grid.grid_spacing = 0.1;
1985 mtop->ffparams.cmap_grid.cmapdata = NULL;
1988 if (file_version >= 57) {
1989 do_groups(fio, &mtop->groups,bRead,&(mtop->symtab),file_version);
1992 if (file_version < 57) {
1993 do_block(fio, &mtop->moltype[0].cgs,bRead,file_version);
1994 if (bRead && gmx_debug_at) {
1995 pr_block(debug,0,"cgs",&mtop->moltype[0].cgs,TRUE);
1997 do_block(fio, &mtop->mols,bRead,file_version);
1998 /* Add the posres coordinates to the molblock */
1999 add_posres_molblock(mtop);
2002 if (file_version >= 57) {
2003 mtop->mols = mtop_mols(mtop);
2006 pr_block(debug,0,"mols",&mtop->mols,TRUE);
2010 if (file_version < 51) {
2011 /* Here used to be the shake blocks */
2012 do_blocka(fio, &dumb,bRead,file_version);
2020 close_symtab(&(mtop->symtab));
2024 /* If TopOnlyOK is TRUE then we can read even future versions
2025 * of tpx files, provided the file_generation hasn't changed.
2026 * If it is FALSE, we need the inputrecord too, and bail out
2027 * if the file is newer than the program.
2029 * The version and generation if the topology (see top of this file)
2030 * are returned in the two last arguments.
2032 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2034 static void do_tpxheader(t_fileio *fio,gmx_bool bRead,t_tpxheader *tpx,
2035 gmx_bool TopOnlyOK, int *file_version,
2036 int *file_generation)
2045 gmx_fio_checktype(fio);
2046 gmx_fio_setdebug(fio,bDebugMode());
2048 /* NEW! XDR tpb file */
2049 precision = sizeof(real);
2051 gmx_fio_do_string(fio,buf);
2052 if (strncmp(buf,"VERSION",7))
2053 gmx_fatal(FARGS,"Can not read file %s,\n"
2054 " this file is from a Gromacs version which is older than 2.0\n"
2055 " Make a new one with grompp or use a gro or pdb file, if possible",
2056 gmx_fio_getname(fio));
2057 gmx_fio_do_int(fio,precision);
2058 bDouble = (precision == sizeof(double));
2059 if ((precision != sizeof(float)) && !bDouble)
2060 gmx_fatal(FARGS,"Unknown precision in file %s: real is %d bytes "
2061 "instead of %d or %d",
2062 gmx_fio_getname(fio),precision,sizeof(float),sizeof(double));
2063 gmx_fio_setprecision(fio,bDouble);
2064 fprintf(stderr,"Reading file %s, %s (%s precision)\n",
2065 gmx_fio_getname(fio),buf,bDouble ? "double" : "single");
2068 gmx_fio_write_string(fio,GromacsVersion());
2069 bDouble = (precision == sizeof(double));
2070 gmx_fio_setprecision(fio,bDouble);
2071 gmx_fio_do_int(fio,precision);
2073 fgen = tpx_generation;
2076 /* Check versions! */
2077 gmx_fio_do_int(fio,fver);
2080 gmx_fio_do_int(fio,fgen);
2084 if(file_version!=NULL)
2085 *file_version = fver;
2086 if(file_generation!=NULL)
2087 *file_generation = fgen;
2090 if ((fver <= tpx_incompatible_version) ||
2091 ((fver > tpx_version) && !TopOnlyOK) ||
2092 (fgen > tpx_generation))
2093 gmx_fatal(FARGS,"reading tpx file (%s) version %d with version %d program",
2094 gmx_fio_getname(fio),fver,tpx_version);
2096 do_section(fio,eitemHEADER,bRead);
2097 gmx_fio_do_int(fio,tpx->natoms);
2099 gmx_fio_do_int(fio,tpx->ngtc);
2103 gmx_fio_do_int(fio,idum);
2104 gmx_fio_do_real(fio,rdum);
2106 gmx_fio_do_real(fio,tpx->lambda);
2107 gmx_fio_do_int(fio,tpx->bIr);
2108 gmx_fio_do_int(fio,tpx->bTop);
2109 gmx_fio_do_int(fio,tpx->bX);
2110 gmx_fio_do_int(fio,tpx->bV);
2111 gmx_fio_do_int(fio,tpx->bF);
2112 gmx_fio_do_int(fio,tpx->bBox);
2114 if((fgen > tpx_generation)) {
2115 /* This can only happen if TopOnlyOK=TRUE */
2120 static int do_tpx(t_fileio *fio, gmx_bool bRead,
2121 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop,
2122 gmx_bool bXVallocated)
2127 gmx_bool TopOnlyOK,bDum=TRUE;
2128 int file_version,file_generation;
2132 gmx_bool bPeriodicMols;
2135 tpx.natoms = state->natoms;
2136 tpx.ngtc = state->ngtc;
2137 tpx.lambda = state->lambda;
2138 tpx.bIr = (ir != NULL);
2139 tpx.bTop = (mtop != NULL);
2140 tpx.bX = (state->x != NULL);
2141 tpx.bV = (state->v != NULL);
2142 tpx.bF = (f != NULL);
2146 TopOnlyOK = (ir==NULL);
2148 do_tpxheader(fio,bRead,&tpx,TopOnlyOK,&file_version,&file_generation);
2152 state->lambda = tpx.lambda;
2153 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
2157 init_state(state,0,tpx.ngtc,0,0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
2158 state->natoms = tpx.natoms;
2159 state->nalloc = tpx.natoms;
2163 init_state(state,tpx.natoms,tpx.ngtc,0,0); /* nose-hoover chains */
2167 #define do_test(fio,b,p) if (bRead && (p!=NULL) && !b) gmx_fatal(FARGS,"No %s in %s",#p,gmx_fio_getname(fio))
2169 do_test(fio,tpx.bBox,state->box);
2170 do_section(fio,eitemBOX,bRead);
2172 gmx_fio_ndo_rvec(fio,state->box,DIM);
2173 if (file_version >= 51) {
2174 gmx_fio_ndo_rvec(fio,state->box_rel,DIM);
2176 /* We initialize box_rel after reading the inputrec */
2177 clear_mat(state->box_rel);
2179 if (file_version >= 28) {
2180 gmx_fio_ndo_rvec(fio,state->boxv,DIM);
2181 if (file_version < 56) {
2183 gmx_fio_ndo_rvec(fio,mdum,DIM);
2188 if (state->ngtc > 0 && file_version >= 28) {
2190 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
2191 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
2192 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
2193 snew(dumv,state->ngtc);
2194 if (file_version < 69) {
2195 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2197 /* These used to be the Berendsen tcoupl_lambda's */
2198 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2202 /* Prior to tpx version 26, the inputrec was here.
2203 * I moved it to enable partial forward-compatibility
2204 * for analysis/viewer programs.
2206 if(file_version<26) {
2207 do_test(fio,tpx.bIr,ir);
2208 do_section(fio,eitemIR,bRead);
2211 do_inputrec(fio, ir,bRead,file_version,
2212 mtop ? &mtop->ffparams.fudgeQQ : NULL);
2214 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2217 do_inputrec(fio, &dum_ir,bRead,file_version,
2218 mtop ? &mtop->ffparams.fudgeQQ :NULL);
2220 pr_inputrec(debug,0,"inputrec",&dum_ir,FALSE);
2221 done_inputrec(&dum_ir);
2227 do_test(fio,tpx.bTop,mtop);
2228 do_section(fio,eitemTOP,bRead);
2231 do_mtop(fio,mtop,bRead, file_version);
2233 do_mtop(fio,&dum_top,bRead,file_version);
2234 done_mtop(&dum_top,TRUE);
2237 do_test(fio,tpx.bX,state->x);
2238 do_section(fio,eitemX,bRead);
2241 state->flags |= (1<<estX);
2243 gmx_fio_ndo_rvec(fio,state->x,state->natoms);
2246 do_test(fio,tpx.bV,state->v);
2247 do_section(fio,eitemV,bRead);
2250 state->flags |= (1<<estV);
2252 gmx_fio_ndo_rvec(fio,state->v,state->natoms);
2255 do_test(fio,tpx.bF,f);
2256 do_section(fio,eitemF,bRead);
2257 if (tpx.bF) gmx_fio_ndo_rvec(fio,f,state->natoms);
2259 /* Starting with tpx version 26, we have the inputrec
2260 * at the end of the file, so we can ignore it
2261 * if the file is never than the software (but still the
2262 * same generation - see comments at the top of this file.
2267 bPeriodicMols = FALSE;
2268 if (file_version >= 26) {
2269 do_test(fio,tpx.bIr,ir);
2270 do_section(fio,eitemIR,bRead);
2272 if (file_version >= 53) {
2273 /* Removed the pbc info from do_inputrec, since we always want it */
2276 bPeriodicMols = ir->bPeriodicMols;
2278 gmx_fio_do_int(fio,ePBC);
2279 gmx_fio_do_gmx_bool(fio,bPeriodicMols);
2281 if (file_generation <= tpx_generation && ir) {
2282 do_inputrec(fio, ir,bRead,file_version,mtop ? &mtop->ffparams.fudgeQQ : NULL);
2284 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2285 if (file_version < 51)
2286 set_box_rel(ir,state);
2287 if (file_version < 53) {
2289 bPeriodicMols = ir->bPeriodicMols;
2292 if (bRead && ir && file_version >= 53) {
2293 /* We need to do this after do_inputrec, since that initializes ir */
2295 ir->bPeriodicMols = bPeriodicMols;
2304 if (state->ngtc == 0)
2306 /* Reading old version without tcoupl state data: set it */
2307 init_gtc_state(state,ir->opts.ngtc,0,ir->opts.nhchainlength);
2309 if (tpx.bTop && mtop)
2311 if (file_version < 57)
2313 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
2315 ir->eDisre = edrSimple;
2319 ir->eDisre = edrNone;
2322 set_disres_npair(mtop);
2326 if (tpx.bTop && mtop)
2328 gmx_mtop_finalize(mtop);
2331 if (file_version >= 57)
2335 env = getenv("GMX_NOCHARGEGROUPS");
2338 sscanf(env,"%d",&ienv);
2339 fprintf(stderr,"\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
2344 "Will make single atomic charge groups in non-solvent%s\n",
2345 ienv > 1 ? " and solvent" : "");
2346 gmx_mtop_make_atomic_charge_groups(mtop,ienv==1);
2348 fprintf(stderr,"\n");
2356 /************************************************************
2358 * The following routines are the exported ones
2360 ************************************************************/
2362 t_fileio *open_tpx(const char *fn,const char *mode)
2364 return gmx_fio_open(fn,mode);
2367 void close_tpx(t_fileio *fio)
2372 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
2373 int *file_version, int *file_generation)
2377 fio = open_tpx(fn,"r");
2378 do_tpxheader(fio,TRUE,tpx,TopOnlyOK,file_version,file_generation);
2382 void write_tpx_state(const char *fn,
2383 t_inputrec *ir,t_state *state,gmx_mtop_t *mtop)
2387 fio = open_tpx(fn,"w");
2388 do_tpx(fio,FALSE,ir,state,NULL,mtop,FALSE);
2392 void read_tpx_state(const char *fn,
2393 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop)
2397 fio = open_tpx(fn,"r");
2398 do_tpx(fio,TRUE,ir,state,f,mtop,FALSE);
2402 int read_tpx(const char *fn,
2403 t_inputrec *ir, matrix box,int *natoms,
2404 rvec *x,rvec *v,rvec *f,gmx_mtop_t *mtop)
2412 fio = open_tpx(fn,"r");
2413 ePBC = do_tpx(fio,TRUE,ir,&state,f,mtop,TRUE);
2415 *natoms = state.natoms;
2417 copy_mat(state.box,box);
2425 int read_tpx_top(const char *fn,
2426 t_inputrec *ir, matrix box,int *natoms,
2427 rvec *x,rvec *v,rvec *f,t_topology *top)
2433 ePBC = read_tpx(fn,ir,box,natoms,x,v,f,&mtop);
2435 *top = gmx_mtop_t_to_t_topology(&mtop);
2440 gmx_bool fn2bTPX(const char *file)
2442 switch (fn2ftp(file)) {
2452 gmx_bool read_tps_conf(const char *infile,char *title,t_topology *top,int *ePBC,
2453 rvec **x,rvec **v,matrix box,gmx_bool bMass)
2456 int natoms,i,version,generation;
2457 gmx_bool bTop,bXNULL;
2459 t_topology *topconv;
2462 bTop = fn2bTPX(infile);
2465 read_tpxheader(infile,&header,TRUE,&version,&generation);
2467 snew(*x,header.natoms);
2469 snew(*v,header.natoms);
2471 *ePBC = read_tpx(infile,NULL,box,&natoms,
2472 (x==NULL) ? NULL : *x,(v==NULL) ? NULL : *v,NULL,mtop);
2473 *top = gmx_mtop_t_to_t_topology(mtop);
2475 strcpy(title,*top->name);
2476 tpx_make_chain_identifiers(&top->atoms,&top->mols);
2479 get_stx_coordnum(infile,&natoms);
2480 init_t_atoms(&top->atoms,natoms,FALSE);
2481 bXNULL = (x == NULL);
2485 read_stx_conf(infile,title,&top->atoms,*x,(v==NULL) ? NULL : *v,ePBC,box);
2491 aps = gmx_atomprop_init();
2492 for(i=0; (i<natoms); i++)
2493 if (!gmx_atomprop_query(aps,epropMass,
2494 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2495 *top->atoms.atomname[i],
2496 &(top->atoms.atom[i].m))) {
2498 fprintf(debug,"Can not find mass for atom %s %d %s, setting to 1\n",
2499 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2500 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
2501 *top->atoms.atomname[i]);
2503 gmx_atomprop_destroy(aps);
2505 top->idef.ntypes=-1;